!################################################################# ! ! this module declares and deallocates the actual data ! dimensions are defined in module dims ! MK december 2002 ! !### 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" ! !################################################################# module global_data use GO, only : gol, goPr, goErr use GO , only : TDate use dims , only : region_data, wind_data use dims , only : nregions, im, jm, lm, nlon360, nlat180 use dims , only : okdebug, lmax_conv use chem_param , only : ntracet use global_types, only : emis_data use global_types, only : d2_data, d23_data, d3_data #ifndef without_convection use global_types, only : conv_data #endif use tracer_data , only : mass_data, mass_dat implicit none ! --- in/out ----------------------------------- private ! variables ! NOTE: ifort compiler requires that both type and variable are public .. #ifdef with_pycasso public :: rcfile #else public :: rcfile, rcfile_runtime #endif public :: inputdir public :: outdir public :: region_data, region_dat public :: wind_data, wind_dat public :: mass_data, mass_dat #ifndef without_convection public :: conv_data, conv_dat #endif public :: emis_data public :: d2_data, d23_data, d3_data public :: fcmode, tfcday0 ! subroutines public :: declare_fields, free_fields ! --- const ---------------------------------- ! module name character(len=*), parameter :: mname = 'global_data' ! --- var --------------------------------------------- #ifdef with_pycasso ! name of rc file: character(len=256) :: rcfile #else ! defintion of parameter 'rcfile' with rcfile name: include "rcfilename.h" #endif ! path to input files: character(len=256) :: inputdir = './' ! path to output files: character(len=256) :: outdir = './' 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 ! forecast mode ? logical :: fcmode type(TDate), save :: tfcday0 contains ! ========================================================== ! ! subroutine to allocate the memory for all the regions ! should be called at start of program ! subroutine declare_fields use tracer_data, only : assign_massdat use partools, only : lmloc, allocate_mass ! local integer :: region integer :: imr,jmr,lmr ! start do region=1,nregions imr = im(region) jmr = jm(region) lmr = lm(region) allocate ( region_dat(region)%zoomed(imr,jmr)) allocate ( region_dat(region)%edge(imr,jmr)) allocate ( region_dat(region)%dxyp(jmr)) ! fluxes through boundaries: allocate ( wind_dat(region)%am_t(0:imr+1,0:jmr+1,0:lmr+1) ) allocate ( wind_dat(region)%bm_t(0:imr+1,0:jmr+1,0:lmr+1) ) allocate ( wind_dat(region)%cm_t(0:imr+1,0:jmr+1,0:lmr+1) ) ! fill with zero for safety ... wind_dat(region)%am_t = 0.0 wind_dat(region)%bm_t = 0.0 wind_dat(region)%cm_t = 0.0 #ifndef without_convection #ifndef without_diffusion allocate ( conv_dat(region)%dkg(imr,jmr,lmax_conv) ) ; conv_dat(region)%dkg = 0.0 allocate ( conv_dat(region)%blh(imr,jmr) ) ; conv_dat(region)%blh = 0.0 #endif allocate ( conv_dat(region)%cloud_base(imr,jmr)) allocate ( conv_dat(region)%cloud_top (imr,jmr)) allocate ( conv_dat(region)%cloud_lfs (imr,jmr)) #endif ! WP! initialize mass_dat on tracer domain call assign_massdat( region, 1 ) ! WP! initialize mass_dat on level domain if ( .not. allocate_mass ) call assign_massdat( region, 0 ) if ( okdebug ) print *, 'declare_fields: Fields for grid, mass, ', & 'and winds allocated, region:', region end do end subroutine declare_fields subroutine free_fields ! ! subroutine to free the memory for all the regions ! should be called at end of program ! use tracer_data, only : free_massdat use partools , only : myid, allocate_mass implicit none ! local integer :: region ! start do region=1,nregions #ifdef MPI if ( okdebug ) print *, 'free_fields: deallocate... ', & myid, region,im(region),jm(region) #else if ( okdebug ) print *, 'free_fields: deallocate...', & region,im(region),jm(region) #endif deallocate ( region_dat(region)%zoomed ) deallocate ( region_dat(region)%edge ) deallocate ( region_dat(region)%dxyp ) deallocate ( wind_dat(region)%am_t ) deallocate ( wind_dat(region)%bm_t ) deallocate ( wind_dat(region)%cm_t ) #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 ! WP! free mass_dat on tracer domain call free_massdat( region, 0 ) ! WP! free mass_dat on level domain if ( .not. allocate_mass ) call free_massdat( region, 1 ) if ( okdebug ) print *, 'free_fields: Fields for grid, mass, ', & 'and winds deallocated, region:', myid, region end do end subroutine free_fields end module global_data