! #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr ! #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if #define IF_ERROR_RETURN(action) if (status >0) then; TRACEBACK; action; return; end if ! #define PRISM_ERR write (gol,'("from call to PRISM routine")'); call goErr #define IF_PRISM_NOTOK_RETURN(action) if (status/=OASIS_OK) then; PRISM_ERR; TRACEBACK; action; return; end if ! #include "tm5.inc" ! !----------------------------------------------------------------------------- ! TM5 ! !----------------------------------------------------------------------------- !BOP ! ! !MODULE: TM5_PRISM ! ! !DESCRIPTION: !\\ !\\ ! !INTERFACE: ! MODULE TM5_PRISM ! ! !USES: ! USE MOD_OASIS USE GO, ONLY : gol, goPr, goErr, TDate IMPLICIT NONE PRIVATE ! ! !PUBLIC MEMBER FUNCTIONS: ! public :: TM5_Prism_Init ! read control parameters from rc file public :: TM5_Prism_Init2 ! prism init: grids, partition, coupled variables public :: TM5_Prism_Done ! deallocate public :: SetPrismTime ! current time/date into prism format public :: InqCplVar ! check if a coupled variable exists public :: Init, Done, Setup, Remap ! methods for TshRemap object ! ! !PUBLIC DATA MEMBERS: ! character(len=6), public, parameter :: comp_name = 'ctm5mp' integer, public :: PRISM_start_date(6) ! prism reference start time integer, public :: comp_id ! tm5 id for prism integer, public :: ifs_cpl_freq ! ifs coupling frequency in hours integer, public :: lpj_cpl_freq ! lpjg coupling frequency in hours integer, public :: pis_cpl_freq ! pisces coupling frequency in hours integer, public :: ifs_cpl_nlev ! number of levels for coupling with IFS integer, public :: ifs_cpl_nlev_cutoff ! reduced number of levels for coupling with IFS, ! applied for aerosols integer, public :: ifs_shT ! ifs spectral fields : resolution integer, public :: ifs_shn ! ifs spectral fields : number of coeff ! logical, public :: refine_levels ! do we send/receive all IFS levels? (then we need to "refine them" here) logical, public :: coupled_to_lpj ! is TM5 coupled to LPJ-Guess? logical, public :: coupled_to_pis ! is TM5 coupled to PISCES? ! ! !PRIVATE DATA MEMBERS: ! character(len=*), parameter :: mname = 'TM5_Prism' integer, parameter :: len_grid_name = 4 ! oasis3: fixed length for gridname ! integer :: ifs_npes ! ifs number of procs integer, parameter :: wp = SELECTED_REAL_KIND(12,307) ! working precision = double character(len=256) :: error_message ! ! !PUBLIC TYPES: ! TYPE, PUBLIC :: TshRemap ! ! remapping: for each element in received grid, identify target indices: ! ! remap(:,:,1) : 1 = real part, 2 = imag part ! remap(:,:,2) : index in triangle : 1=(0,0) 2=(0,1) 3=(0,2) ... np=(m,n) ! remap(:,:,3) : level ! ! Info array has values : m*1000 + n + lev/100 and negative for imaginary part ! ! receive spectral info every timestep to avoid memory growth type(TDate) :: t ! store time of current info integer, pointer :: remap(:,:,:) ! remapping indices integer :: shT ! truncation END TYPE TshRemap TYPE, PUBLIC :: TCplVar ! --- TM Coupled Variable --- character(len=128) :: name ! tm5 internal name character(len=128) :: cpl_name ! short name used in rcfile integer, pointer :: var_id(:) ! list of id return by oasis3 (one per level) character(len=128), pointer :: var_name(:)! list of names used in the namcouple (one per level) logical :: serial ! serial transfer character(len=3) :: intent ! in or out character(len=2) :: grid_type ! spectral or gridpoint integer :: region ! integer :: nlev ! real :: west_deg, dlon_deg, south_deg, dlat_deg ! lon/lat grids integer :: nlon, nlat integer :: shT, shn, shn_recv ! spectral info integer :: itr ! tracer index if any logical :: cache ! cache type(TDate) :: cache_tmid real, pointer :: cache_data(:,:,:) => null() END TYPE TCplVar ! --- var ----------------------------------- integer, parameter :: maxcplvar = 165 ! max nb of coupled fields type(TCplVar), public, save :: CplVar(maxcplvar) ! array of coupled fields integer, public :: ncplvar ! actual nb of coupled fields integer :: region_glb, region_sfc character(len=1024) :: prism_get_list ! comma seperated lists of coupled fields character(len=1024) :: prism_put_list integer, allocatable :: part_id(:) integer, allocatable :: var_shape(:,:) integer :: sp_part_id integer, allocatable :: sp_var_shape(:) ! ! !INTERFACE: ! interface Init module procedure shremap_Init end interface interface Done module procedure shremap_Done end interface interface Setup module procedure shremap_Setup end interface interface Remap module procedure shremap_Remap end interface interface SetPrismTime module procedure SetPrismTime_date end interface ! ! !REVISION HISTORY: ! 10 Sep 2013 - Ph. Le Sager - cleanup, document, remove oasis4 stuff (can ! always be retrieved with svn if needed) ! 8 Oct 2013 - Ph. Le Sager - dummy CO2 exchange with LPJ-Guess ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ CONTAINS !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: TM5_PRISM_INIT ! ! !DESCRIPTION: read coupling information from rc file. !\\ !\\ ! !INTERFACE: ! SUBROUTINE TM5_Prism_Init( rcfile, status ) ! ! !USES: ! use GO, only : TrcFile, Init, Done, ReadRc ! ! !INPUT PARAMETERS: ! character(len=*), intent(in) :: rcfile ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/TM5_Prism_Init' type(TrcFile) :: rcF ! --- begin ----------------------------------- ! * settings from rcfile call Init( rcF, rcfile, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'ifs.cpl.nlev', ifs_cpl_nlev, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcf, 'ifs.cpl.nlev.cutoff', ifs_cpl_nlev_cutoff, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'ifs.shresol', ifs_shT, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'cpl.ifs.period', ifs_cpl_freq, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'cpl.lpj.period', lpj_cpl_freq, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'cpl.pis.period', pis_cpl_freq, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'prism.get', prism_get_list, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'prism.put', prism_put_list, status ) IF_NOTOK_RETURN(status=1) call ReadRc( rcF, 'coupled_to_lpjguess', coupled_to_lpj, status, default=.false. ) IF_ERROR_RETURN(status=1) call ReadRc( rcF, 'coupled_to_pisces', coupled_to_pis, status, default=.false. ) IF_ERROR_RETURN(status=1) select case (ifs_cpl_nlev) case (91,62) refine_levels=.true. case(34,31,10,4) refine_levels=.false. case default write(gol,*) " Wrong (sub)set of levels is exchanged between IFS and TM5 " ; call goErr write(gol,*) " Either (4, 10 or 34 out of) 91, or (31 out of) 62" ; call goErr status=1 TRACEBACK; return end select call Done( rcF, status ) IF_NOTOK_RETURN(status=1) ! * spectral grids ifs_shn = (ifs_shT+1)*(ifs_shT+2)/2 ! number of coeff ! check that compilation was perform as expected with respect to optical feedback #ifndef with_ecearth_optics if (index(prism_put_list,'AOD_01') /=0) then write(gol,*) "Feedback of aerosols optical properties requires 'with_ecearth_optics' cpp"; call goErr write(gol,*) "You must recompile TM5-MP with cpp defs 'with_ecearth_optics'"; call goErr write(gol,*) "This can be done in the config-build.xml file with the TM5_MDEFS_FFLAGS key."; call goErr status=1 TRACEBACK; return endif #endif status = 0 END SUBROUTINE TM5_PRISM_INIT !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: TM5_PRISM_INIT2 ! ! !DESCRIPTION: prism grid writing, prism partition, init coupled variables !\\ !\\ ! !INTERFACE: ! subroutine TM5_Prism_Init2( nreg_all, nreg, ireg_sfc, lli, levi, status ) ! ! !USES: ! use Grid, only : TllGridInfo, TLevelInfo use Grid, only : TshGridInfo, Init, Done #ifdef parallel_cplng use tm5_distgrid, only : dgrid, Get_DistGrid #endif use GO, only : NewDate, goReadFromLine use dims, only : lm use chem_param, only : names, io3, ich4, ico2 use partools, only : isRoot #ifdef with_m7 use chem_param, only : inus_n, iso4nus, iais_n, iso4ais, ibcais, ipomais use chem_param, only : iacs_n, iso4acs, ibcacs, ipomacs, issacs, iduacs use chem_param, only : icos_n, iso4cos, ibccos, ipomcos, isscos, iducos use chem_param, only : iaii_n, ibcaii, ipomaii, iaci_n, iduaci, icoi_n, iducoi use chem_param, only : ino3_a, imsa #endif ! ! !INPUT PARAMETERS: ! integer, intent(in) :: nreg_all integer, intent(in) :: nreg integer, intent(in) :: ireg_sfc type(TllGridInfo), intent(in) :: lli(1:nreg_all) type(TLevelInfo), intent(in) :: levi ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 8 Oct 2013 - Ph. Le Sager - coupling with LPJG ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/TM5_Prism_Init2' integer :: ireg, region, id, ivar, ilev, i1, j1, i2, j2 integer, parameter :: nc = 4 ! corners integer :: nx, ny integer :: i, j, k, m, n, z character(len=len_grid_name) :: point_name real(ip_realwp_p), allocatable :: lons(:,:), lats(:,:) real(ip_realwp_p), allocatable :: clons(:,:,:), clats(:,:,:) real(ip_realwp_p), allocatable :: area(:,:) integer, allocatable :: mask(:,:) character(len=len_grid_name) :: sp_point_name real(ip_realwp_p), allocatable :: sp_lons(:,:), sp_lats(:,:) integer, allocatable :: sp_mask(:,:) real(ip_realwp_p) :: sp_dlon, sp_dlat #ifdef parallel_cplng integer :: part_val(1:5) #else integer :: part_val(1:3) #endif integer :: sp_part_val(1:3) character(len=128) :: cpl_name integer :: var_nodims(2) integer :: var_type integer(kind=ip_intwp_p) :: var_intent logical :: write_grid #ifdef parallel_cplng type(TllGridInfo) :: local_lli ! local Lat-Lon grid info #endif ! --- begin ----------------------------------- write (gol,'("initialize prism coupling:")'); call goPr write (gol,'(" component : ",a)') trim(comp_name); call goPr ! store in module variables region_glb = 1 region_sfc = ireg_sfc ! storage for variable shape: allocate( part_id(nreg_all) ) allocate( var_shape(4,nreg_all) ) allocate( sp_var_shape(4) ) ! init to zero on all pe's part_id = 0 var_shape = 0 sp_part_id = 0 sp_var_shape = 0 ! ============ oasis3 grid writing ================= write_grid=.false. !! HARDCODED !! ! Define the grids by master proc if ( isroot .and. write_grid ) then call oasis_start_grids_writing( status ) ! **** lon/lat grid do region = 1, nreg_all ! name of grid points if ( region == region_glb ) then ! global region point_name = 'CTM3' else if ( region == region_sfc ) then ! global surface grid: point_name = 'CTM1' else ! global grids only ... cycle end if write (gol,'(" define points ",a," ...")') point_name; call goPr ! grid size: nx = lli(region)%nlon ny = lli(region)%nlat write (gol,'(" region : ",i6)') region; call goPr write (gol,'(" size : ",2i6)') nx, ny; call goPr allocate( lons(nx,ny) ) allocate( lats(nx,ny) ) allocate( clons(nx,ny,nc) ) allocate( clats(nx,ny,nc) ) allocate( area(nx,ny) ) allocate( mask(nx,ny) ) ! set lon/lat grid (grid cell centers): do i = 1, nx lons(i,:) = lli(region)%lon_deg(i) end do do j = 1, ny lats(:,j) = lli(region)%lat_deg(j) end do call oasis_write_grid( point_name, nx, ny, lons, lats ) ! set corner lon/lat: ! 3 o o 2 ! 4 o o 1 do i = 1, nx clons(i,:,1) = lli(region)%blon_deg(i ) clons(i,:,2) = lli(region)%blon_deg(i ) clons(i,:,3) = lli(region)%blon_deg(i-1) clons(i,:,4) = lli(region)%blon_deg(i-1) end do do j = 1, ny clats(:,j,1) = lli(region)%blat_deg(j-1) clats(:,j,2) = lli(region)%blat_deg(j ) clats(:,j,3) = lli(region)%blat_deg(j ) clats(:,j,4) = lli(region)%blat_deg(j-1) end do write (gol,'(" write corners ...")'); call goPr call oasis_write_corner( point_name, nx, ny, nc, clons, clats ) ! land/sea mask mask = 0 ! not masked; gives warnings about 'sea-world' cplout ... write (gol,'(" write mask ...")'); call goPr call oasis_write_mask( point_name, nx, ny, mask ) do j = 1, ny area(:,j) = lli(region)%area_m2(j) end do write (gol,'(" write area ...")'); call goPr call oasis_write_area( point_name, nx, ny, area ) deallocate( lons ) deallocate( lats ) deallocate( clons ) deallocate( clats ) deallocate( area ) deallocate( mask ) end do ! regions ! *** SPECTRAL GRID write(sp_point_name, '("C",i3.3)') ifs_shT write (gol,'(" define points ",a," ...")') trim(sp_point_name); call goPr allocate( sp_lons(2*ifs_shn,1) ) allocate( sp_lats(2*ifs_shn,1) ) allocate( sp_mask(2*ifs_shn,1) ) ! Triangular storage: ! ! NSMAX * * .. * ! : ! ! 1 * * ! n 0 * ! ! 0 1 .. NSMAX ! m "wavenumber" ! ! dummy locations : (n*2+z+0.5)*dlon, (m+0.5)*dlat ! where z=0 is real part and z=1 is imaginary part ! dummy spacing: sp_dlon = 0.1 ! degree sp_dlat = 0.1 ! degree ! index in coeff array: k = 0 ! loop over global wavenumbers: do m = 0, ifs_shT ! loop over complex coeff: do n = m, ifs_shT ! loop over real/complex do z = 0, 1 ! next coeff: k = k + 1 ! cell centers: sp_lons(k,1) = -180.0 + (n*2+z+0.5) * sp_dlon sp_lats(k,1) = -90.0 + (m +0.5) * sp_dlat ! no mask: sp_mask(k,1) = 0 ! not masked end do ! re,im end do ! n end do ! m call oasis_write_grid( sp_point_name, 2*ifs_shn, 1, sp_lons, sp_lats ) call oasis_write_mask( sp_point_name, 2*ifs_shn, 1, sp_mask ) deallocate( sp_lons ) deallocate( sp_lats ) deallocate( sp_mask ) call oasis_terminate_grids_writing() else write (gol,'(" not necessary to write grids (not root) ...")'); call goPr end if ! root and write_grid ! ============ Partition ================= write (gol,'(" define partitions ...")'); call goPr ! *** LAT/LON reg: do region = 1, nreg_all if ( (region /= region_glb) .and. (region /= region_sfc) ) cycle ! global grids only nx = lli(region)%nlon ny = lli(region)%nlat #ifdef parallel_cplng call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2 , J_STRT=j1, J_STOP=j2 ) ! store shape: var_shape(1:2,region) = (/1,i2-i1+1/) var_shape(3:4,region) = (/1,j2-j1+1/) part_val(1) = 2 part_val(2) = i1-1+(j1-1)*nx part_val(3) = i2-i1+1 part_val(4) = j2-j1+1 part_val(5) = nx #else ! store shape: var_shape(1:2,region) = (/1,nx/) var_shape(3:4,region) = (/1,ny/) part_val(1) = 0 ! serial partition part_val(2) = 0 part_val(3) = 0 if ( isroot ) part_val(3) = nx*ny ! total grid size #endif status = OASIS_OK ! <-- status was not set by PRISM_Def_Partition_Proto (is it still true in OASIS3-MCT?) call oasis_def_partition( part_id(region), part_val, status ) if (status/=OASIS_OK) then write (error_message,'("from OASIS_DEF_PARTITION : ",i6)') status TRACEBACK call oasis_abort( comp_id, rname, error_message ) end if #ifdef parallel_cplng write (gol,'(" partition : ",i6," ; ",5i6)') part_id(region), part_val; call goPr #else write (gol,'(" partition : ",i6," ; ",3i6)') part_id(region), part_val; call goPr #endif end do reg ! *** SPECTRAL sp_part_val(1) = 0 sp_part_val(2) = 0 sp_part_val(3) = 0 if ( isroot ) sp_part_val(3) = ifs_shn*2 ! total grid size !status = OASIS_OK ! <-- status was not set by PRISM_Def_Partition_Proto (is it still true in OASIS3-MCT?) call oasis_def_partition( sp_part_id, sp_part_val, status ) if (status/=OASIS_OK) then write (error_message,'("from OASIS_DEF_PARTITION : ",i6)') status TRACEBACK call oasis_abort( comp_id, rname, error_message ) end if write (gol,'(" partition : ",i6," ; ",3i6)') sp_part_id, sp_part_val; call goPr ! ------------------------------------------------------------------------- ! * CONFIGURE COUPLING FIELDS ! ------------------------------------------------------------------------- ! ! (0) DEFAULT ! write (gol,'(" init cplvar list ...")'); call goPr var_nodims(1) = 2 ! rank of coupling field var_nodims(2) = 1 ! number of bundles in coupling field (always 1) var_type = OASIS_Real do ivar = 1, size(CplVar) CplVar(ivar)%cpl_name = '' nullify( CplVar(ivar)%var_id ) nullify( CplVar(ivar)%var_name ) CplVar(ivar)%grid_type = 'll' ! lon/lat CplVar(ivar)%region = -1 CplVar(ivar)%nlev = -1 #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'xxx' CplVar(ivar)%cache = .false. ! at init, flag means "will use cache" CplVar(ivar)%cache_tmid = NewDate(year=9999) nullify( CplVar(ivar)%cache_data ) end do ! Above init should be same as: ! CplVar = TCplVar('','',null(),null(),.true.,'xxx','ll',-1,-1,0.,0.,0.,0.,0,0,0,0,0,.false.,NewDate(year=9999),null()) ! We could also just set pointers to => null() in declaration l.117 ncplvar = 0 ! no fields defined yet ! ! (1) METEO VARIABLES ! write (gol,'(" set meteo cplvars ...")'); call goPr write (gol,'(" list : ",a)') trim(prism_get_list); call goPr ivar = 0 GET: DO if ( len_trim(prism_get_list) == 0 ) exit ! leave if empty ! extract next name from list call goReadFromLine( prism_get_list, cpl_name, status, sep=',' ) IF_NOTOK_RETURN(status=1) write (gol,'(" extracted ",a," ...")') trim(cpl_name); call goPr ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = trim(cpl_name) CplVar(ivar)%intent = 'in' select case ( trim(cpl_name) ) case ( 'LNSP' ) CplVar(ivar)%name = 'LNSP' CplVar(ivar)%grid_type = 'sh' ! spherical harmonics CplVar(ivar)%serial = .true. CplVar(ivar)%nlev = 1 CplVar(ivar)%cache = .true. case ( 'VOR' ) CplVar(ivar)%name = 'VO' CplVar(ivar)%grid_type = 'sh' ! spherical harmonics CplVar(ivar)%serial = .true. CplVar(ivar)%nlev = ifs_cpl_nlev CplVar(ivar)%cache = .true. case ( 'DIV' ) CplVar(ivar)%name = 'D' CplVar(ivar)%grid_type = 'sh' ! spherical harmonics CplVar(ivar)%serial = .true. CplVar(ivar)%nlev = ifs_cpl_nlev CplVar(ivar)%cache = .true. case ( 'SPRES' ) CplVar(ivar)%name = 'sp' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 CplVar(ivar)%cache = .true. case ( 'TMP' ) CplVar(ivar)%name = 'T' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'HUM' ) CplVar(ivar)%name = 'Q' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'OROG' ) CplVar(ivar)%name = 'oro' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'CLW' ) CplVar(ivar)%name = 'CLWC' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'CIW') CplVar(ivar)%name = 'CIWC' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'CC' ) CplVar(ivar)%name = 'CC' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'CCO' ) CplVar(ivar)%name = 'CCO' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'CCU' ) CplVar(ivar)%name = 'CCU' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'UMF' ) CplVar(ivar)%name = 'UDMF' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'DMF' ) CplVar(ivar)%name = 'DDMF' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'UDR' ) CplVar(ivar)%name = 'UDDR' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'DDR' ) CplVar(ivar)%name = 'DDDR' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = ifs_cpl_nlev case ( 'LSMSK' ) CplVar(ivar)%name = 'lsm' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'ALB' ) CplVar(ivar)%name = 'albedo' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SR' ) CplVar(ivar)%name = 'sr' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'CI' ) CplVar(ivar)%name = 'ci' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SST' ) CplVar(ivar)%name = 'sst' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'U10M' ) CplVar(ivar)%name = 'u10m' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'V10M' ) CplVar(ivar)%name = 'v10m' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'WSPD' ) CplVar(ivar)%name = 'wspd' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SRC' ) CplVar(ivar)%name = 'src' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'D2M' ) CplVar(ivar)%name = 'd2m' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'T2M' ) CplVar(ivar)%name = 't2m' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SSHF' ) CplVar(ivar)%name = 'sshf' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SLHF' ) CplVar(ivar)%name = 'slhf' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'EWSS' ) CplVar(ivar)%name = 'ewss' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'NSSS' ) CplVar(ivar)%name = 'nsss' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'CP' ) CplVar(ivar)%name = 'cp' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'LSP' ) CplVar(ivar)%name = 'lsp' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SSR' ) CplVar(ivar)%name = 'ssr' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SKT__') CplVar(ivar)%name = 'skt' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SF___' ) CplVar(ivar)%name = 'sf' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SD' ) CplVar(ivar)%name = 'sd' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'SWVL1' ) CplVar(ivar)%name = 'swvl1' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV01' ) CplVar(ivar)%name = 'tv01' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV02' ) CplVar(ivar)%name = 'tv02' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV03' ) CplVar(ivar)%name = 'tv03' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV04' ) CplVar(ivar)%name = 'tv04' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV05' ) CplVar(ivar)%name = 'tv05' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV06' ) CplVar(ivar)%name = 'tv06' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV07' ) CplVar(ivar)%name = 'tv07' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV09' ) CplVar(ivar)%name = 'tv09' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV10' ) CplVar(ivar)%name = 'tv10' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV11' ) CplVar(ivar)%name = 'tv11' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV13' ) CplVar(ivar)%name = 'tv13' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV16' ) CplVar(ivar)%name = 'tv16' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV17' ) CplVar(ivar)%name = 'tv17' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV18' ) CplVar(ivar)%name = 'tv18' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'TV19' ) CplVar(ivar)%name = 'tv19' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'CVL' ) CplVar(ivar)%name = 'cvl' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case ( 'CVH' ) CplVar(ivar)%name = 'cvh' CplVar(ivar)%region = region_sfc CplVar(ivar)%nlev = 1 case default write (gol,'("unsupported cpl_name : ",a)') trim(cpl_name); call goErr TRACEBACK; status=1; return end select ! check if ( CplVar(ivar)%nlev < 1 ) then write (gol,'("found nlev ",i6," in cplvar ",i4," (",a,")")') CplVar(ivar)%nlev, ivar, trim(cpl_name); call goErr TRACEBACK; status=1; return end if ! storage per level allocate( CplVar(ivar)%var_id (CplVar(ivar)%nlev) ) allocate( CplVar(ivar)%var_name(CplVar(ivar)%nlev) ) ! name used in namcouple if ( CplVar(ivar)%nlev == 1 ) then ilev = 1 write (CplVar(ivar)%var_name(ilev),'("C_",a)') trim(CplVar(ivar)%cpl_name) write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr else do ilev = 1, CplVar(ivar)%nlev write (CplVar(ivar)%var_name(ilev),'("C_",a,".L",i3.3)') trim(CplVar(ivar)%cpl_name), ilev write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr end do end if ! store latest number: ncplvar = ivar END DO GET ! list coupled var to get ! ! (2) VARIABLES SENT TO IFS ! write (gol,'(" set concentration/optics coupled vars ...")'); call goPr write (gol,'(" list : ",a)') trim(prism_put_list); call goPr ivar = ncplvar PUT: DO if ( len_trim(prism_put_list) == 0 ) exit ! extract next name from list call goReadFromLine( prism_put_list, cpl_name, status, sep=',' ) IF_NOTOK_RETURN(status=1) write (gol,'(" extracted ",a," ...")') trim(cpl_name); call goPr ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = trim(cpl_name) #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'out' CplVar(ivar)%region = region_glb select case ( trim(cpl_name) ) case ( 'O3', 'CH4', 'CO2' ) ! send whole atmosphere for ozone, methane and CO2 if (.not.refine_levels) then ! in both cases this should be ifs_cpl_nlev! CplVar(ivar)%nlev = lm(region_glb) else CplVar(ivar)%nlev = ifs_cpl_nlev endif case default ! for aerosols, do not send all levels in stratosphere ! this works, for refine_levels true or false CplVar(ivar)%nlev = ifs_cpl_nlev_cutoff end select select case ( trim(cpl_name) ) case ( 'CO2' ) CplVar(ivar)%itr = ico2 case ( 'O3' ) CplVar(ivar)%itr = io3 case ( 'CH4' ) CplVar(ivar)%itr = ich4 #ifdef with_m7 case ( 'N1' ) CplVar(ivar)%itr = inus_n case ( 'SU1' ) CplVar(ivar)%itr = iso4nus case ( 'N2' ) CplVar(ivar)%itr = iais_n case ( 'SU2' ) CplVar(ivar)%itr = iso4ais case ( 'BC2' ) CplVar(ivar)%itr = ibcais case ( 'OC2' ) CplVar(ivar)%itr = ipomais case ( 'N3' ) CplVar(ivar)%itr = iacs_n case ( 'SU3' ) CplVar(ivar)%itr = iso4acs case ( 'BC3' ) CplVar(ivar)%itr = ibcacs case ( 'OC3' ) CplVar(ivar)%itr = ipomacs case ( 'SS3' ) CplVar(ivar)%itr = issacs case ( 'DU3' ) CplVar(ivar)%itr = iduacs case ( 'N4' ) CplVar(ivar)%itr = icos_n case ( 'SU4' ) CplVar(ivar)%itr = iso4cos case ( 'BC4' ) CplVar(ivar)%itr = ibccos case ( 'OC4' ) CplVar(ivar)%itr = ipomcos case ( 'SS4' ) CplVar(ivar)%itr = isscos case ( 'DU4' ) CplVar(ivar)%itr = iducos case ( 'N5' ) CplVar(ivar)%itr = iaii_n case ( 'BC5' ) CplVar(ivar)%itr = ibcaii case ( 'OC5' ) CplVar(ivar)%itr = ipomaii case ( 'N6' ) CplVar(ivar)%itr = iaci_n case ( 'DU6' ) CplVar(ivar)%itr = iduaci case ( 'N7' ) CplVar(ivar)%itr = icoi_n case ( 'DU7' ) CplVar(ivar)%itr = iducoi case ( 'NO3' ) CplVar(ivar)%itr = ino3_a case ( 'MSA' ) CplVar(ivar)%itr = imsa #endif #ifdef with_ecearth_optics case ( 'AOD_01', 'AOD_02', 'AOD_03', 'AOD_04', 'AOD_05', 'AOD_06', 'AOD_07', & 'AOD_08', 'AOD_09', 'AOD_10', 'AOD_11', 'AOD_12', 'AOD_13', 'AOD_14', & 'SSA_01', 'SSA_02', 'SSA_03', 'SSA_04', 'SSA_05', 'SSA_06', 'SSA_07', & 'SSA_08', 'SSA_09', 'SSA_10', 'SSA_11', 'SSA_12', 'SSA_13', 'SSA_14', & 'ASF_01', 'ASF_02', 'ASF_03', 'ASF_04', 'ASF_05', 'ASF_06', 'ASF_07', & 'ASF_08', 'ASF_09', 'ASF_10', 'ASF_11', 'ASF_12', 'ASF_13', 'ASF_14' ) #endif case default write (gol,'("unsupported cpl_name : ",a)') trim(cpl_name); call goErr TRACEBACK; status=1; return end select ! set name: select case ( trim(cpl_name) ) #ifdef with_ecearth_optics case ( 'AOD_01', 'AOD_02', 'AOD_03', 'AOD_04', 'AOD_05', 'AOD_06', 'AOD_07', & 'AOD_08', 'AOD_09', 'AOD_10', 'AOD_11', 'AOD_12', 'AOD_13', 'AOD_14', & 'SSA_01', 'SSA_02', 'SSA_03', 'SSA_04', 'SSA_05', 'SSA_06', 'SSA_07', & 'SSA_08', 'SSA_09', 'SSA_10', 'SSA_11', 'SSA_12', 'SSA_13', 'SSA_14', & 'ASF_01', 'ASF_02', 'ASF_03', 'ASF_04', 'ASF_05', 'ASF_06', 'ASF_07', & 'ASF_08', 'ASF_09', 'ASF_10', 'ASF_11', 'ASF_12', 'ASF_13', 'ASF_14' ) write(CplVar(ivar)%name,'(a)') trim(cpl_name) #endif case default CplVar(ivar)%name = names(CplVar(ivar)%itr) end select ! check .. if ( CplVar(ivar)%nlev < 1 ) then write (gol,'("found nlev ",i6," in cplvar ",i4," (",a,")")') CplVar(ivar)%nlev, ivar, trim(cpl_name); call goErr TRACEBACK; status=1; return end if ! storage per level allocate( CplVar(ivar)%var_id (CplVar(ivar)%nlev) ) allocate( CplVar(ivar)%var_name(CplVar(ivar)%nlev) ) ! name used in namcouple if ( CplVar(ivar)%nlev == 1 ) then ilev = 1 write (CplVar(ivar)%var_name(ilev),'(a,"TM5")') trim(CplVar(ivar)%cpl_name) write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr else do ilev = 1, CplVar(ivar)%nlev write (CplVar(ivar)%var_name(ilev),'("C_",a,".L",i3.3)') trim(CplVar(ivar)%cpl_name), ilev end do write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev-1, trim(CplVar(ivar)%var_name(ilev-1)); call goPr end if ncplvar = ivar END DO PUT ! list with coupled names of var to send to IFS ! ! (3) COUPLING WITH LPJG ! if (coupled_to_lpj) then ! Sent to LPJ-Guess ivar = ncplvar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'co2_for_lpjg' CplVar(ivar)%itr = ico2 CplVar(ivar)%name = 'ppmCO2' ! Reserve names(CplVar(ivar)%itr) for tracer mass #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'out' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "LCO2_TM5" ! Land CO2 write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ! Received from LPJ-Guess ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'land_c_flux_nat' CplVar(ivar)%name = 'land_c_flux_nat' #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'in' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 CplVar(ivar)%itr = ico2 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "TM5_LandCNAT" write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'land_c_flux_ant' CplVar(ivar)%name = 'land_c_flux_ant' #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'in' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 CplVar(ivar)%itr = ico2 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "TM5_LandCANT" write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'land_c_npp' CplVar(ivar)%name = 'land_c_npp' #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'in' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 CplVar(ivar)%itr = 999 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "TM5_LandCNPP" write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ncplvar = ivar end if ! ! (4) COUPLING WITH PISCES ! if (coupled_to_pis) then ! Sent to PISCES ivar = ncplvar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'co2_for_pis' CplVar(ivar)%itr = ico2 CplVar(ivar)%name = 'ppmCO2' ! Reserve names(CplVar(ivar)%itr) for tracer mass #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'out' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "OCO2_TM5" ! Ocean CO2 write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ! Received from PISCES ivar = ivar + 1 if ( ivar > maxcplvar ) then write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr TRACEBACK; status=1; return end if CplVar(ivar)%cpl_name = 'oce_c_flux' ! C fluxes from Ocean CplVar(ivar)%name = 'oce_c_flux' #ifdef parallel_cplng CplVar(ivar)%serial = .false. #else CplVar(ivar)%serial = .true. #endif CplVar(ivar)%intent = 'in' CplVar(ivar)%region = region_glb CplVar(ivar)%nlev = 1 CplVar(ivar)%itr = 999 allocate( CplVar(ivar)%var_id (1) ) allocate( CplVar(ivar)%var_name(1) ) CplVar(ivar)%var_name(1) = "TM5_OceCFLX" write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr ncplvar = ivar end if ! ! (5) DEFINE OASIS VARIABLES ! write (gol,'(" define oasis variables ...")'); call goPr do ivar = 1, ncplvar ireg = CplVar(ivar)%region select case ( CplVar(ivar)%intent ) case ( 'in' ) var_intent = OASIS_In case ( 'out' ) var_intent = OASIS_Out case default write (gol,'("unsupported intent : ",a)') trim(CplVar(ivar)%intent); call goErr TRACEBACK; status=1; return end select do ilev = 1, CplVar(ivar)%nlev select case ( CplVar(ivar)%grid_type ) case ( 'sh' ) !DBG write (gol,'(" ",i4," (",i3,") spectral variable ",a," ...")') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr !DBG write (gol,'(" region : ",i6)' ) ireg ; call goPr !DBG write (gol,'(" partition : ",i6)' ) sp_part_id ; call goPr !DBG write (gol,'(" nodims : ",2i6)' ) var_nodims ; call goPr !DBG write (gol,'(" shape : ",4i6)' ) sp_var_shape(1:4) ; call goPr call oasis_def_var( & CplVar(ivar)%var_id(ilev), & trim(CplVar(ivar)%var_name(ilev)), & sp_part_id, var_nodims, var_intent, & sp_var_shape(1:4), var_type, status ) IF_PRISM_NOTOK_RETURN(status=1) case ( 'll' ) !DBG write (gol,'(" ",i4," (",i3,") gridded variable ",a," ...")') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr !DBG write (gol,'(" region : ",i6)' ) ireg ; call goPr !DBG write (gol,'(" partition : ",i6)' ) part_id(ireg) ; call goPr !DBG write (gol,'(" nodims : ",2i6)' ) var_nodims ; call goPr !DBG write (gol,'(" shape : ",4i6)' ) var_shape(1:4,ireg) ; call goPr call oasis_def_var( & CplVar(ivar)%var_id(ilev), & trim(CplVar(ivar)%var_name(ilev)), & part_id(ireg), var_nodims, var_intent, & var_shape(1:4,ireg), var_type, status ) IF_PRISM_NOTOK_RETURN(status=1) case default write (gol,'("unsupported grid_type : ",a)') trim(CplVar(ivar)%grid_type); call goErr TRACEBACK; status=1; return end select end do ! levels end do ! var ! ! (6) STORE GRID INFO ! write (gol,'("add grid info to cplvars ...")'); call goPr do ivar = 1, ncplvar do ilev = 1, CplVar(ivar)%nlev select case ( CplVar(ivar)%grid_type ) case ( 'sh' ) CplVar(ivar)%shT = ifs_shT CplVar(ivar)%shn = ifs_shn CplVar(ivar)%shn_recv = ifs_shn if ( CplVar(ivar)%cache ) allocate( CplVar(ivar)%cache_data(2,CplVar(ivar)%shn,CplVar(ivar)%nlev) ) case ( 'll' ) #ifdef parallel_cplng call Get_DistGrid( dgrid(CplVar(ivar)%region), lli=local_lli) CplVar(ivar)%west_deg = local_lli%lon_deg(1) CplVar(ivar)%south_deg = local_lli%lat_deg(1) CplVar(ivar)%dlon_deg = local_lli%dlon_deg CplVar(ivar)%dlat_deg = local_lli%dlat_deg CplVar(ivar)%nlon = local_lli%nlon CplVar(ivar)%nlat = local_lli%nlat if ( CplVar(ivar)%cache ) & allocate( CplVar(ivar)%cache_data(CplVar(ivar)%nlon, CplVar(ivar)%nlat, CplVar(ivar)%nlev) ) #else CplVar(ivar)%west_deg = lli(CplVar(ivar)%region)%lon_deg(1) CplVar(ivar)%south_deg = lli(CplVar(ivar)%region)%lat_deg(1) CplVar(ivar)%dlon_deg = lli(CplVar(ivar)%region)%dlon_deg CplVar(ivar)%dlat_deg = lli(CplVar(ivar)%region)%dlat_deg CplVar(ivar)%nlon = lli(CplVar(ivar)%region)%nlon CplVar(ivar)%nlat = lli(CplVar(ivar)%region)%nlat if ( CplVar(ivar)%cache ) & allocate( CplVar(ivar)%cache_data(CplVar(ivar)%nlon,CplVar(ivar)%nlat,CplVar(ivar)%nlev) ) #endif case default write (gol,'("unsupported grid_type : ",a)') trim(CplVar(ivar)%grid_type); call goErr TRACEBACK; status=1; return end select end do ! levels end do ! var ! ! (7) FINALISE ! call oasis_enddef( status ) if (status/=OASIS_OK) then write (error_message,'("from OASIS_ENDDEF : ",i6)') status TRACEBACK call oasis_abort( comp_id, rname, error_message ) end if if (isRoot) then write (gol,'(" ")' ) ; call goPr write (gol,'("initialized oasis coupling:")' ) ; call goPr write (gol,'(" component : ",a)' ) trim(comp_name) ; call goPr write (gol,'(" real kind : ",i4)' ) wp ; call goPr write (gol,'(" ")' ) ; call goPr end if status = 0 END SUBROUTINE TM5_PRISM_INIT2 !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: TM5_PRISM_DONE ! ! !DESCRIPTION: cleanup (ie deallocate). !\\ !\\ ! !INTERFACE: ! SUBROUTINE TM5_Prism_Done( status ) ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 10 Sep 2013 - Ph. Le Sager - ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/TM5_Prism_Done' integer :: ireg integer :: ivar ! --- begin ----------------------------------- deallocate( part_id ) deallocate( var_shape ) deallocate( sp_var_shape ) ! clear descriptions: do ivar = 1, ncplvar deallocate( CplVar(ivar)%var_id ) deallocate( CplVar(ivar)%var_name ) if ( associated(CplVar(ivar)%cache_data) ) deallocate( CplVar(ivar)%cache_data ) end do status = 0 END SUBROUTINE TM5_PRISM_DONE !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: InqCplVar ! ! !DESCRIPTION: Inquire if there is a coupled variable with 'name'. !\\ !\\ ! !INTERFACE: ! SUBROUTINE InqCplVar( name, status, ivar, var_id, var_name, nlev ) ! ! !INPUT PARAMETERS: ! character(len=*), intent(in) :: name ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status integer, intent(out), optional :: ivar integer, intent(out), optional :: var_id(:) character(len=*), intent(out), optional :: var_name(:) integer, intent(out), optional :: nlev ! ! !REVISION HISTORY: ! 10 Sep 2013 - Ph. Le Sager - ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/InqCplVar' integer :: i, iv ! --- begin ----------------------------------- ! loop over defined variables: iv = -1 do i = 1, ncplvar ! check name: if ( trim(name) == trim(CplVar(i)%name) ) then iv = i exit end if end do ! not found ? if ( iv < 0 ) then write (gol,'("name of cplvar not found : ",a)') trim(name) ; call goErr write (gol,'(" available names : ")' ) ; call goErr do i = 1, ncplvar write (gol,'(" ",i4," ",a)') i, trim(CplVar(i)%name) ; call goErr end do end if ! fill requested arguments: if ( present(ivar ) ) ivar = iv if ( present(var_id ) ) var_id = CplVar(iv)%var_id if ( present(var_name) ) var_name = CplVar(iv)%var_name if ( present(nlev ) ) nlev = CplVar(iv)%nlev ! ok status = 0 END SUBROUTINE InqCplVar !EOC ! ************************************************************************** ! *** ! *** spectral field remapping routines ! *** ! ************************************************************************** !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: SHREMAP_INIT ! ! !DESCRIPTION: Init TshRemap object !\\ !\\ ! !INTERFACE: ! SUBROUTINE SHREMAP_INIT( shR, status ) ! ! !USES: ! use GO, only : NewDate ! ! !OUTPUT PARAMETERS: ! type(TshRemap), intent(out) :: shR integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/shremap_Init' ! --- begin --------------------------------------- ! no time stored yet: shR%t = NewDate(9999,9,9) ! safety: nullify( shR%remap ) ! nu truncation determined yet: shR%shT = 0 status = 0 END SUBROUTINE SHREMAP_INIT !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: SHREMAP_DONE ! ! !DESCRIPTION: deallocate var !\\ !\\ ! !INTERFACE: ! SUBROUTINE SHREMAP_DONE( shR, status ) ! ! !INPUT/OUTPUT PARAMETERS: ! type(TshRemap), intent(inout) :: shR ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/shremap_Done' ! --- begin --------------------------------------- if ( associated(shR%remap) ) deallocate( shR%remap ) status = 0 END SUBROUTINE SHREMAP_DONE !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: SHREMAP_SETUP ! ! !DESCRIPTION: !\\ !\\ ! !INTERFACE: ! SUBROUTINE SHREMAP_SETUP( shR, spinf, spinf_nan, status ) ! ! !INPUT/OUTPUT PARAMETERS: ! type(TshRemap), intent(inout) :: shR ! ! !INPUT PARAMETERS: ! real, intent(in) :: spinf(:,:,:) ! spectral info field real, intent(in) :: spinf_nan ! not-a-number in spinf ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/shremap_Setup' integer :: nlev integer :: sh_tripos(0:ifs_shT,0:ifs_shT) integer :: vri, vm, vn, vp, vk logical, allocatable :: sh_field(:,:,:) integer :: i, j, k real :: val integer :: nzero, nerr ! --- begin --------------------------------------- ! number of levels: nlev = size(spinf,2) ! triangle position: sh_tripos = 0 vp = 0 do vm = 0, ifs_shT do vn = vm, ifs_shT vp = vp + 1 sh_tripos(vm,vn) = vp end do end do ! storage for mapping indices: if ( .not. associated(shR%remap) ) then allocate( shR%remap(ifs_shn*2,nlev,3) ) end if shR%remap = -1 ! flags for target values; not filled remains negative: allocate( sh_field(2,ifs_shn,nlev) ) sh_field = .false. ! loop over levels: do k = 1, nlev ! no zero's detected yet ... nzero = 0 ! loop over spectral elements in layer: do i = 1, ifs_shn*2 !if (k==1) then !write (gol,'(" k, i, spinf : ",2i6,f16.4)') k, i, spinf(i,k,1); call goPr !endif ! not a number ? then this is an extra element due to the partitioning if ( spinf(i,k,1) == spinf_nan ) cycle ! extract m, n, and level: ! ! OLD : mmmnnnkk.0 real part ! mmmnnnkk.5 imag part ! !vri = 1 ! real part !if ( spinf(i,k,1)-floor(spinf(i,k,1)) == 0.5 ) vri = 2 ! imaginary part !vk = modulo( floor(spinf(i,k,1)), 100 ) ! level !vn = modulo( floor((spinf(i,k,1)-vk)/100.0), 1000 ) ! n !vm = floor(spinf(i,k,1)/100000.0) ! m ! ! NEW : mmmnnn.kk real part ! -mmmnnn.kk imag part ! ! Note that real and imag for m=0,n=0 are both 000000.00 for nlev=1; ! for nlev > 1, the values are both 000000.01 ! val = spinf(i,k,1) if ( val > 0.0 ) then ! positive value means real part: vri = 1 else if ( val < 0.0 ) then ! negative value means imag part: vri = 2 else ! zero values for both real and imag part of (0,0) nzero = nzero + 1 ! counter for number of zero values found !-- !! test number of zero values: !select case ( nzero ) ! case ( 1 ) ; vri = 1 ! real part of (0,0) ! case ( 2 ) ; vri = 2 ! imag part of (0,0) ! case default ! write (gol,'("found more than 2 zero values in spectral info, ")'); call goErr ! write (gol,'("while only expected for real and imag part of (0,0)")'); call goErr ! TRACEBACK; status=1; return ! cycle ! next value from received array !end select !-- ! assume that the extra zero's are the imaginary part for m=0, ! which is zero anyway and does not need to be received: if ( (nzero == 1) .and. (nlev == 1) ) then vri = 1 ! real part of (0,0) in spinf2d else cycle ! next value from received array end if end if val = abs(val) vk = nint( ( val - floor(val) )*100.0 ) ! level is fractional part vn = modulo( floor(val), 1000 ) ! last 3 numbers is n vm = floor( val/1000.0 ) ! first 3 numbers is m ! trap surface level ... if ( nlev == 1 ) vk = vk + 1 ! check ... if ( (vri < 1) .or. (vri > 2) .or. & (vm < 0) .or. (vm > ifs_shT) .or. (vn < vm) .or. (vn > ifs_shT) .or. & ((nlev==1) .and. (vk/=1)) .or. & ((nlev>1) .and. ((vk < 1) .or. (vk > nlev))) ) then write (gol,'("strange values extracted from spectral info value:")') ; call goErr write (gol,'(" spinf : ",f16.4)' ) spinf(i,k,1) ; call goErr write (gol,'(" ri : ",i4," (1=real,2=imag)")' ) vri ; call goErr write (gol,'(" m : ",i4," (0 .. ",i4,")")' ) vm, ifs_shT ; call goErr write (gol,'(" n : ",i4," (m .. ",i4,")")' ) vn, ifs_shT ; call goErr write (gol,'(" k : ",i4," (1 .. ",i4,")")' ) vk, nlev ; call goErr write (gol,'(" nzero : ",i4)' ) nzero ; call goErr TRACEBACK; status=1; return end if ! position in triangle: vp = sh_tripos(vm,vn) ! check ... if ( (vp < 1) .or. (vp > ifs_shn) ) then write (gol,'("strange triangle position:")' ) ; call goErr write (gol,'(" ifs T : ",i4)' ) ifs_shT ; call goErr write (gol,'(" m : ",i4)' ) vm ; call goErr write (gol,'(" n : ",i4)' ) vm ; call goErr write (gol,'(" p : ",i8)' ) vp ; call goErr TRACEBACK; status=1; return end if ! store: shR%remap(i,k,1) = vri shR%remap(i,k,2) = vp shR%remap(i,k,3) = vk ! maximum truncation: shR%shT = max( shR%shT, max( vm, vn ) ) ! flag ... sh_field(shR%remap(i,k,1),shR%remap(i,k,2),shR%remap(i,k,3)) = .true. end do ! received coeff i end do ! level k ! check on missing target values: if ( .not. all(sh_field) ) then ! error counter: nerr = 0 ! loop over levels: do k = 1, nlev ! init triangle position counter: vp = 0 ! loop over spectral triangle: do vm = 0, ifs_shT do vn = vm, ifs_shT ! increase triangle position counter: vp = vp + 1 ! negative values at unexpected places ? if ( ( (vm==0) .and. (.not. sh_field(1,vp,k) ) ) .or. & ( (vm> 0) .and. (.not. all(sh_field(:,vp,k))) ) ) then ! increase error counter: nerr = nerr + 1 ! intro message ? if ( nerr == 1 ) then write (gol,'("not all sh target values filled :")'); call goErr write (gol,'(" ifs T : ",i4)') ifs_shT; call goErr end if ! show error: write (gol,'(" (m, n) p, ; k ; real imag : (",2i5,") ",i8," ; ",i4," ; ",2l2)') vm, vn, vp, k, sh_field(:,vp,k); call goErr end if ! negatives ? end do ! n end do ! m end do ! level ! leave ? if (nerr>0) then TRACEBACK; status=1; return end if end if ! some negatives ? ! done deallocate( sh_field ) status = 0 END SUBROUTINE SHREMAP_SETUP !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: SHREMAP_REMAP ! ! !DESCRIPTION: !\\ !\\ ! !INTERFACE: ! SUBROUTINE SHREMAP_REMAP( shR, sh_recv, shi, sh_ri, status ) ! ! !USES: ! use grid, only : TshGridInfo ! ! !INPUT/OUTPUT PARAMETERS: ! type(TshRemap), intent(inout) :: shR ! ! !INPUT PARAMETERS: ! real, intent(in) :: sh_recv(:,:,:) type(TshGridInfo), intent(in) :: shi ! ! !OUTPUT PARAMETERS: ! real, intent(out) :: sh_ri(:,:,:) integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/shremap_Remap' integer :: nlev, i, k ! --- begin --------------------------------------- ! number of levels: nlev = size(sh_recv,2) ! check shape of input array ... if ( any( shape(sh_recv) /= (/ifs_shn*2,nlev,1/) ) ) then write (gol,'("strange input shape :")' ) ; call goErr write (gol,'(" sh_recv : ",3i6)' ) shape(sh_recv) ; call goErr write (gol,'(" expected : ",3i6)' ) ifs_shn*2,nlev,1 ; call goErr TRACEBACK; status=1; return end if ! check shape of output array ... if ( any( shape(sh_ri) /= (/2,ifs_shn,nlev/) ) ) then write (gol,'("strange input shape :")' ) ; call goErr write (gol,'(" sh : ",3i6)' ) shape(sh_ri) ; call goErr write (gol,'(" expected : ",3i6)' ) 2,ifs_shn,nlev ; call goErr TRACEBACK; status=1; return end if ! initial zero: sh_ri = 0.0 ! loop over all elements of received array: do k = 1, nlev do i = 1, ifs_shn*2 ! the triplet shR%remap(i,k,:) defines (/ 1=real/2=imag, traingle-position, level /) ! all negative ? ! o this is a dummy element due to partitioning ! o this is the imaginary part for m=0, which should remain zero if ( all( shR%remap(i,k,:) < 0 ) ) cycle ! any negative ? should not be possible... if ( any( shR%remap(i,k,:) < 0 ) ) then write (gol,'("found strange mapping:")' ) ; call goErr write (gol,'(" triangle point : ",i6)' ) i ; call goErr write (gol,'(" level : ",i6)' ) k ; call goErr write (gol,'(" mapping : ",3i6)' ) shR%remap(i,k,:) ; call goErr end if ! copy value from received array into spectral field: sh_ri(shR%remap(i,k,1),shR%remap(i,k,2),shR%remap(i,k,3)) = sh_recv(i,k,1) end do end do status = 0 END SUBROUTINE SHREMAP_REMAP !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: SetPrismTime_date ! ! !DESCRIPTION: returns current time/date into prism format (seconds from ! prism reference start). !\\ !\\ ! !INTERFACE: ! subroutine SetPrismTime_date( prism_t, t, status ) ! ! !USES: ! use GO, only : TDate, NewDate, iTotal, operator(-) ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: prism_t ! seconds from start ! ! !INPUT PARAMETERS: ! type(TDate), intent(in) :: t integer, intent(out) :: status ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/SetPrismTime_date' ! --- begin ---------------------------------------- ! seconds since start prism_t = iTotal( t - NewDate(time6=PRISM_start_date), 'sec' ) status = 0 end subroutine SetPrismTime_date !EOC END MODULE TM5_PRISM