!### macro's ##################################################### ! #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 ! #include "tm5.inc" ! !----------------------------------------------------------------------------- ! TM5 ! !----------------------------------------------------------------------------- !BOP ! ! !MODULE: GLOBAL_DATA ! ! !DESCRIPTION: hold global data (region_dat, wind_dat, conv_dat) and routines ! to allocate (DECLARE_FIELDS) and deallocate (FREE_FIELDS) them. ! Also holds (by inheritance) tracers masses (mass_dat, ! chem_dat) and allocate/deallocate them through their init/done methods. !\\ !\\ ! !INTERFACE: ! MODULE GLOBAL_DATA ! ! !USES: ! use GO, only : gol, goErr, TDate use dims, only : nregions, lm, lmax_conv use tracer_data, only : mass_dat, chem_dat, par_check_mass, tracer_print use global_types IMPLICIT NONE ! ! !PUBLIC DATA MEMBERS: ! character(len=256) :: rcfile ! name of rc file character(len=256) :: inputdir = './' ! path to input files character(len=256) :: outdir = './' ! path to output files logical :: fcmode ! is in forecast mode type(TDate), save :: tfcday0 ! day0 of forecast mode type(region_data), dimension(nregions), target :: region_dat type(wind_data) , dimension(nregions), target :: wind_dat #ifndef without_convection type(conv_data) , dimension(nregions), target :: conv_dat #endif ! ! !PRIVATE DATA MEMBERS: ! character(len=*), parameter, private :: mname = 'global_data' ! module name ! ! !REVISION HISTORY: ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ CONTAINS !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: DECLARE_FIELDS ! ! !DESCRIPTION: Allocate memory for data that are not handled with the meteo ! module (i.e. tracers, region_dat, wind_dat, conv_dat) !\\ !\\ ! !INTERFACE: ! SUBROUTINE DECLARE_FIELDS ! ! !USES: ! USE tracer_data, ONLY : tracer_init USE tm5_distgrid, ONLY : dgrid, Get_DistGrid ! ! !REVISION HISTORY: ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! !EOP !------------------------------------------------------------------------ !BOC INTEGER :: i0, i1, j0, j1, halo, region, lmr ! start do region=1,nregions lmr=lm(region) CALL Get_DistGrid( dgrid(region), I_STRT=i0, I_STOP=i1, J_STRT=j0, J_STOP=j1 ) halo=1 ALLOCATE ( region_dat(region)%zoomed(i0-halo:i1+halo,j0-halo:j1+halo)) ALLOCATE ( region_dat(region)%edge(i0-halo:i1+halo,j0-halo:j1+halo)) ALLOCATE ( region_dat(region)%dxyp(j0-halo:j1+halo)) region_dat(region)%halo = halo ! fluxes through boundaries halo=2 ALLOCATE ( wind_dat(region)%am( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) ) ALLOCATE ( wind_dat(region)%bm( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) ) ALLOCATE ( wind_dat(region)%cm( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) ) wind_dat(region)%halo = halo ! fill with zero for safety wind_dat(region)%am = 0.0 wind_dat(region)%bm = 0.0 wind_dat(region)%cm = 0.0 #ifndef without_convection #ifndef without_diffusion ALLOCATE ( conv_dat(region)%dkg(i0:i1,j0:j1,lmax_conv) ) ; conv_dat(region)%dkg = 0.0 ALLOCATE ( conv_dat(region)%blh(i0:i1,j0:j1) ) ; conv_dat(region)%blh = 0.0 #endif ALLOCATE ( conv_dat(region)%cloud_base( i0:i1, j0:j1 ) ) ALLOCATE ( conv_dat(region)%cloud_top ( i0:i1, j0:j1 ) ) ALLOCATE ( conv_dat(region)%cloud_lfs ( i0:i1, j0:j1 ) ) #endif end do ! allocate tracers CALL TRACER_INIT END SUBROUTINE DECLARE_FIELDS !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: FREE_FIELDS ! ! !DESCRIPTION: Deallocate data allocated through DECLARE_FIELDS !\\ !\\ ! !INTERFACE: ! SUBROUTINE FREE_FIELDS ! ! !USES: ! USE tracer_data, ONLY : tracer_done ! ! !REVISION HISTORY: ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! !EOP !------------------------------------------------------------------------ !BOC INTEGER :: region ! start do region=1,nregions DEALLOCATE ( region_dat(region)%zoomed ) DEALLOCATE ( region_dat(region)%edge ) DEALLOCATE ( region_dat(region)%dxyp ) DEALLOCATE ( wind_dat(region)%am ) DEALLOCATE ( wind_dat(region)%bm ) DEALLOCATE ( wind_dat(region)%cm ) #ifndef without_convection #ifndef without_diffusion DEALLOCATE ( conv_dat(region)%dkg ) DEALLOCATE ( conv_dat(region)%blh ) #endif DEALLOCATE ( conv_dat(region)%cloud_base ) DEALLOCATE ( conv_dat(region)%cloud_top ) DEALLOCATE ( conv_dat(region)%cloud_lfs ) #endif end do CALL TRACER_DONE END SUBROUTINE FREE_FIELDS !EOC END MODULE GLOBAL_DATA