global_data.F90 5.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191
  1. !### macro's #####################################################
  2. !
  3. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  4. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  5. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  6. !
  7. #include "tm5.inc"
  8. !
  9. !-----------------------------------------------------------------------------
  10. ! TM5 !
  11. !-----------------------------------------------------------------------------
  12. !BOP
  13. !
  14. ! !MODULE: GLOBAL_DATA
  15. !
  16. ! !DESCRIPTION: hold global data (region_dat, wind_dat, conv_dat) and routines
  17. ! to allocate (DECLARE_FIELDS) and deallocate (FREE_FIELDS) them.
  18. ! Also holds (by inheritance) tracers masses (mass_dat,
  19. ! chem_dat) and allocate/deallocate them through their init/done methods.
  20. !\\
  21. !\\
  22. ! !INTERFACE:
  23. !
  24. MODULE GLOBAL_DATA
  25. !
  26. ! !USES:
  27. !
  28. use GO, only : gol, goErr, TDate
  29. use dims, only : nregions, lm, lmax_conv
  30. use tracer_data, only : mass_dat, chem_dat, par_check_mass, tracer_print
  31. use global_types
  32. IMPLICIT NONE
  33. !
  34. ! !PUBLIC DATA MEMBERS:
  35. !
  36. character(len=256) :: rcfile ! name of rc file
  37. character(len=256) :: inputdir = './' ! path to input files
  38. character(len=256) :: outdir = './' ! path to output files
  39. logical :: fcmode ! is in forecast mode
  40. type(TDate), save :: tfcday0 ! day0 of forecast mode
  41. type(region_data), dimension(nregions), target :: region_dat
  42. type(wind_data) , dimension(nregions), target :: wind_dat
  43. #ifndef without_convection
  44. type(conv_data) , dimension(nregions), target :: conv_dat
  45. #endif
  46. !
  47. ! !PRIVATE DATA MEMBERS:
  48. !
  49. character(len=*), parameter, private :: mname = 'global_data' ! module name
  50. !
  51. ! !REVISION HISTORY:
  52. ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  53. !
  54. ! !REMARKS:
  55. !
  56. !EOP
  57. !------------------------------------------------------------------------
  58. CONTAINS
  59. !--------------------------------------------------------------------------
  60. ! TM5 !
  61. !--------------------------------------------------------------------------
  62. !BOP
  63. !
  64. ! !IROUTINE: DECLARE_FIELDS
  65. !
  66. ! !DESCRIPTION: Allocate memory for data that are not handled with the meteo
  67. ! module (i.e. tracers, region_dat, wind_dat, conv_dat)
  68. !\\
  69. !\\
  70. ! !INTERFACE:
  71. !
  72. SUBROUTINE DECLARE_FIELDS
  73. !
  74. ! !USES:
  75. !
  76. USE tracer_data, ONLY : tracer_init
  77. USE tm5_distgrid, ONLY : dgrid, Get_DistGrid
  78. !
  79. ! !REVISION HISTORY:
  80. ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  81. !
  82. !EOP
  83. !------------------------------------------------------------------------
  84. !BOC
  85. INTEGER :: i0, i1, j0, j1, halo, region, lmr
  86. ! start
  87. do region=1,nregions
  88. lmr=lm(region)
  89. CALL Get_DistGrid( dgrid(region), I_STRT=i0, I_STOP=i1, J_STRT=j0, J_STOP=j1 )
  90. halo=1
  91. ALLOCATE ( region_dat(region)%zoomed(i0-halo:i1+halo,j0-halo:j1+halo))
  92. ALLOCATE ( region_dat(region)%edge(i0-halo:i1+halo,j0-halo:j1+halo))
  93. ALLOCATE ( region_dat(region)%dxyp(j0-halo:j1+halo))
  94. region_dat(region)%halo = halo
  95. ! fluxes through boundaries
  96. halo=2
  97. ALLOCATE ( wind_dat(region)%am( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) )
  98. ALLOCATE ( wind_dat(region)%bm( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) )
  99. ALLOCATE ( wind_dat(region)%cm( i0-halo:i1+halo, j0-halo:j1+halo, 0:lmr+1) )
  100. wind_dat(region)%halo = halo
  101. ! fill with zero for safety
  102. wind_dat(region)%am = 0.0
  103. wind_dat(region)%bm = 0.0
  104. wind_dat(region)%cm = 0.0
  105. #ifndef without_convection
  106. #ifndef without_diffusion
  107. ALLOCATE ( conv_dat(region)%dkg(i0:i1,j0:j1,lmax_conv) ) ; conv_dat(region)%dkg = 0.0
  108. ALLOCATE ( conv_dat(region)%blh(i0:i1,j0:j1) ) ; conv_dat(region)%blh = 0.0
  109. #endif
  110. ALLOCATE ( conv_dat(region)%cloud_base( i0:i1, j0:j1 ) )
  111. ALLOCATE ( conv_dat(region)%cloud_top ( i0:i1, j0:j1 ) )
  112. ALLOCATE ( conv_dat(region)%cloud_lfs ( i0:i1, j0:j1 ) )
  113. #endif
  114. end do
  115. ! allocate tracers
  116. CALL TRACER_INIT
  117. END SUBROUTINE DECLARE_FIELDS
  118. !EOC
  119. !--------------------------------------------------------------------------
  120. ! TM5 !
  121. !--------------------------------------------------------------------------
  122. !BOP
  123. !
  124. ! !IROUTINE: FREE_FIELDS
  125. !
  126. ! !DESCRIPTION: Deallocate data allocated through DECLARE_FIELDS
  127. !\\
  128. !\\
  129. ! !INTERFACE:
  130. !
  131. SUBROUTINE FREE_FIELDS
  132. !
  133. ! !USES:
  134. !
  135. USE tracer_data, ONLY : tracer_done
  136. !
  137. ! !REVISION HISTORY:
  138. ! 9 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  139. !
  140. !EOP
  141. !------------------------------------------------------------------------
  142. !BOC
  143. INTEGER :: region
  144. ! start
  145. do region=1,nregions
  146. DEALLOCATE ( region_dat(region)%zoomed )
  147. DEALLOCATE ( region_dat(region)%edge )
  148. DEALLOCATE ( region_dat(region)%dxyp )
  149. DEALLOCATE ( wind_dat(region)%am )
  150. DEALLOCATE ( wind_dat(region)%bm )
  151. DEALLOCATE ( wind_dat(region)%cm )
  152. #ifndef without_convection
  153. #ifndef without_diffusion
  154. DEALLOCATE ( conv_dat(region)%dkg )
  155. DEALLOCATE ( conv_dat(region)%blh )
  156. #endif
  157. DEALLOCATE ( conv_dat(region)%cloud_base )
  158. DEALLOCATE ( conv_dat(region)%cloud_top )
  159. DEALLOCATE ( conv_dat(region)%cloud_lfs )
  160. #endif
  161. end do
  162. CALL TRACER_DONE
  163. END SUBROUTINE FREE_FIELDS
  164. !EOC
  165. END MODULE GLOBAL_DATA