MTObsFcInfo.F90 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225
  1. ! First include the set of model-wide compiler flags
  2. #include "tm5.inc"
  3. Module MTObsFcInfo
  4. !
  5. ! This module defines the structure to store the output data of the NRT retrieval
  6. ! Henk Eskes, KNMI, Feb 2015, Aug 2015
  7. !
  8. implicit none
  9. private
  10. public :: TObsFcInfo, ObsFcInfoAllocate, ObsFcInfoDeallocate
  11. ! ----------------------------------------------------
  12. ! structure containing observation, forecast and the observation operator
  13. ! ----------------------------------------------------
  14. !
  15. ! count : number of NO2 observations (ground pixels)
  16. ! copy of the number of observations provided as input
  17. ! flag : Quality flag
  18. !
  19. ! --> retrieval results
  20. ! note: vertical, slant columns and errors have units: 10^15 molecules cm^-2
  21. ! no2vcd,no2vcdsig : vertical column and error
  22. ! Note: this is not equal to the sum of trop and strat
  23. ! no2vcdsum : = no2vcdtrop + no2vcdstrat
  24. ! no2vcdsumsig : vertical column error (for no2vcdtrop + no2vcdstrat)
  25. ! no2vcdsigak : vertical column error, when kernels are used
  26. ! no2vcdtrop : tropospheric vertical column
  27. ! no2vcdtropsig : tropospheric vertical column error (including smoothing)
  28. ! no2vcdtropsigak : tropospheric vertical column error,
  29. ! when kernels are used (excluding smoothing)
  30. ! no2vcdstrat : stratospheric vertical column and error
  31. ! no2vcdstratsig : stratospheric vertical column and error
  32. ! no2slcstrat : (model) stratospheric slant column
  33. ! amf,amftrop,amfstrat : air-mass factors: total, tropospere, stratosphere (unit 1)
  34. ! avkernel : averaging kernel vector (unit 1)
  35. ! dimension ("lm","count")
  36. ! ghostcol : ghost column (model predicted part of column as obscured by clouds)
  37. !
  38. ! --> model derived values
  39. ! icell,jcell : pixel longitude, latitude grid index
  40. ! psurf : model surface pressure at pixel, unit "hPa"
  41. ! no2vcdfc,no2vcdfctrop : model forecast vertical column (10^15 mol/cm^2)
  42. ! no2slcfc,no2slcfctrop : model slant column (10^15 mol/cm^2)
  43. ! no2slcstrat : model stratospheric slant column (10^15 mol/cm^2)
  44. ! levtropopause : index of highest model level in the troposphere
  45. ! modelTerrainHeight : orography at TM5 model resolution, unit "m"
  46. !
  47. ! -- extra fields --
  48. ! cloudradfraction : (optional) computed when not available as input, [0-1]
  49. ! cloudradfraction_modified : .true. when CRF is computed in the observation operator module
  50. ! amfgeo : geometrical airmass factor (used internally)
  51. ! amfclear : clear-sky AMF (cloud fraction neglected)
  52. !
  53. ! The one-dimensional arrays all have dimension "count"
  54. ! The unit of the vertical and slant columns is "10^15 molecules / cm^2"
  55. ! The air mass factors have unit "1"
  56. !
  57. type TObsFcInfo
  58. integer :: count
  59. ! flagging
  60. integer,dimension(:),allocatable :: flag
  61. ! vcd, amf
  62. real,dimension(:),allocatable :: no2vcd
  63. real,dimension(:),allocatable :: no2vcdsig, no2vcdsigak
  64. real,dimension(:),allocatable :: no2vcdsum, no2vcdsumsig
  65. real,dimension(:),allocatable :: no2vcdtrop, no2vcdtropsig
  66. real,dimension(:),allocatable :: no2vcdtropsigak
  67. real,dimension(:),allocatable :: no2vcdstrat, no2vcdstratsig
  68. real,dimension(:),allocatable :: amf, amftrop, amfstrat
  69. real,dimension(:,:),allocatable :: avkernel
  70. real,dimension(:),allocatable :: ghostcol
  71. ! model derived quantities
  72. integer,dimension(:),allocatable :: icell, jcell
  73. real,dimension(:),allocatable :: psurf
  74. real,dimension(:),allocatable :: no2vcdfc, no2vcdfctrop
  75. real,dimension(:),allocatable :: no2slcfc, no2slcfctrop
  76. real,dimension(:),allocatable :: no2slcstrat
  77. integer,dimension(:),allocatable :: levtropopause
  78. real,dimension(:),allocatable :: modelTerrainHeight
  79. ! extra data
  80. real,dimension(:),allocatable :: amfgeo
  81. real,dimension(:),allocatable :: amfclear
  82. real,dimension(:),allocatable :: cloudRadFraction
  83. logical :: cloudRadFraction_computed
  84. end type TObsFcInfo
  85. contains
  86. subroutine ObsFcInfoAllocate ( obsFcInfo, npix, nlev )
  87. !
  88. ! Allocate arrays in structure
  89. !
  90. use pqf_module, only : PQF_E_GENERIC_EXCEPTION
  91. implicit none
  92. !
  93. type(TObsFcInfo),intent(inout) :: obsFcInfo
  94. integer, intent(in) :: npix, nlev
  95. !
  96. real,parameter :: nf_fill_double = 9.9692099683868690d+36
  97. integer, parameter :: nf_fill_int = -2147483647
  98. ! print *, 'Allocating ObsFcInfo, dimensions ',npix,nlev
  99. allocate ( obsFcInfo%icell(npix) )
  100. obsFcInfo%icell(:) = nf_fill_int
  101. allocate ( obsFcInfo%jcell(npix) )
  102. obsFcInfo%jcell(:) = nf_fill_int
  103. allocate ( obsFcInfo%no2vcd(npix) )
  104. obsFcInfo%no2vcd(:) = nf_fill_double
  105. allocate ( obsFcInfo%no2vcdsum(npix) )
  106. obsFcInfo%no2vcdsum(:) = nf_fill_double
  107. allocate ( obsFcInfo%no2vcdsumsig(npix) )
  108. obsFcInfo%no2vcdsumsig(:) = nf_fill_double
  109. allocate ( obsFcInfo%no2vcdsig(npix) )
  110. obsFcInfo%no2vcdsig(:) = nf_fill_double
  111. allocate ( obsFcInfo%no2vcdsigak(npix) )
  112. obsFcInfo%no2vcdsigak(:) = nf_fill_double
  113. allocate ( obsFcInfo%no2vcdtrop(npix) )
  114. obsFcInfo%no2vcdtrop(:) = nf_fill_double
  115. allocate ( obsFcInfo%no2vcdtropsig(npix) )
  116. obsFcInfo%no2vcdtropsig(:) = nf_fill_double
  117. allocate ( obsFcInfo%no2vcdtropsigak(npix) )
  118. obsFcInfo%no2vcdtropsigak(:) = nf_fill_double
  119. allocate ( obsFcInfo%no2vcdstrat(npix) )
  120. obsFcInfo%no2vcdstrat(:) = nf_fill_double
  121. allocate ( obsFcInfo%no2vcdstratsig(npix) )
  122. obsFcInfo%no2vcdstratsig(:) = nf_fill_double
  123. allocate ( obsFcInfo%no2slcstrat(npix) )
  124. obsFcInfo%no2vcdstrat(:) = nf_fill_double
  125. allocate ( obsFcInfo%amf(npix) )
  126. obsFcInfo%amf(:) = nf_fill_double
  127. allocate ( obsFcInfo%amftrop(npix) )
  128. obsFcInfo%amftrop(:) = nf_fill_double
  129. allocate ( obsFcInfo%amfstrat(npix) )
  130. obsFcInfo%amfstrat(:) = nf_fill_double
  131. allocate ( obsFcInfo%psurf(npix) )
  132. obsFcInfo%psurf(:) = nf_fill_double
  133. allocate ( obsFcInfo%no2slcfc(npix) )
  134. obsFcInfo%no2slcfc(:) = nf_fill_double
  135. allocate ( obsFcInfo%no2slcfctrop(npix) )
  136. obsFcInfo%no2slcfctrop(:) = nf_fill_double
  137. allocate ( obsFcInfo%no2vcdfc(npix) )
  138. obsFcInfo%no2vcdfc(:) = nf_fill_double
  139. allocate ( obsFcInfo%no2vcdfctrop(npix) )
  140. obsFcInfo%no2vcdfctrop(:) = nf_fill_double
  141. allocate ( obsFcInfo%flag(npix) )
  142. obsFcInfo%flag(:) = PQF_E_GENERIC_EXCEPTION
  143. allocate ( obsFcInfo%levtropopause(npix) )
  144. obsFcInfo%levtropopause(:) = nf_fill_int
  145. allocate ( obsFcInfo%avkernel(nlev,npix) )
  146. obsFcInfo%avkernel(:,:) = nf_fill_double
  147. allocate ( obsFcInfo%ghostcol(npix) )
  148. obsFcInfo%ghostcol(:) = nf_fill_double
  149. allocate ( obsFcInfo%modelterrainheight(npix) )
  150. obsFcInfo%modelterrainheight(:) = nf_fill_double
  151. ! extra
  152. allocate ( obsFcInfo%cloudradfraction(npix) )
  153. obsFcInfo%cloudradfraction(:) = nf_fill_double
  154. allocate ( obsFcInfo%amfgeo(npix) )
  155. obsFcInfo%amfgeo(:) = nf_fill_double
  156. allocate ( obsFcInfo%amfclear(npix) )
  157. obsFcInfo%amfclear(:) = nf_fill_double
  158. !
  159. obsFcInfo%count = npix
  160. !
  161. end subroutine ObsFcInfoAllocate
  162. subroutine ObsFcInfoDeallocate ( obsFcInfo )
  163. !
  164. ! De-allocate arrays in structure
  165. !
  166. implicit none
  167. !
  168. type(TObsFcInfo),intent(inout) :: obsFcInfo
  169. !
  170. if ( allocated( obsFcInfo%icell ) ) deallocate ( obsFcInfo%icell )
  171. if ( allocated( obsFcInfo%jcell ) ) deallocate ( obsFcInfo%jcell )
  172. if ( allocated( obsFcInfo%no2vcd ) ) deallocate ( obsFcInfo%no2vcd )
  173. if ( allocated( obsFcInfo%no2vcdsum ) ) deallocate ( obsFcInfo%no2vcdsum )
  174. if ( allocated( obsFcInfo%no2vcdsumsig ) ) deallocate ( obsFcInfo%no2vcdsumsig )
  175. if ( allocated( obsFcInfo%no2vcdsig ) ) deallocate ( obsFcInfo%no2vcdsig )
  176. if ( allocated( obsFcInfo%no2vcdsigak ) ) deallocate ( obsFcInfo%no2vcdsigak )
  177. if ( allocated( obsFcInfo%no2vcdtrop ) ) deallocate ( obsFcInfo%no2vcdtrop )
  178. if ( allocated( obsFcInfo%no2vcdtropsig ) ) deallocate ( obsFcInfo%no2vcdtropsig )
  179. if ( allocated( obsFcInfo%no2vcdtropsigak ) ) deallocate ( obsFcInfo%no2vcdtropsigak )
  180. if ( allocated( obsFcInfo%no2vcdstrat ) ) deallocate ( obsFcInfo%no2vcdstrat )
  181. if ( allocated( obsFcInfo%no2vcdstratsig ) ) deallocate ( obsFcInfo%no2vcdstratsig )
  182. if ( allocated( obsFcInfo%no2slcstrat ) ) deallocate ( obsFcInfo%no2slcstrat )
  183. if ( allocated( obsFcInfo%amf ) ) deallocate ( obsFcInfo%amf )
  184. if ( allocated( obsFcInfo%amftrop ) ) deallocate ( obsFcInfo%amftrop )
  185. if ( allocated( obsFcInfo%amfstrat ) ) deallocate ( obsFcInfo%amfstrat )
  186. if ( allocated( obsFcInfo%psurf ) ) deallocate ( obsFcInfo%psurf )
  187. if ( allocated( obsFcInfo%no2slcfc ) ) deallocate ( obsFcInfo%no2slcfc )
  188. if ( allocated( obsFcInfo%no2slcfctrop ) ) deallocate ( obsFcInfo%no2slcfctrop )
  189. if ( allocated( obsFcInfo%no2vcdfc ) ) deallocate ( obsFcInfo%no2vcdfc )
  190. if ( allocated( obsFcInfo%no2vcdfctrop ) ) deallocate ( obsFcInfo%no2vcdfctrop )
  191. if ( allocated( obsFcInfo%flag ) ) deallocate ( obsFcInfo%flag )
  192. if ( allocated( obsFcInfo%levtropopause ) ) deallocate ( obsFcInfo%levtropopause )
  193. if ( allocated( obsFcInfo%avkernel ) ) deallocate ( obsFcInfo%avkernel )
  194. if ( allocated( obsFcInfo%ghostcol ) ) deallocate ( obsFcInfo%ghostcol )
  195. if ( allocated( obsFcInfo%modelterrainheight ) ) deallocate ( obsFcInfo%modelterrainheight )
  196. ! extra
  197. if ( allocated( obsFcInfo%cloudradfraction ) ) deallocate ( obsFcInfo%cloudradfraction )
  198. if ( allocated( obsFcInfo%amfgeo ) ) deallocate ( obsFcInfo%amfgeo )
  199. if ( allocated( obsFcInfo%amfclear ) ) deallocate ( obsFcInfo%amfclear )
  200. !
  201. obsFcInfo%count = -1
  202. !
  203. end subroutine ObsFcInfoDeallocate
  204. end Module MTObsFcInfo