[2404] | 1 | program simu_MCS |
---|
| 2 | ! program for a MCS observer simulator of GCM data |
---|
| 3 | ! author : Antoine Bierjon, 2019-2020 |
---|
| 4 | ! contact : antoine.bierjon@lmd.jussieu.fr |
---|
| 5 | ! |
---|
| 6 | !=================================================================================================== |
---|
| 7 | ! PREFACE |
---|
| 8 | !=================================================================================================== |
---|
| 9 | ! This program loads on one hand, a GCM file, which can be : |
---|
| 10 | ! - a GCM diagfi.nc or concat.nc that covers the observation period, |
---|
| 11 | ! - a GCM stats file (stats.nc) which gets extended to cover the observation period |
---|
| 12 | ! |
---|
| 13 | ! On the other hand, it loads a MRO/MCS data file (binned by Luca Montabone). |
---|
| 14 | ! This MCS file serves then as a reference for the interpolation and the binning of the |
---|
| 15 | ! GCM variables contained in the GCM file and given as inputs in the simu_MCS.def |
---|
| 16 | ! |
---|
| 17 | ! Since the MCS data is binned by intervals of Ls, the GCM simulation's data is interpolated at |
---|
| 18 | ! the MCS spatial coordinates, at every GCM sol value contained in each interval of Ls, before |
---|
| 19 | ! being averaged to create output bins for each interval. The binning in sols also takes into |
---|
| 20 | ! account the variability of Local Time in the bin and tries to represent it (with a number of LT |
---|
| 21 | ! values in each sol bin equal to the number of values in each associated MCS bin : numbintemp/dust/wice, |
---|
| 22 | ! and centered around the MCS bin LT average : timeave). |
---|
| 23 | ! |
---|
| 24 | ! There are also specific GCM variables that the program looks for in order to make them |
---|
| 25 | ! comparable with their equivalent in MCS files. These variables are : |
---|
| 26 | ! GCM MCS |
---|
| 27 | ! # temp ---> dtemp,ntemp |
---|
| 28 | ! # dso(/dsodust/qdust)+rho ---> ddust,ndust (dust opacity/km) |
---|
| 29 | ! # h2o_ice+rho ---> dwice,nwice (water ice opacity/km) |
---|
| 30 | ! |
---|
[2558] | 31 | ! For dust especially, if a dust opacity variable has been created in the output file, the program also |
---|
| 32 | ! computes the ratio of the Column-integrated Dust Optical Depths from the binned GCM file and the |
---|
| 33 | ! MCS file, that can then be used to normalize the dust opacity profiles when doing comparisons. |
---|
| 34 | ! |
---|
[2404] | 35 | ! Eventually, the program outputs a netcdf file, filled with the GCM binned data for dayside |
---|
| 36 | ! and nightside and edited with the same format than the MCS file (and with the same |
---|
| 37 | ! dimensions' declaration order). |
---|
| 38 | ! |
---|
| 39 | ! Minimal requirements : |
---|
| 40 | ! - the MCS file must contain the variables : |
---|
| 41 | ! dtimeave,dtimemax,dtimemin, |
---|
| 42 | ! dtemp,dnumbintemp, |
---|
| 43 | ! ddust,dnumbindust, |
---|
| 44 | ! dwice,dnumbinwice, |
---|
| 45 | ! ntimeave,ntimemax,ntimemin, |
---|
| 46 | ! ntemp,nnumbintemp, |
---|
| 47 | ! ndust,nnumbindust, |
---|
| 48 | ! nwice,nnumbinwice |
---|
| 49 | ! - the GCM file must : |
---|
| 50 | ! # have the same altitude type as the MCS file ; |
---|
| 51 | ! # cover completely the observation period |
---|
| 52 | ! |
---|
| 53 | ! See also NOTA BENE in section 2.2 |
---|
| 54 | ! |
---|
| 55 | ! |
---|
| 56 | ! |
---|
| 57 | ! Algorithm : |
---|
| 58 | ! 0. Variable declarations |
---|
| 59 | ! 1. OPENING OF THE FILES AND INITIALIZATION OF THE DIMENSIONS |
---|
| 60 | ! 1.1 MCS data file : obsfile |
---|
| 61 | ! 1.1.1 Open the Observer data file |
---|
| 62 | ! 1.1.2 Get dimensions lon,lat,alt,time from the observation file |
---|
| 63 | ! 1.2. GCM simulation file : gcmfile |
---|
| 64 | ! 1.2.1 Open the GCM simulation file |
---|
| 65 | ! 1.2.2 Get dimensions lon,lat,alt,time from the GCM file |
---|
| 66 | ! 1.3 Create the output file & initialize the coordinates |
---|
| 67 | ! 2. VARIABLES MANAGEMENT |
---|
| 68 | ! 2.1 List of the GCM variables to be interpolated |
---|
| 69 | ! 2.1.1 Read the GCM variables |
---|
| 70 | ! 2.1.2 Handle dust and wice opacities (first set-ups) |
---|
| 71 | ! 2.2 Definition of LT variables from obsfile to outfile |
---|
| 72 | ![day/night loop begins... |
---|
| 73 | ! 2.2.1 Average of Local Time in the OBS bins |
---|
| 74 | ! 2.2.2 Maximum of Local Time in the OBS bins |
---|
| 75 | ! 2.2.3 Minimum of Local Time in the OBS bins |
---|
| 76 | ! 2.3 Definition of numbin variables from obsfile to outfile |
---|
| 77 | ! 2.3.1 Number of values in the OBS "temp" bins |
---|
| 78 | ! 2.3.2 Number of values in the OBS "dust" bins |
---|
| 79 | ! 2.3.3 Number of values in the OBS "wice" bins |
---|
| 80 | ! 2.4 Opening of the GCM variables |
---|
| 81 | ! [var loop begins... |
---|
| 82 | ! 2.4.1 Generic reading of the variable |
---|
| 83 | ! 2.4.2 Handle dust and wice opacities (second part) |
---|
| 84 | ! 2.5 Opening of the associated MCS variables |
---|
| 85 | ! 2.5.1 MCS reference variable (for the missing values) |
---|
| 86 | ! 2.5.2 Number of values in the OBS bin (for the sol binning) |
---|
| 87 | ! 2.6 Definition of GCM variables in outfile |
---|
| 88 | ! 3. EXTRACTION OF THE VARIABLE |
---|
| 89 | ! [coordinates loop begins... |
---|
| 90 | ! 3.1 Do some checks and preparations before the extraction |
---|
| 91 | ! 3.2 Compute GCM sol date corresponding to Observer Ls (via m_(min/max)sol) |
---|
| 92 | ! and LT (via OBSLT(min/max)) |
---|
| 93 | ! 3.3 Do the interpolation and binning for the given location |
---|
| 94 | ! ..coordinates loop ends] |
---|
| 95 | ! 3.4 Write the data in the netcdf output file |
---|
| 96 | ! ..var loop ends] |
---|
[2558] | 97 | ! 4. CDOD RATIO (tau_GCM/tau_MCS) |
---|
| 98 | ! 4.1 Preliminary checks |
---|
| 99 | ! 4.2 tau_ratio computation |
---|
| 100 | ! 4.3 Write tau_ratio in the netcdf output file |
---|
| 101 | ! 5. END OF THE DAY/NIGHT LOOP |
---|
[2404] | 102 | !..day/night loop ends] |
---|
[2558] | 103 | ! 6. CLOSE THE FILES |
---|
[2404] | 104 | ! |
---|
| 105 | ! Subroutines |
---|
| 106 | ! extraction |
---|
| 107 | ! inidim |
---|
| 108 | ! ls2sol |
---|
| 109 | ! gen_sol_list |
---|
| 110 | ! status_check |
---|
| 111 | ! LTmod |
---|
| 112 | !=================================================================================================== |
---|
| 113 | |
---|
| 114 | |
---|
| 115 | use netcdf |
---|
| 116 | |
---|
| 117 | |
---|
| 118 | !=================================================================================== |
---|
| 119 | ! 0. Variable declarations |
---|
| 120 | !=================================================================================== |
---|
| 121 | |
---|
| 122 | implicit none ! for no implicitly typed variables |
---|
| 123 | |
---|
| 124 | !------------------------ |
---|
| 125 | ! Files: |
---|
| 126 | character(len=256) :: gcmfile ! GCM simulation file |
---|
| 127 | logical :: is_stats = .false. ! to check if the GCM file is a stats.nc or a diagfi.nc file |
---|
| 128 | character(len=256) :: obsfile ! observation file = MCS/Observer data file |
---|
| 129 | character(len=256) :: outfile ! output file |
---|
| 130 | |
---|
| 131 | !------------------------ |
---|
| 132 | ! NetCDF stuff |
---|
| 133 | integer :: status ! NetCDF routines return code |
---|
| 134 | character (len=256) :: error_text ! to store the error text |
---|
| 135 | integer :: gcmfid ! NetCDF gcm file ID |
---|
| 136 | integer :: obsfid ! NetCDF observation file ID |
---|
| 137 | integer :: outfid ! NetCDF output file ID |
---|
| 138 | integer :: GCMvarid, OBSvarid, LT_id, numbin_id, outvarid ! to store the ID of a variable |
---|
| 139 | integer :: lat_dimid_obs,lon_dimid_obs,alt_dimid_obs,time_dimid_obs ! dimensions' ID in OBS and output files |
---|
| 140 | integer :: lat_dimid_gcm,lon_dimid_gcm,alt_dimid_gcm,time_dimid_gcm ! dimensions' ID in GCM file |
---|
| 141 | integer :: GCMvarshape(4), OBSvarshape(4), LTshape(3),numbinshape(4) ! to store a variable's coordinates order |
---|
| 142 | |
---|
| 143 | !------------------------ |
---|
| 144 | ! Dimensions |
---|
| 145 | real,dimension(:),allocatable :: GCMlon, OBSlon ! longitude in the GCM & the Observer files |
---|
| 146 | integer GCMlonlen, OBSlonlen ! # of grid points along GCMlon & OBSlon |
---|
| 147 | real,dimension(:),allocatable :: GCMlat, OBSlat ! latitude in the GCM & the Observer files |
---|
| 148 | integer GCMlatlen, OBSlatlen ! # of grid points along GCMlat & OBSlat |
---|
| 149 | real,dimension(:),allocatable :: GCMalt, OBSalt ! altitude/pressure in the GCM & the Observer files |
---|
| 150 | integer GCMaltlen, OBSaltlen ! # of grid point along GCMalt & OBSalt |
---|
| 151 | character(len=1) :: GCMalttype, OBSalttype ! altitude coord. type:'z' (altitude, m) 'p' (pressure, Pa) |
---|
| 152 | real,dimension(:),allocatable :: GCMtime, OBSLs ! time in the GCM diagfi (sols) & the Observer files (Ls) |
---|
| 153 | real,dimension(:),allocatable :: GCMstatstime ! time in the GCM stats file (LT at lon 0°) |
---|
| 154 | integer :: GCMtimelen, GCMstatstimelen, OBSLslen ! # of points along GCMtime, GCMstatstime, OBSLs |
---|
| 155 | real :: starttimeoffset=0. ! offset (in sols) wrt Ls=0 of sol 0 in GCM file |
---|
| 156 | |
---|
| 157 | !------------------------ |
---|
| 158 | ! Variables |
---|
| 159 | character(len=10) :: dayornight ! are we in the "dayside" or "nightside" loop? |
---|
| 160 | character(len=64),dimension(:),allocatable :: gcm_vars ! list of GCM variables to interpolate |
---|
| 161 | character(len=10),dimension(15) :: notprocessed ! names of the (15) variables that won't be processed |
---|
| 162 | integer :: nbvarfile,Nnotprocessed ! nbs of variables to deal with the non-processed ones |
---|
| 163 | integer :: nbvar ! nb of variables that will be processed |
---|
| 164 | logical :: var_ok ! is this variable to be processed? |
---|
[2525] | 165 | logical :: dustok1,dustok2,dustok3,wiceok ! is it possible to compute opacities and how? |
---|
[2404] | 166 | character(len=64) :: GCMvarname,OBSvarname,outvarname ! name of the variables |
---|
| 167 | integer :: nbdim ! nb of dimensions of a variable |
---|
| 168 | real,dimension(:,:,:,:),allocatable :: GCM_var,OBS_var ! the 4D variable extracted from GCM & OBS files |
---|
[2558] | 169 | real,dimension(:,:,:,:),allocatable :: GCM_rho ! atmospheric density for opacities' computation |
---|
[2404] | 170 | character(len=64) :: long_name,units,comment ! netcdf attributes of the variable |
---|
| 171 | |
---|
| 172 | character(len=64) :: OBSLTave_name,OBSLTmax_name,OBSLTmin_name |
---|
| 173 | ! names of the average, max and min of LT |
---|
| 174 | real,dimension(:,:,:),allocatable :: OBSLT,OBSLTmax,OBSLTmin |
---|
| 175 | ! 3D variables extracted from obsfile (ave, max and min of LT in a bin) |
---|
| 176 | |
---|
| 177 | character (len=64) :: numbin_name ! name of the nb of values in an OBS bin |
---|
| 178 | real,dimension(:,:,:,:),allocatable :: numbin ! nb of values in an OBS temp/dust/wice bin |
---|
| 179 | |
---|
| 180 | real :: GCMmiss_val, OBSmiss_val, LTmiss_val ! value to denote non-existant data |
---|
| 181 | real :: extr_value ! result of the extraction subroutine |
---|
| 182 | real, dimension(:,:,:,:), allocatable :: outvar ! outvar(,,,): 4D array to store the output variable's data |
---|
| 183 | |
---|
| 184 | !------------------------ |
---|
| 185 | ! Time binning management |
---|
| 186 | real :: OBSdeltaLs ! difference of Ls between each observation bin |
---|
| 187 | real :: sol, maxsol, minsol ! sol date corresponding to Observed Ls and LT |
---|
| 188 | integer :: m_maxsol, m_minsol ! indexes of the maximum and minimum GCM sol taken in a bin for interpolation |
---|
| 189 | |
---|
| 190 | external LTmod ! declaration of the function LTmod |
---|
| 191 | real :: LTmod ! declaration of the type of the function LTmod |
---|
| 192 | integer :: LTcount ! nb of LT samples on which the interpolation is performed, |
---|
| 193 | ! for every LT interval (= numbin[lon,lat,alt,Ls]) |
---|
| 194 | |
---|
| 195 | integer :: solcount ! number of GCM sol integer values in one Ls interval |
---|
| 196 | real,dimension(:),allocatable :: int_sol_list ! list of the integer values of GCM sol |
---|
| 197 | real,dimension(:),allocatable :: sol_list ! list of the sol values used for the interpolation |
---|
| 198 | integer :: solerrcount ! nb of GCM missing values during interpolation, removed for the binning |
---|
| 199 | integer :: errcount = 0 ! total number of GCM missing values |
---|
| 200 | real :: solbinned_value ! extracted value averaged on sol samples, which is finally put in the output bin |
---|
| 201 | |
---|
| 202 | !------------------------ |
---|
| 203 | ! Extraction & loop indices |
---|
| 204 | real :: lon_val, lat_val, alt_val, Ls_val, LT_val ! where and when the output file is written at |
---|
| 205 | |
---|
| 206 | integer :: i,j,k,l ! loop iteration indices |
---|
| 207 | integer :: m ! sol binning loops iteration index |
---|
| 208 | integer :: v,vnot ! variable loops indices |
---|
| 209 | |
---|
[2558] | 210 | !------------------------ |
---|
| 211 | ! CDOD ratio (tau_GCM/tau_MCS) |
---|
| 212 | logical :: tau_ratio_ok ! true if it is possible to compute the CDOD ratio |
---|
[2404] | 213 | |
---|
[2558] | 214 | character(len=64)::out_dztauname,OBS_dztauname,OBS_tempname ! input variable names in files |
---|
| 215 | integer :: out_dztauid,OBS_dztauid,OBS_tempid ! input variable ID in files |
---|
| 216 | |
---|
| 217 | real,dimension(:,:,:,:),allocatable :: out_dztau,OBS_dztau ! the 4D dust opacities from outfile and obsfile |
---|
| 218 | real,dimension(:,:,:,:),allocatable :: OBS_temp ! temperature from obsfile |
---|
| 219 | |
---|
| 220 | real,dimension(:),allocatable :: OBS_plev ! pressure levels encompassing each Observer pressure levels |
---|
| 221 | |
---|
| 222 | real :: OBS_rho ! temporary OBS atm density computed from OBS temp and pressure |
---|
| 223 | real,parameter :: r_atm=191. ! Mars atmosphere specific gas constant |
---|
| 224 | real,parameter :: g=3.72 ! Mars gravity |
---|
| 225 | real :: out_tau,OBS_tau ! integrated dust columns from GCM (outfile) and OBS |
---|
| 226 | |
---|
| 227 | real,dimension(:,:,:),allocatable :: tau_ratio ! ratio out_tau/OBS_tau |
---|
| 228 | |
---|
| 229 | |
---|
[2404] | 230 | !=================================================================================== |
---|
| 231 | ! 1. OPENING OF THE FILES AND INITIALIZATION OF THE DIMENSIONS |
---|
| 232 | !=================================================================================== |
---|
| 233 | write(*,*) "Welcome in the MRO/MCS Observer Simulator program !" |
---|
| 234 | !=============================================================================== |
---|
| 235 | ! 1.1 MCS data file : obsfile |
---|
| 236 | !=============================================================================== |
---|
| 237 | !================================================================ |
---|
| 238 | ! 1.1.1 Open the Observer data file |
---|
| 239 | !================================================================ |
---|
| 240 | ! Ask the user to give a netcdf observation file |
---|
| 241 | WRITE(*,*) "-> Enter observation file name :" |
---|
| 242 | READ(*,*) obsfile |
---|
| 243 | |
---|
| 244 | status=NF90_OPEN(obsfile,nf90_nowrite,obsfid) ! nowrite mode=the program can only read the opened file |
---|
| 245 | error_text="Error: could not open file "//trim(obsfile) |
---|
| 246 | call status_check(status,error_text) |
---|
| 247 | |
---|
| 248 | !================================================================ |
---|
| 249 | ! 1.1.2 Get dimensions lon,lat,alt,time from the observation file |
---|
| 250 | !================================================================ |
---|
| 251 | ! OBS Latitude |
---|
| 252 | !-------------- |
---|
| 253 | status=nf90_inq_dimid(obsfid,"latitude",lat_dimid_obs) |
---|
| 254 | error_text="Failed to find Observer latitude dimension" |
---|
| 255 | call status_check(status,error_text) |
---|
| 256 | |
---|
| 257 | status=nf90_inquire_dimension(obsfid,lat_dimid_obs,len=OBSlatlen) |
---|
| 258 | error_text="Failed to find Observer latitude length" |
---|
| 259 | call status_check(status,error_text) |
---|
| 260 | allocate(OBSlat(OBSlatlen)) |
---|
| 261 | |
---|
| 262 | status=nf90_inq_varid(obsfid,"latitude",OBSvarid) |
---|
| 263 | error_text="Failed to find Observer latitude ID" |
---|
| 264 | call status_check(status,error_text) |
---|
| 265 | |
---|
| 266 | ! Read OBSlat |
---|
| 267 | status=NF90_GET_VAR(obsfid,OBSvarid,OBSlat) |
---|
| 268 | error_text="Failed to load OBSlat" |
---|
| 269 | call status_check(status,error_text) |
---|
| 270 | |
---|
| 271 | |
---|
| 272 | ! OBS Longitude |
---|
| 273 | !-------------- |
---|
| 274 | status=nf90_inq_dimid(obsfid,"longitude",lon_dimid_obs) |
---|
| 275 | error_text="Failed to find Observer longitude dimension" |
---|
| 276 | call status_check(status,error_text) |
---|
| 277 | |
---|
| 278 | status=nf90_inquire_dimension(obsfid,lon_dimid_obs,len=OBSlonlen) |
---|
| 279 | error_text="Failed to find Observer longitude length" |
---|
| 280 | call status_check(status,error_text) |
---|
| 281 | allocate(OBSlon(OBSlonlen)) |
---|
| 282 | |
---|
| 283 | status=nf90_inq_varid(obsfid,"longitude",OBSvarid) |
---|
| 284 | error_text="Failed to find Observer longitude ID" |
---|
| 285 | call status_check(status,error_text) |
---|
| 286 | |
---|
| 287 | ! Read OBSlon |
---|
| 288 | status=NF90_GET_VAR(obsfid,OBSvarid,OBSlon) |
---|
| 289 | error_text="Failed to load OBSlon" |
---|
| 290 | call status_check(status,error_text) |
---|
| 291 | |
---|
| 292 | |
---|
| 293 | ! OBS Time (Ls) |
---|
| 294 | !-------------- |
---|
| 295 | status=nf90_inq_dimid(obsfid,"time",time_dimid_obs) |
---|
| 296 | error_text="Failed to find Observer time (Ls) dimension" |
---|
| 297 | call status_check(status,error_text) |
---|
| 298 | |
---|
| 299 | status=nf90_inquire_dimension(obsfid,time_dimid_obs,len=OBSLslen) |
---|
| 300 | error_text="Failed to find Observer time (Ls) length" |
---|
| 301 | call status_check(status,error_text) |
---|
| 302 | allocate(OBSLs(OBSLslen)) |
---|
| 303 | |
---|
| 304 | status=nf90_inq_varid(obsfid,"time",OBSvarid) |
---|
| 305 | error_text="Failed to find Observer time (Ls) ID" |
---|
| 306 | call status_check(status,error_text) |
---|
| 307 | |
---|
| 308 | ! Read OBSLs |
---|
| 309 | status=NF90_GET_VAR(obsfid,OBSvarid,OBSLs) |
---|
| 310 | error_text="Failed to load OBSLs" |
---|
| 311 | call status_check(status,error_text) |
---|
| 312 | |
---|
| 313 | ! Get the observation timestep between bins |
---|
| 314 | OBSdeltaLs=OBSLs(2)-OBSLs(1) |
---|
| 315 | |
---|
| 316 | |
---|
| 317 | ! OBS Altitude |
---|
| 318 | !-------------- |
---|
| 319 | status=nf90_inq_dimid(obsfid,"altitude",alt_dimid_obs) |
---|
| 320 | error_text="Failed to find Observer altitude dimension" |
---|
| 321 | call status_check(status,error_text) |
---|
| 322 | |
---|
| 323 | status=nf90_inquire_dimension(obsfid,alt_dimid_obs,len=OBSaltlen) |
---|
| 324 | error_text="Failed to find Observer altitude length" |
---|
| 325 | call status_check(status,error_text) |
---|
| 326 | allocate(OBSalt(OBSaltlen)) |
---|
| 327 | |
---|
| 328 | status=nf90_inq_varid(obsfid,"altitude",OBSvarid) |
---|
| 329 | error_text="Failed to find Observer altitude ID" |
---|
| 330 | call status_check(status,error_text) |
---|
| 331 | |
---|
| 332 | ! Read OBSalt |
---|
| 333 | status=NF90_GET_VAR(obsfid,OBSvarid,OBSalt) |
---|
| 334 | error_text="Failed to load OBSalt" |
---|
| 335 | call status_check(status,error_text) |
---|
| 336 | |
---|
| 337 | ! Check altitude attribute "units" to find out altitude type and compare with the GCM file |
---|
| 338 | status=nf90_get_att(obsfid,OBSvarid,"units",units) |
---|
| 339 | error_text="Failed to load Observer altitude units attribute" |
---|
| 340 | call status_check(status,error_text) |
---|
| 341 | ! an unknown and invisible character is placed just after the unit's |
---|
| 342 | ! characters in the Observer file so we only take the first characters |
---|
| 343 | ! corresponding to the sought unit |
---|
| 344 | if (trim(units(1:2)).eq."Pa") then |
---|
| 345 | units="Pa" |
---|
| 346 | OBSalttype='p' ! pressure coordinate |
---|
| 347 | else if (trim(units(1:2)).eq."m") then |
---|
| 348 | units="m" |
---|
| 349 | OBSalttype='z' ! altitude coordinate |
---|
| 350 | else |
---|
| 351 | write(*,*)" I do not understand this unit ",trim(units)," for Observer altitude!" |
---|
| 352 | stop |
---|
| 353 | endif |
---|
| 354 | |
---|
| 355 | !=============================================================================== |
---|
| 356 | ! 1.2. GCM simulation file : gcmfile |
---|
| 357 | !=============================================================================== |
---|
| 358 | !================================================================ |
---|
| 359 | ! 1.2.1 Open the GCM simulation file |
---|
| 360 | !================================================================ |
---|
| 361 | ! Ask the user to give a netcdf input file |
---|
| 362 | write(*,*)"";WRITE(*,*) "-> Enter input file name (GCM simulation) :" |
---|
| 363 | READ(*,*) gcmfile |
---|
| 364 | |
---|
| 365 | ! Open GCM file |
---|
| 366 | status=NF90_OPEN(gcmfile,nf90_nowrite,gcmfid) |
---|
| 367 | ! nowrite mode=the program can only read the opened file |
---|
| 368 | error_text="Failed to open datafile "//trim(gcmfile) |
---|
| 369 | call status_check(status,error_text) |
---|
| 370 | |
---|
| 371 | !================================================================ |
---|
| 372 | ! 1.2.2 Get dimensions lon,lat,alt,time from the GCM file |
---|
| 373 | !================================================================ |
---|
| 374 | ! GCM Latitude |
---|
| 375 | !-------------- |
---|
| 376 | status=nf90_inq_dimid(gcmfid,"latitude",lat_dimid_gcm) |
---|
| 377 | error_text="Failed to find GCM latitude dimension" |
---|
| 378 | call status_check(status,error_text) |
---|
| 379 | |
---|
| 380 | status=nf90_inquire_dimension(gcmfid,lat_dimid_gcm,len=GCMlatlen) |
---|
| 381 | error_text="Failed to find GCM latitude length" |
---|
| 382 | call status_check(status,error_text) |
---|
| 383 | allocate(GCMlat(GCMlatlen)) |
---|
| 384 | |
---|
| 385 | status=nf90_inq_varid(gcmfid,"latitude",GCMvarid) |
---|
| 386 | error_text="Failed to find GCM latitude ID" |
---|
| 387 | call status_check(status,error_text) |
---|
| 388 | |
---|
| 389 | ! Read GCMlat |
---|
| 390 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCMlat) |
---|
| 391 | error_text="Failed to load GCMlat" |
---|
| 392 | call status_check(status,error_text) |
---|
| 393 | |
---|
| 394 | |
---|
| 395 | ! GCM Longitude |
---|
| 396 | !-------------- |
---|
| 397 | status=nf90_inq_dimid(gcmfid,"longitude",lon_dimid_gcm) |
---|
| 398 | error_text="Failed to find GCM longitude dimension" |
---|
| 399 | call status_check(status,error_text) |
---|
| 400 | |
---|
| 401 | status=nf90_inquire_dimension(gcmfid,lon_dimid_gcm,len=GCMlonlen) |
---|
| 402 | error_text="Failed to find GCM longitude length" |
---|
| 403 | call status_check(status,error_text) |
---|
| 404 | allocate(GCMlon(GCMlonlen)) |
---|
| 405 | |
---|
| 406 | status=nf90_inq_varid(gcmfid,"longitude",GCMvarid) |
---|
| 407 | error_text="Failed to find GCM longitude ID" |
---|
| 408 | call status_check(status,error_text) |
---|
| 409 | |
---|
| 410 | ! Read GCMlon |
---|
| 411 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCMlon) |
---|
| 412 | error_text="Failed to load GCMlon" |
---|
| 413 | call status_check(status,error_text) |
---|
| 414 | |
---|
| 415 | |
---|
| 416 | ! GCM Time |
---|
| 417 | !-------------- |
---|
| 418 | status=nf90_inq_dimid(gcmfid,"Time",time_dimid_gcm) |
---|
| 419 | error_text="Failed to find GCM time dimension" |
---|
| 420 | call status_check(status,error_text) |
---|
| 421 | |
---|
| 422 | status=nf90_inquire_dimension(gcmfid,time_dimid_gcm,len=GCMtimelen) |
---|
| 423 | error_text="Failed to find GCM time length" |
---|
| 424 | call status_check(status,error_text) |
---|
| 425 | allocate(GCMtime(GCMtimelen)) |
---|
| 426 | |
---|
| 427 | status=nf90_inq_varid(gcmfid,"Time",GCMvarid) |
---|
| 428 | error_text="Failed to find GCM time ID" |
---|
| 429 | call status_check(status,error_text) |
---|
| 430 | |
---|
| 431 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCMtime) |
---|
| 432 | error_text="Failed to load GCMtime" |
---|
| 433 | call status_check(status,error_text) |
---|
| 434 | |
---|
| 435 | ! is_stats ? |
---|
| 436 | IF ((GCMtimelen.eq.12).and.(GCMtime(1).eq.2.).and.(GCMtime(GCMtimelen).eq.24.)) then |
---|
| 437 | ! if GCM file is a stats, time is in LT at longitude 0° and not in sols at longitude 0° |
---|
| 438 | write(*,*)"The GCM file is recognized as a stats file." |
---|
| 439 | is_stats = .true. |
---|
| 440 | deallocate(GCMtime) |
---|
| 441 | GCMstatstimelen = GCMtimelen |
---|
| 442 | allocate(GCMstatstime(GCMstatstimelen)) |
---|
| 443 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCMstatstime) |
---|
| 444 | error_text="Failed to load GCMstatstime (LT at lon 0)" |
---|
| 445 | call status_check(status,error_text) |
---|
| 446 | ELSE |
---|
| 447 | write(*,*)"The GCM file is recognized as a diagfi/concatnc file." |
---|
| 448 | ENDIF |
---|
| 449 | |
---|
| 450 | ! Simulation time offset management |
---|
| 451 | WRITE(*,*) "Beginning date of the simulation file?" |
---|
| 452 | WRITE(*,*) "(i.e. number of sols since Ls=0 at the Time=0.0 in the GCM file)" |
---|
| 453 | READ(*,*) starttimeoffset |
---|
| 454 | if (.not.is_stats) then |
---|
| 455 | ! Add the offset to GCMtime(:) if the file is not a stats file |
---|
| 456 | GCMtime(:)=GCMtime(:)+starttimeoffset |
---|
| 457 | endif |
---|
| 458 | |
---|
| 459 | ! Check of temporal coherence between gcmfile & obsfile |
---|
| 460 | call ls2sol(OBSLs(OBSLslen),maxsol) ! maximum date considered |
---|
| 461 | call ls2sol(OBSLs(1),minsol) ! minimum date considered |
---|
| 462 | |
---|
| 463 | IF (.not.is_stats) then ! if it is a diagfi, we check the time coherence between the 2 files |
---|
| 464 | if ((maxsol.gt.maxval(GCMtime)).or.(minsol.lt.minval(GCMtime))) then |
---|
| 465 | write(*,*)"Error : obsfile temporal bounds exceed the GCM simulation bounds." |
---|
| 466 | write(*,*)"Please use a GCM file whose time interval contains the observation period." |
---|
| 467 | stop |
---|
| 468 | else |
---|
| 469 | write(*,*)"Both files are temporally coherent. The program continues..." |
---|
| 470 | endif |
---|
| 471 | |
---|
| 472 | ELSE ! if it is a stats, we create the array GCMtime array (in sols) covering the observation period |
---|
| 473 | ! and filled with the mean GCM day stored in stats.nc |
---|
| 474 | |
---|
| 475 | GCMtimelen = ((ceiling(maxsol)-floor(minsol)+1)+2) ! we add 2 days in the beginning and the end |
---|
| 476 | ! to be sure we cover the observation period |
---|
| 477 | allocate(GCMtime(GCMstatstimelen * GCMtimelen)) |
---|
| 478 | do l=1,GCMtimelen |
---|
| 479 | do m=1,GCMstatstimelen |
---|
| 480 | GCMtime(m+(l-1)*GCMstatstimelen) = (floor(minsol)-1) + (l-1) + GCMstatstime(m)/24. |
---|
| 481 | enddo |
---|
| 482 | enddo |
---|
| 483 | GCMtimelen = GCMstatstimelen * GCMtimelen |
---|
| 484 | write(*,*)"GCMtime has been created from the stats.nc time and the observation period. The program continues..." |
---|
| 485 | ENDIF |
---|
| 486 | |
---|
| 487 | |
---|
| 488 | ! GCM Altitude |
---|
| 489 | !-------------- |
---|
| 490 | status=nf90_inq_dimid(gcmfid,"altitude",alt_dimid_gcm) |
---|
| 491 | error_text="Failed to find GCM altitude dimension" |
---|
| 492 | call status_check(status,error_text) |
---|
| 493 | |
---|
| 494 | status=nf90_inquire_dimension(gcmfid,alt_dimid_gcm,len=GCMaltlen) |
---|
| 495 | error_text="Failed to find GCM altitude length" |
---|
| 496 | call status_check(status,error_text) |
---|
| 497 | allocate(GCMalt(GCMaltlen)) |
---|
| 498 | |
---|
| 499 | status=nf90_inq_varid(gcmfid,"altitude",GCMvarid) |
---|
| 500 | error_text="Failed to find GCM altitude ID" |
---|
| 501 | call status_check(status,error_text) |
---|
| 502 | |
---|
| 503 | ! Read GCMalt |
---|
| 504 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCMalt) |
---|
| 505 | error_text="Failed to load GCMalt" |
---|
| 506 | call status_check(status,error_text) |
---|
| 507 | |
---|
| 508 | ! Check altitude attribute "units" to find out altitude type |
---|
| 509 | status=nf90_get_att(gcmfid,GCMvarid,"units",units) |
---|
| 510 | error_text="Failed to load GCM altitude units attribute" |
---|
| 511 | call status_check(status,error_text) |
---|
| 512 | if (trim(units).eq."Pa") then |
---|
| 513 | GCMalttype='p' ! pressure coordinate |
---|
| 514 | else if (trim(units).eq."m") then |
---|
| 515 | GCMalttype='z' ! altitude coordinate |
---|
| 516 | else |
---|
| 517 | write(*,*)"I do not understand this unit ",trim(units)," for GCM altitude!" |
---|
| 518 | if (OBSalttype.eq.'p') then |
---|
| 519 | write(*,*)"Please use zrecast to put the altitude in the same type as the MCS file (pressure in Pa)" |
---|
| 520 | else if (OBSalttype.eq.'z') then |
---|
| 521 | write(*,*)"Please use zrecast to put the altitude in the same type as the MCS file (altitude in m)" |
---|
| 522 | endif |
---|
| 523 | stop |
---|
| 524 | endif |
---|
| 525 | IF(OBSalttype.ne.GCMalttype) then |
---|
| 526 | write(*,*)"Observer altitude type (", OBSalttype,") and ", & |
---|
| 527 | "GCM altitude type (",GCMalttype,") don't match!" |
---|
| 528 | stop |
---|
| 529 | ENDIF |
---|
| 530 | |
---|
| 531 | !=============================================================================== |
---|
| 532 | ! 1.3 Create the output file & initialize the coordinates |
---|
| 533 | !=============================================================================== |
---|
| 534 | ! Name of the outfile |
---|
| 535 | IF (.not.is_stats) then |
---|
| 536 | outfile=obsfile(1:index(obsfile, ".nc")-1)//"_GCMdiagfi.nc" |
---|
| 537 | ELSE |
---|
| 538 | outfile=obsfile(1:index(obsfile, ".nc")-1)//"_GCMstats.nc" |
---|
| 539 | ENDIF |
---|
| 540 | |
---|
| 541 | ! Creation of the outfile |
---|
| 542 | status=NF90_CREATE(outfile,nf90_clobber,outfid)!NB: clobber mode=overwrite any dataset with the same file name ! |
---|
| 543 | error_text="Error: could not create file "//trim(outfile) |
---|
| 544 | call status_check(status,error_text) |
---|
| 545 | write(*,*)"";WRITE(*,*)"-> Output file is: ",trim(outfile) |
---|
| 546 | |
---|
| 547 | ! Creation of the dimensions |
---|
| 548 | call inidim(outfid,OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen,OBSlon,OBSlat,OBSalt,OBSLs,OBSalttype,& |
---|
| 549 | lon_dimid_obs,lat_dimid_obs,alt_dimid_obs,time_dimid_obs) |
---|
| 550 | |
---|
| 551 | !write(*,*)"Dimensions ID in the outfile are :" |
---|
| 552 | !write(*,*)"lon_dimid=",lon_dimid_obs |
---|
| 553 | !write(*,*)"lat_dimid=",lat_dimid_obs |
---|
| 554 | !write(*,*)"alt_dimid=",alt_dimid_obs |
---|
| 555 | !write(*,*)"time_dimid=",time_dimid_obs |
---|
| 556 | |
---|
| 557 | |
---|
| 558 | |
---|
| 559 | !=================================================================================== |
---|
| 560 | ! 2. VARIABLES MANAGEMENT |
---|
| 561 | !=================================================================================== |
---|
| 562 | !=============================================================================== |
---|
| 563 | ! 2.1 List of the GCM variables to be interpolated |
---|
| 564 | !=============================================================================== |
---|
| 565 | !================================================================ |
---|
| 566 | ! 2.1.1 Read the GCM variables |
---|
| 567 | !================================================================ |
---|
[2525] | 568 | ! Initialize logicals |
---|
| 569 | dustok1 = .false. |
---|
| 570 | dustok2 = .false. |
---|
| 571 | dustok3 = .false. |
---|
| 572 | wiceok = .false. |
---|
| 573 | |
---|
[2404] | 574 | ! Get nbvarfile (total number of variables in the GCM file) |
---|
| 575 | status=NF90_INQUIRE(gcmfid,nVariables=nbvarfile) |
---|
| 576 | error_text="Error : Pb with nf90_inquire(gcmfid,nVariables=nbvarfile)" |
---|
| 577 | call status_check(status,error_text) |
---|
| 578 | |
---|
| 579 | ! List of variables that should not be processed |
---|
| 580 | notprocessed(1)='Time' |
---|
| 581 | notprocessed(2)='controle' |
---|
| 582 | notprocessed(3)='rlonu' |
---|
| 583 | notprocessed(4)='latitude' |
---|
| 584 | notprocessed(5)='longitude' |
---|
| 585 | notprocessed(6)='altitude' |
---|
| 586 | notprocessed(7)='rlatv' |
---|
| 587 | notprocessed(8)='aps' |
---|
| 588 | notprocessed(9)='bps' |
---|
| 589 | notprocessed(10)='ap' |
---|
| 590 | notprocessed(11)='bp' |
---|
| 591 | notprocessed(12)='cu' |
---|
| 592 | notprocessed(13)='cv' |
---|
| 593 | notprocessed(14)='aire' |
---|
| 594 | notprocessed(15)='phisinit' |
---|
| 595 | |
---|
| 596 | |
---|
| 597 | ! List of variables in the GCM file |
---|
| 598 | write(*,*)"" |
---|
| 599 | write(*,*)"List of variables in the GCM file :" |
---|
| 600 | Nnotprocessed=0 |
---|
| 601 | do v=1,nbvarfile |
---|
| 602 | status=NF90_INQUIRE_VARIABLE(gcmfid,v,name=GCMvarname) |
---|
| 603 | ! GCMvarname now contains the "name" of variable of ID # v |
---|
| 604 | var_ok=.true. |
---|
| 605 | do vnot=1,15 |
---|
| 606 | if (GCMvarname.eq.notprocessed(vnot)) then |
---|
| 607 | var_ok=.false. |
---|
| 608 | Nnotprocessed=Nnotprocessed+1 |
---|
| 609 | endif |
---|
| 610 | enddo |
---|
| 611 | if (var_ok) write(*,*) trim(GCMvarname) |
---|
| 612 | |
---|
| 613 | ! Detect if we can compute dust and wice opacities |
---|
| 614 | if (trim(GCMvarname).eq."dso") then |
---|
| 615 | dustok1 = .true. |
---|
| 616 | else if (trim(GCMvarname).eq."dsodust") then |
---|
| 617 | dustok2 = .true. |
---|
| 618 | else if (trim(GCMvarname).eq."dustq") then |
---|
| 619 | dustok3 = .true. |
---|
| 620 | else if (trim(GCMvarname).eq."h2o_ice") then |
---|
| 621 | wiceok = .true. |
---|
| 622 | endif |
---|
| 623 | enddo |
---|
| 624 | |
---|
| 625 | ! Nnotprocessed: # of variables that won't be processed |
---|
| 626 | ! nbvarfile: total # of variables in file |
---|
| 627 | ! +2: the dust and wice opacities |
---|
| 628 | allocate(gcm_vars(nbvarfile-Nnotprocessed+2),stat=status) |
---|
| 629 | if (status.ne.0) then |
---|
| 630 | write(*,*) "Error: failed allocation of gcm_vars(nbvarfile-Nnotprocessed+2)" |
---|
| 631 | write(*,*) " nbvarfile=",nbvarfile |
---|
| 632 | write(*,*) " Nnotprocessed=",Nnotprocessed |
---|
| 633 | stop |
---|
| 634 | endif |
---|
| 635 | |
---|
| 636 | ! List of variables to process |
---|
| 637 | write(*,*) |
---|
| 638 | write(*,*) "Which variables do you want to redistribute ?" |
---|
| 639 | write(*,*) "list of <variables> (separated by <Return>s)" |
---|
| 640 | write(*,*) "(an empty line , i.e: just <Return>, implies end of list)" |
---|
| 641 | write(*,*) "NB: this program handles only 4D netcdf variables for now" |
---|
| 642 | nbvar=0 |
---|
| 643 | read(*,'(a50)') GCMvarname |
---|
| 644 | do while ((GCMvarname.ne.' ').AND.(trim(GCMvarname).ne."all")) |
---|
| 645 | nbvar=nbvar+1 |
---|
| 646 | gcm_vars(nbvar)=GCMvarname |
---|
| 647 | read(*,'(a50)') GCMvarname |
---|
| 648 | enddo |
---|
| 649 | |
---|
| 650 | if (GCMvarname.eq."all") then |
---|
| 651 | nbvar=nbvarfile-Nnotprocessed |
---|
| 652 | do v=Nnotprocessed+1,nbvarfile |
---|
| 653 | status=nf90_inquire_variable(gcmfid,v,name=gcm_vars(v-Nnotprocessed)) |
---|
| 654 | enddo |
---|
| 655 | ! Variables names from the file are stored in gcm_vars() |
---|
| 656 | nbvar=nbvarfile-Nnotprocessed |
---|
| 657 | do v=1,nbvar |
---|
| 658 | status=nf90_inquire_variable(gcmfid,v+Nnotprocessed,name=gcm_vars(v)) |
---|
| 659 | write(*,'(a9,1x,i2,1x,a1,1x,a64)') "variable ",v,":",gcm_vars(v) |
---|
| 660 | enddo |
---|
| 661 | else if(nbvar==0) then |
---|
| 662 | write(*,*) "No variables to process in the GCM file... program stopped" |
---|
| 663 | stop |
---|
| 664 | endif ! of if (GCMvarname.eq."all") |
---|
| 665 | |
---|
| 666 | !================================================================ |
---|
| 667 | ! 2.1.2 Handle dust and wice opacities (first set-ups) |
---|
| 668 | !================================================================ |
---|
| 669 | ! 2nd part is in section 2.4.2 |
---|
| 670 | write(*,*) |
---|
| 671 | ! Load atmospheric density "rho" |
---|
| 672 | if (dustok1.or.dustok2.or.dustok3.or.wiceok) then |
---|
| 673 | ! Check that the GCM file contains that variable |
---|
| 674 | status=nf90_inq_varid(gcmfid,"rho",GCMvarid) |
---|
| 675 | if (status.ne.nf90_noerr) then |
---|
| 676 | write(*,*) "Failed to find variable rho in "//trim(gcmfile) |
---|
| 677 | write(*,*) "No computation of opacities..." |
---|
| 678 | dustok1 =.false. |
---|
| 679 | dustok2 =.false. |
---|
| 680 | dustok3 =.false. |
---|
| 681 | wiceok =.false. |
---|
| 682 | else |
---|
| 683 | ! Length allocation for each dimension of the 4D variable |
---|
[2558] | 684 | allocate(GCM_rho(GCMlonlen,GCMlatlen,GCMaltlen,GCMtimelen)) |
---|
[2404] | 685 | |
---|
| 686 | ! Load datasets |
---|
| 687 | if (.not.is_stats) then |
---|
[2558] | 688 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCM_rho) |
---|
[2404] | 689 | error_text="Failed to load rho" |
---|
| 690 | call status_check(status,error_text) |
---|
| 691 | else |
---|
| 692 | ! if it is a stats file, we load only the first sol, and then copy it to all the other sols |
---|
[2558] | 693 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCM_rho(:,:,:,1:GCMstatstimelen)) |
---|
[2404] | 694 | error_text="Failed to load rho" |
---|
| 695 | call status_check(status,error_text) |
---|
| 696 | ! write(*,*)"GCMstatstimelen = ", GCMstatstimelen |
---|
| 697 | ! write(*,*)"GCMtimelen = ", GCMtimelen |
---|
| 698 | do l=(GCMstatstimelen+1),GCMtimelen |
---|
| 699 | if (modulo(l,GCMstatstimelen).ne.0) then |
---|
[2558] | 700 | GCM_rho(:,:,:,l) = GCM_rho(:,:,:,modulo(l,GCMstatstimelen)) |
---|
[2404] | 701 | else ! if l is a multiple of GCMstatstimelen, since the index modulo(l,GCMstatstimelen)=0 |
---|
| 702 | ! doesn't exist, we make a special case |
---|
[2558] | 703 | GCM_rho(:,:,:,l) = GCM_rho(:,:,:,GCMstatstimelen) |
---|
[2404] | 704 | endif |
---|
| 705 | enddo |
---|
| 706 | endif |
---|
| 707 | write(*,*) "Variable rho loaded from the GCM file" |
---|
| 708 | endif |
---|
| 709 | endif ! dustok1.or.dustok2.or.dustok3.or.wiceok |
---|
| 710 | |
---|
| 711 | ! Dust and wice opacity booleans |
---|
| 712 | if (dustok1.or.dustok2.or.dustok3) then |
---|
| 713 | nbvar=nbvar+1 |
---|
| 714 | gcm_vars(nbvar)="dust" |
---|
| 715 | endif |
---|
| 716 | |
---|
| 717 | if (wiceok) then |
---|
| 718 | nbvar=nbvar+1 |
---|
| 719 | gcm_vars(nbvar)="wice" |
---|
| 720 | endif |
---|
| 721 | |
---|
[2540] | 722 | !write(*,*) "gcm_vars retrieved : ",gcm_vars(1:nbvar) |
---|
| 723 | |
---|
[2404] | 724 | !=============================================================================== |
---|
| 725 | ! 2.2 Definition of LT variables from obsfile to outfile |
---|
| 726 | !=============================================================================== |
---|
| 727 | ! --> the day/night loop begins here |
---|
| 728 | |
---|
| 729 | !******************** NOTA BENE (cf sections 2.2 and 4)************************* |
---|
| 730 | ! We execute the program a first time with the daytime values, and then a second |
---|
| 731 | ! time with the nighttime values. |
---|
| 732 | !******************************************************************************* |
---|
| 733 | |
---|
| 734 | dayornight = "dayside" ! we begin with daytime temperature |
---|
| 735 | write(*,*)"" ; write(*,*) "Beginning the 1st loop, on daytime values"; write(*,*)"" |
---|
| 736 | DAY_OR_NIGHT: DO ! (the end of the loop is in section 4.) |
---|
| 737 | |
---|
| 738 | SELECT CASE (dayornight) |
---|
| 739 | CASE ("dayside") |
---|
| 740 | OBSLTave_name = "dtimeave" |
---|
| 741 | OBSLTmax_name = "dtimemax" |
---|
| 742 | OBSLTmin_name = "dtimemin" |
---|
| 743 | CASE ("nightside") |
---|
| 744 | OBSLTave_name = "ntimeave" |
---|
| 745 | OBSLTmax_name = "ntimemax" |
---|
| 746 | OBSLTmin_name = "ntimemin" |
---|
| 747 | END SELECT |
---|
| 748 | |
---|
| 749 | !================================================================ |
---|
| 750 | ! 2.2.1 Average of Local Time in the OBS bins |
---|
| 751 | !================================================================ |
---|
| 752 | ! Read the OBS file |
---|
| 753 | !------------------ |
---|
| 754 | status=nf90_inq_varid(obsfid,trim(OBSLTave_name),LT_id) |
---|
| 755 | error_text="Failed to find Observer local time ("//trim(OBSLTave_name)//") ID in "//trim(obsfile) |
---|
| 756 | call status_check(status,error_text) |
---|
| 757 | status=nf90_inquire_variable(obsfid,LT_id,dimids=LTshape) |
---|
| 758 | error_text="Failed to get the dim shape of variable "//trim(OBSLTave_name) |
---|
| 759 | call status_check(status,error_text) |
---|
| 760 | |
---|
| 761 | ! Length allocation for each dimension of the 3D variable |
---|
| 762 | allocate(OBSLT(OBSlonlen,OBSlatlen,OBSLslen)) |
---|
| 763 | |
---|
| 764 | ! Load datasets |
---|
| 765 | status=NF90_GET_VAR(obsfid,LT_id,OBSLT) |
---|
| 766 | error_text="Failed to load "//trim(OBSLTave_name)//" from the obsfile" |
---|
| 767 | call status_check(status,error_text) |
---|
| 768 | write(*,*) trim(OBSLTave_name)," loaded from the obsfile" |
---|
| 769 | |
---|
| 770 | ! Get LT missing_value attribute |
---|
| 771 | status=nf90_get_att(obsfid,LT_id,"_FillValue",LTmiss_val) |
---|
| 772 | error_text="Failed to load missing_value attribute" |
---|
| 773 | call status_check(status,error_text) |
---|
| 774 | |
---|
| 775 | ! Create the variable in the outfile |
---|
| 776 | !----------------------------------- |
---|
| 777 | ! Switch to netcdf define mode |
---|
| 778 | status=nf90_redef(outfid) |
---|
| 779 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 780 | call status_check(status,error_text) |
---|
| 781 | |
---|
| 782 | ! Definition of the variable |
---|
| 783 | status=NF90_DEF_VAR(outfid,trim(OBSLTave_name),nf90_float,LTshape,LT_id) |
---|
| 784 | error_text="Error: could not define the variable "//trim(OBSLTave_name)//" in the outfile" |
---|
| 785 | call status_check(status,error_text) |
---|
| 786 | |
---|
| 787 | ! Write the attributes |
---|
| 788 | select case (dayornight) |
---|
| 789 | case ("dayside") |
---|
| 790 | status=nf90_put_att(outfid,LT_id,"long_name","Average local time in bin - day side [6h, 18h]") |
---|
| 791 | case ("nightside") |
---|
| 792 | status=nf90_put_att(outfid,LT_id,"long_name","Average local time in bin - night side [18h, 6h]") |
---|
| 793 | end select |
---|
| 794 | status=nf90_put_att(outfid,LT_id,"units","hours") |
---|
| 795 | status=nf90_put_att(outfid,LT_id,"_FillValue",LTmiss_val) |
---|
| 796 | |
---|
| 797 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 798 | status=nf90_enddef(outfid) |
---|
| 799 | error_text="Error: could not close the define mode of the outfile" |
---|
| 800 | call status_check(status,error_text) |
---|
| 801 | |
---|
| 802 | ! Write the data in the output file |
---|
| 803 | status = NF90_PUT_VAR(outfid, LT_id, OBSLT) ! write the MCS d/ntimeave as the output Local Time |
---|
| 804 | error_text="Error: could not write "//trim(OBSLTave_name)//" data in the outfile" |
---|
| 805 | call status_check(status,error_text) |
---|
| 806 | |
---|
| 807 | write(*,*)"Local Time (",trim(OBSLTave_name),") has been created in the outfile" |
---|
| 808 | write(*,'(" with missing_value attribute : ",1pe12.5)')LTmiss_val |
---|
| 809 | |
---|
| 810 | !================================================================ |
---|
| 811 | ! 2.2.2 Maximum of Local Time in the OBS bins |
---|
| 812 | !================================================================ |
---|
| 813 | ! Read the OBS file |
---|
| 814 | !------------------ |
---|
| 815 | status=nf90_inq_varid(obsfid,trim(OBSLTmax_name),LT_id) |
---|
| 816 | error_text="Failed to find Observer max local time ("//trim(OBSLTmax_name)//") ID in "//trim(obsfile) |
---|
| 817 | call status_check(status,error_text) |
---|
| 818 | status=nf90_inquire_variable(obsfid,LT_id,dimids=LTshape) |
---|
| 819 | error_text="Failed to get the dim shape of variable "//trim(OBSLTmax_name) |
---|
| 820 | call status_check(status,error_text) |
---|
| 821 | |
---|
| 822 | ! Length allocation for each dimension of the 3D variable |
---|
| 823 | allocate(OBSLTmax(OBSlonlen,OBSlatlen,OBSLslen)) |
---|
| 824 | |
---|
| 825 | ! Load datasets |
---|
| 826 | status=NF90_GET_VAR(obsfid,LT_id,OBSLTmax) |
---|
| 827 | error_text="Failed to load "//trim(OBSLTmax_name)//" from the obsfile" |
---|
| 828 | call status_check(status,error_text) |
---|
| 829 | write(*,*) trim(OBSLTmax_name)," loaded from the obsfile" |
---|
| 830 | |
---|
| 831 | !================================================================ |
---|
| 832 | ! 2.2.3 Minimum of Local Time in the OBS bins |
---|
| 833 | !================================================================ |
---|
| 834 | ! Read the OBS file |
---|
| 835 | !------------------ |
---|
| 836 | status=nf90_inq_varid(obsfid,trim(OBSLTmin_name),LT_id) |
---|
| 837 | error_text="Failed to find Observer min local time ("//trim(OBSLTmin_name)//") ID in "//trim(obsfile) |
---|
| 838 | call status_check(status,error_text) |
---|
| 839 | status=nf90_inquire_variable(obsfid,LT_id,dimids=LTshape) |
---|
| 840 | error_text="Failed to obtain information on variable "//trim(OBSLTmin_name) |
---|
| 841 | call status_check(status,error_text) |
---|
| 842 | |
---|
| 843 | ! Length allocation for each dimension of the 3D variable |
---|
| 844 | allocate(OBSLTmin(OBSlonlen,OBSlatlen,OBSLslen)) |
---|
| 845 | |
---|
| 846 | ! Load datasets |
---|
| 847 | status=NF90_GET_VAR(obsfid,LT_id,OBSLTmin) |
---|
| 848 | error_text="Failed to load "//trim(OBSLTmin_name)//" from the obsfile" |
---|
| 849 | call status_check(status,error_text) |
---|
| 850 | write(*,*) trim(OBSLTmin_name)," loaded from the obsfile" |
---|
| 851 | write(*,*)"" |
---|
| 852 | |
---|
| 853 | !=============================================================================== |
---|
| 854 | ! 2.3 Definition of numbin variables from obsfile to outfile |
---|
| 855 | !=============================================================================== |
---|
| 856 | !================================================================ |
---|
| 857 | ! 2.3.1 Number of values in the OBS "temp" bins |
---|
| 858 | !================================================================ |
---|
| 859 | SELECT CASE (dayornight) |
---|
| 860 | CASE ("dayside") |
---|
| 861 | numbin_name = "dnumbintemp" |
---|
| 862 | CASE ("nightside") |
---|
| 863 | numbin_name = "nnumbintemp" |
---|
| 864 | END SELECT |
---|
| 865 | |
---|
| 866 | ! Read the OBS file |
---|
| 867 | !------------------ |
---|
| 868 | status=nf90_inq_varid(obsfid,trim(numbin_name),numbin_id) |
---|
| 869 | error_text="Failed to find Observer nb of temp values in bin ("//trim(numbin_name)//")'s ID in "//trim(obsfile) |
---|
| 870 | call status_check(status,error_text) |
---|
| 871 | status=nf90_inquire_variable(obsfid,numbin_id,dimids=numbinshape) |
---|
| 872 | error_text="Failed to obtain information on variable "//trim(numbin_name) |
---|
| 873 | call status_check(status,error_text) |
---|
| 874 | |
---|
| 875 | ! Length allocation for each dimension of the 4D variable |
---|
| 876 | allocate(numbin(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 877 | |
---|
| 878 | ! Load datasets |
---|
| 879 | status=NF90_GET_VAR(obsfid,numbin_id,numbin) |
---|
| 880 | error_text="Failed to load "//trim(numbin_name)//" from the obsfile" |
---|
| 881 | call status_check(status,error_text) |
---|
| 882 | write(*,*) trim(numbin_name)," loaded from the obsfile" |
---|
| 883 | |
---|
| 884 | ! Create the variable in the outfile |
---|
| 885 | !----------------------------------- |
---|
| 886 | ! Switch to netcdf define mode |
---|
| 887 | status=nf90_redef(outfid) |
---|
| 888 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 889 | call status_check(status,error_text) |
---|
| 890 | |
---|
| 891 | ! Definition of the variable |
---|
| 892 | status=NF90_DEF_VAR(outfid,trim(numbin_name),nf90_float,numbinshape,numbin_id) |
---|
| 893 | error_text="Error: could not define the variable "//trim(numbin_name)//" in the outfile" |
---|
| 894 | call status_check(status,error_text) |
---|
| 895 | |
---|
| 896 | ! Write the attributes |
---|
| 897 | select case (dayornight) |
---|
| 898 | case ("dayside") |
---|
| 899 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of temp values in bin - day side") |
---|
| 900 | case ("nightside") |
---|
| 901 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of temp values in bin - night side") |
---|
| 902 | end select |
---|
| 903 | |
---|
| 904 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 905 | status=nf90_enddef(outfid) |
---|
| 906 | error_text="Error: could not close the define mode of the outfile" |
---|
| 907 | call status_check(status,error_text) |
---|
| 908 | |
---|
| 909 | ! Write the data in the output file |
---|
| 910 | status = NF90_PUT_VAR(outfid, numbin_id, numbin) |
---|
| 911 | error_text="Error: could not write "//trim(numbin_name)//" data in the outfile" |
---|
| 912 | call status_check(status,error_text) |
---|
| 913 | |
---|
| 914 | write(*,*)"Number of temp values in bin (",trim(numbin_name),") has been created in the outfile" |
---|
| 915 | write(*,*)"" |
---|
| 916 | deallocate(numbin) |
---|
| 917 | |
---|
| 918 | !================================================================ |
---|
| 919 | ! 2.3.2 Number of values in the OBS "dust" bins |
---|
| 920 | !================================================================ |
---|
| 921 | SELECT CASE (dayornight) |
---|
| 922 | CASE ("dayside") |
---|
| 923 | numbin_name = "dnumbindust" |
---|
| 924 | CASE ("nightside") |
---|
| 925 | numbin_name = "nnumbindust" |
---|
| 926 | END SELECT |
---|
| 927 | |
---|
| 928 | ! Read the OBS file |
---|
| 929 | !------------------ |
---|
| 930 | status=nf90_inq_varid(obsfid,trim(numbin_name),numbin_id) |
---|
| 931 | error_text="Failed to find Observer nb of dust values in bin ("//trim(numbin_name)//")'s ID in "//trim(obsfile) |
---|
| 932 | call status_check(status,error_text) |
---|
| 933 | status=nf90_inquire_variable(obsfid,numbin_id,dimids=numbinshape) |
---|
| 934 | error_text="Failed to obtain information on variable "//trim(numbin_name) |
---|
| 935 | call status_check(status,error_text) |
---|
| 936 | |
---|
| 937 | ! Length allocation for each dimension of the 4D variable |
---|
| 938 | allocate(numbin(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 939 | |
---|
| 940 | ! Load datasets |
---|
| 941 | status=NF90_GET_VAR(obsfid,numbin_id,numbin) |
---|
| 942 | error_text="Failed to load "//trim(numbin_name)//" from the obsfile" |
---|
| 943 | call status_check(status,error_text) |
---|
| 944 | write(*,*) trim(numbin_name)," loaded from the obsfile" |
---|
| 945 | |
---|
| 946 | ! Create the variable in the outfile |
---|
| 947 | !----------------------------------- |
---|
| 948 | ! Switch to netcdf define mode |
---|
| 949 | status=nf90_redef(outfid) |
---|
| 950 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 951 | call status_check(status,error_text) |
---|
| 952 | |
---|
| 953 | ! Definition of the variable |
---|
| 954 | status=NF90_DEF_VAR(outfid,trim(numbin_name),nf90_float,numbinshape,numbin_id) |
---|
| 955 | error_text="Error: could not define the variable "//trim(numbin_name)//" in the outfile" |
---|
| 956 | call status_check(status,error_text) |
---|
| 957 | |
---|
| 958 | ! Write the attributes |
---|
| 959 | select case (dayornight) |
---|
| 960 | case ("dayside") |
---|
| 961 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of dust values in bin - day side") |
---|
| 962 | case ("nightside") |
---|
| 963 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of dust values in bin - night side") |
---|
| 964 | end select |
---|
| 965 | |
---|
| 966 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 967 | status=nf90_enddef(outfid) |
---|
| 968 | error_text="Error: could not close the define mode of the outfile" |
---|
| 969 | call status_check(status,error_text) |
---|
| 970 | |
---|
| 971 | ! Write the data in the output file |
---|
| 972 | status = NF90_PUT_VAR(outfid, numbin_id, numbin) |
---|
| 973 | error_text="Error: could not write "//trim(numbin_name)//" data in the outfile" |
---|
| 974 | call status_check(status,error_text) |
---|
| 975 | |
---|
| 976 | write(*,*)"Number of dust values in bin (",trim(numbin_name),") has been created in the outfile" |
---|
| 977 | write(*,*)"" |
---|
| 978 | deallocate(numbin) |
---|
| 979 | |
---|
| 980 | !================================================================ |
---|
| 981 | ! 2.3.3 Number of values in the OBS "wice" bins |
---|
| 982 | !================================================================ |
---|
| 983 | SELECT CASE (dayornight) |
---|
| 984 | CASE ("dayside") |
---|
| 985 | numbin_name = "dnumbinwice" |
---|
| 986 | CASE ("nightside") |
---|
| 987 | numbin_name = "nnumbinwice" |
---|
| 988 | END SELECT |
---|
| 989 | |
---|
| 990 | ! Read the OBS file |
---|
| 991 | !------------------ |
---|
| 992 | status=nf90_inq_varid(obsfid,trim(numbin_name),numbin_id) |
---|
| 993 | error_text="Failed to find Observer nb of wice values in bin ("//trim(numbin_name)//")'s ID in "//trim(obsfile) |
---|
| 994 | call status_check(status,error_text) |
---|
| 995 | status=nf90_inquire_variable(obsfid,numbin_id,dimids=numbinshape) |
---|
| 996 | error_text="Failed to obtain information on variable "//trim(numbin_name) |
---|
| 997 | call status_check(status,error_text) |
---|
| 998 | |
---|
| 999 | ! Length allocation for each dimension of the 4D variable |
---|
| 1000 | allocate(numbin(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1001 | |
---|
| 1002 | ! Load datasets |
---|
| 1003 | status=NF90_GET_VAR(obsfid,numbin_id,numbin) |
---|
| 1004 | error_text="Failed to load "//trim(numbin_name)//" from the obsfile" |
---|
| 1005 | call status_check(status,error_text) |
---|
| 1006 | write(*,*) trim(numbin_name)," loaded from the obsfile" |
---|
| 1007 | |
---|
| 1008 | ! Create the variable in the outfile |
---|
| 1009 | !----------------------------------- |
---|
| 1010 | ! Switch to netcdf define mode |
---|
| 1011 | status=nf90_redef(outfid) |
---|
| 1012 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 1013 | call status_check(status,error_text) |
---|
| 1014 | |
---|
| 1015 | ! Definition of the variable |
---|
| 1016 | status=NF90_DEF_VAR(outfid,trim(numbin_name),nf90_float,numbinshape,numbin_id) |
---|
| 1017 | error_text="Error: could not define the variable "//trim(numbin_name)//" in the outfile" |
---|
| 1018 | call status_check(status,error_text) |
---|
| 1019 | |
---|
| 1020 | ! Write the attributes |
---|
| 1021 | select case (dayornight) |
---|
| 1022 | case ("dayside") |
---|
| 1023 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of wice values in bin - day side") |
---|
| 1024 | case ("nightside") |
---|
| 1025 | status=nf90_put_att(outfid,numbin_id,"long_name","Number of wice values in bin - night side") |
---|
| 1026 | end select |
---|
| 1027 | |
---|
| 1028 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 1029 | status=nf90_enddef(outfid) |
---|
| 1030 | error_text="Error: could not close the define mode of the outfile" |
---|
| 1031 | call status_check(status,error_text) |
---|
| 1032 | |
---|
| 1033 | ! Write the data in the output file |
---|
| 1034 | status = NF90_PUT_VAR(outfid, numbin_id, numbin) |
---|
| 1035 | error_text="Error: could not write "//trim(numbin_name)//" data in the outfile" |
---|
| 1036 | call status_check(status,error_text) |
---|
| 1037 | |
---|
| 1038 | write(*,*)"Number of wice values in bin (",trim(numbin_name),") has been created in the outfile" |
---|
| 1039 | write(*,*)"" |
---|
| 1040 | deallocate(numbin) |
---|
| 1041 | |
---|
| 1042 | |
---|
| 1043 | !=============================================================================== |
---|
| 1044 | ! 2.4 Opening of the GCM variables |
---|
| 1045 | !=============================================================================== |
---|
| 1046 | ! --> the var loop begins here |
---|
| 1047 | |
---|
| 1048 | VAR: DO v=1,nbvar ! LOOP ON ALL THE GCM VARIABLES TO PROCESS |
---|
| 1049 | ! (the end of the loop is in section 3.4) |
---|
| 1050 | |
---|
| 1051 | GCMvarname = gcm_vars(v) |
---|
| 1052 | |
---|
[2558] | 1053 | ! Detect the dust and wice opacities special cases |
---|
[2404] | 1054 | if (trim(GCMvarname).eq."dust") then |
---|
[2558] | 1055 | ! Default : dso > dsodust > dustq |
---|
[2404] | 1056 | if (dustok1) then ! "dso" is detected in gcmfile |
---|
| 1057 | GCMvarname="dso" |
---|
| 1058 | dustok2=.false. |
---|
| 1059 | dustok3=.false. |
---|
| 1060 | else if (dustok2) then ! "dsodust" is detected in gcmfile |
---|
| 1061 | GCMvarname="dsodust" |
---|
| 1062 | dustok3=.false. |
---|
| 1063 | else if (dustok3) then ! "dustq" is detected in gcmfile |
---|
| 1064 | GCMvarname="dustq" |
---|
| 1065 | endif |
---|
| 1066 | write(*,*) "Computing dust opacity..." |
---|
| 1067 | endif |
---|
[2558] | 1068 | |
---|
[2404] | 1069 | if (trim(GCMvarname).eq."wice") then ! "h2o_ice" detected in gcmfile |
---|
| 1070 | GCMvarname="h2o_ice" |
---|
| 1071 | write(*,*) "Computing water ice opacity..." |
---|
| 1072 | endif |
---|
| 1073 | |
---|
| 1074 | !================================================================ |
---|
| 1075 | ! 2.4.1 Generic reading of the variable |
---|
| 1076 | !================================================================ |
---|
| 1077 | ! Check that the GCM file contains that variable |
---|
| 1078 | status=nf90_inq_varid(gcmfid,trim(GCMvarname),GCMvarid) |
---|
| 1079 | if (status.ne.nf90_noerr) then |
---|
| 1080 | write(*,*) "Failed to find variable "//trim(GCMvarname)//" in "//trim(gcmfile) |
---|
| 1081 | write(*,*) "We'll skip that variable..." |
---|
| 1082 | CYCLE VAR ! go directly to the next variable |
---|
| 1083 | endif |
---|
| 1084 | |
---|
| 1085 | ! Sanity checks on the variable |
---|
| 1086 | status=nf90_inquire_variable(gcmfid,GCMvarid,ndims=nbdim,dimids=GCMvarshape) |
---|
| 1087 | error_text="Failed to obtain information on variable "//trim(GCMvarname) |
---|
| 1088 | call status_check(status,error_text) |
---|
| 1089 | |
---|
| 1090 | ! Check that it is a 4D variable |
---|
| 1091 | if (nbdim.ne.4) then |
---|
| 1092 | write(*,*) "Error:",trim(GCMvarname)," is not a 4D variable" |
---|
| 1093 | write(*,*) "We'll skip that variable...";write(*,*)"" |
---|
| 1094 | CYCLE VAR ! go directly to the next variable |
---|
| 1095 | endif |
---|
| 1096 | ! Check that its dimensions are indeed lon,lat,alt,time (in the right order) |
---|
| 1097 | if (GCMvarshape(1).ne.lon_dimid_gcm) then |
---|
| 1098 | write(*,*) "Error, expected first dimension of ",trim(GCMvarname)," to be longitude!" |
---|
| 1099 | write(*,*) "We'll skip that variable..." |
---|
| 1100 | CYCLE VAR ! go directly to the next variable |
---|
| 1101 | endif |
---|
| 1102 | if (GCMvarshape(2).ne.lat_dimid_gcm) then |
---|
| 1103 | write(*,*) "Error, expected second dimension of ",trim(GCMvarname)," to be latitude!" |
---|
| 1104 | write(*,*) "We'll skip that variable..." |
---|
| 1105 | CYCLE VAR ! go directly to the next variable |
---|
| 1106 | endif |
---|
| 1107 | if (GCMvarshape(3).ne.alt_dimid_gcm) then |
---|
| 1108 | write(*,*) "Error, expected third dimension of ",trim(GCMvarname)," to be altitude!" |
---|
| 1109 | write(*,*) "We'll skip that variable..." |
---|
| 1110 | CYCLE VAR ! go directly to the next variable |
---|
| 1111 | endif |
---|
| 1112 | if (GCMvarshape(4).ne.time_dimid_gcm) then |
---|
| 1113 | write(*,*) "Error, expected fourth dimension of ",trim(GCMvarname)," to be time!" |
---|
| 1114 | write(*,*) "We'll skip that variable..." |
---|
| 1115 | CYCLE VAR ! go directly to the next variable |
---|
| 1116 | endif |
---|
| 1117 | |
---|
| 1118 | ! Length allocation for each dimension of the 4D variable |
---|
| 1119 | allocate(GCM_var(GCMlonlen,GCMlatlen,GCMaltlen,GCMtimelen)) |
---|
| 1120 | |
---|
| 1121 | ! Load datasets |
---|
| 1122 | if (.not.is_stats) then |
---|
| 1123 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCM_var) |
---|
| 1124 | error_text="Failed to load "//trim(GCMvarname) |
---|
| 1125 | call status_check(status,error_text) |
---|
| 1126 | else |
---|
| 1127 | ! if it is a stats file, we load only the first sol, and then copy it to all the other sols |
---|
| 1128 | status=NF90_GET_VAR(gcmfid,GCMvarid,GCM_var(:,:,:,1:GCMstatstimelen)) |
---|
| 1129 | error_text="Failed to load "//trim(GCMvarname) |
---|
| 1130 | call status_check(status,error_text) |
---|
| 1131 | ! write(*,*)"GCMstatstimelen = ", GCMstatstimelen |
---|
| 1132 | ! write(*,*)"GCMtimelen = ", GCMtimelen |
---|
| 1133 | do l=(GCMstatstimelen+1),GCMtimelen |
---|
| 1134 | if (modulo(l,GCMstatstimelen).ne.0) then |
---|
| 1135 | GCM_var(:,:,:,l) = GCM_var(:,:,:,modulo(l,GCMstatstimelen)) |
---|
| 1136 | else ! if l is a multiple of GCMstatstimelen, since the index modulo(l,GCMstatstimelen)=0 |
---|
| 1137 | ! doesn't exist, we make a special case |
---|
| 1138 | GCM_var(:,:,:,l) = GCM_var(:,:,:,GCMstatstimelen) |
---|
| 1139 | endif |
---|
| 1140 | enddo |
---|
| 1141 | endif |
---|
| 1142 | write(*,*) "Variable ",trim(GCMvarname)," loaded from the GCM file" |
---|
| 1143 | |
---|
| 1144 | ! Get dataset's missing_value attribute |
---|
| 1145 | status=nf90_get_att(gcmfid,GCMvarid,"missing_value",GCMmiss_val) |
---|
| 1146 | error_text="Failed to load missing_value attribute" |
---|
| 1147 | call status_check(status,error_text) |
---|
| 1148 | |
---|
| 1149 | ! Get other variable's attributes |
---|
| 1150 | status=nf90_get_att(gcmfid,GCMvarid,"long_name",long_name) |
---|
| 1151 | if (status.ne.nf90_noerr) then |
---|
| 1152 | ! if no attribute "long_name", try "title" |
---|
| 1153 | status=nf90_get_att(gcmfid,GCMvarid,"title",long_name) |
---|
| 1154 | endif |
---|
| 1155 | status=nf90_get_att(gcmfid,GCMvarid,"units",units) |
---|
| 1156 | |
---|
| 1157 | !================================================================ |
---|
| 1158 | ! 2.4.2 Handle dust and wice opacities (second part) |
---|
| 1159 | !================================================================ |
---|
| 1160 | ! DUST |
---|
| 1161 | !----- |
---|
| 1162 | if (trim(gcm_vars(v)).eq."dust") then |
---|
| 1163 | |
---|
| 1164 | IF (dustok1.or.dustok2) THEN |
---|
| 1165 | ! Dust opacity computed from its density-scaled opacity |
---|
| 1166 | ! write(*,*)long_name(index(long_name,"(")+1:index(long_name,")")-1) |
---|
| 1167 | do i=1,GCMlonlen |
---|
| 1168 | do j=1,GCMlatlen |
---|
| 1169 | do k=1,GCMaltlen |
---|
| 1170 | do l=1,GCMtimelen |
---|
| 1171 | if (GCM_var(i,j,k,l).ne.GCMmiss_val) then |
---|
| 1172 | ! Multiply by rho to have opacity [1/km] |
---|
[2558] | 1173 | GCM_var(i,j,k,l) = GCM_var(i,j,k,l) * GCM_rho(i,j,k,l) *1000. |
---|
[2404] | 1174 | |
---|
| 1175 | if (long_name(index(long_name,"(")+1:index(long_name,")")-1).eq."TES") then |
---|
| 1176 | ! The density-scaled opacity was calibrated on TES wavelength (9.3um) |
---|
| 1177 | ! so we must recalibrate it to MCS wavelength (21.6um) using recalibration |
---|
| 1178 | ! coefficients from Montabone et al. 2015, section 2.3 |
---|
| 1179 | GCM_var(i,j,k,l) = 1.3/2.7 * GCM_var(i,j,k,l) |
---|
| 1180 | endif |
---|
| 1181 | endif |
---|
| 1182 | enddo |
---|
| 1183 | enddo |
---|
| 1184 | enddo |
---|
| 1185 | enddo |
---|
| 1186 | |
---|
| 1187 | long_name = "IR Dust opacity (from DSO)" |
---|
| 1188 | |
---|
| 1189 | ELSE IF (dustok3) THEN |
---|
| 1190 | ! Dust opacity computed from its mass mixing ratio |
---|
| 1191 | do i=1,GCMlonlen |
---|
| 1192 | do j=1,GCMlatlen |
---|
| 1193 | do k=1,GCMaltlen |
---|
| 1194 | do l=1,GCMtimelen |
---|
| 1195 | if (GCM_var(i,j,k,l).ne.GCMmiss_val) then |
---|
| 1196 | ! Opacity is computed from the equation of Heavens et al. 2014, section 2.3, |
---|
| 1197 | ! assuming a rho_dust=3000 kg/m3 and an effective radius reff=1.06 microns |
---|
| 1198 | ! + the opacity is here in 1/km and the MMR in kg/kg |
---|
[2558] | 1199 | GCM_var(i,j,k,l) = GCM_var(i,j,k,l) * GCM_rho(i,j,k,l) / 0.012 * 1000 |
---|
[2404] | 1200 | endif |
---|
| 1201 | enddo |
---|
| 1202 | enddo |
---|
| 1203 | enddo |
---|
| 1204 | enddo |
---|
| 1205 | |
---|
| 1206 | long_name = "IR Dust opacity (from MMR)" |
---|
| 1207 | |
---|
| 1208 | ENDIF |
---|
| 1209 | |
---|
| 1210 | GCMvarname = gcm_vars(v) ! reput the right name in GCMvarname |
---|
| 1211 | units = "opacity/km" |
---|
| 1212 | endif ! trim(gcm_vars(v)).eq."dust" |
---|
| 1213 | |
---|
| 1214 | |
---|
| 1215 | ! WICE |
---|
| 1216 | !----- |
---|
| 1217 | if (trim(gcm_vars(v)).eq."wice") then |
---|
| 1218 | ! Water ice opacity computed from its mass mixing ratio |
---|
| 1219 | do i=1,GCMlonlen |
---|
| 1220 | do j=1,GCMlatlen |
---|
| 1221 | do k=1,GCMaltlen |
---|
| 1222 | do l=1,GCMtimelen |
---|
| 1223 | if (GCM_var(i,j,k,l).ne.GCMmiss_val) then |
---|
| 1224 | ! Opacity at MCS wavelength (11.9um) is computed from an equation |
---|
| 1225 | ! similar to the one of Heavens et al. 2014, section 2.3. |
---|
| 1226 | ! We assume a rho_wice=920 kg/m3, an effective radius reff=3um, |
---|
| 1227 | ! an extinction coefficient Qext(wvl,reff)=1.54471 |
---|
[2558] | 1228 | GCM_var(i,j,k,l) = 750*1.54471* GCM_var(i,j,k,l) * GCM_rho(i,j,k,l) / (920*3e-6) |
---|
[2404] | 1229 | endif |
---|
| 1230 | enddo |
---|
| 1231 | enddo |
---|
| 1232 | enddo |
---|
| 1233 | enddo |
---|
| 1234 | |
---|
| 1235 | long_name = "IR Water ice opacity (from MMR)" |
---|
| 1236 | |
---|
| 1237 | GCMvarname = gcm_vars(v) ! reput the right name in GCMvarname |
---|
| 1238 | units = "opacity/km" |
---|
| 1239 | endif ! trim(gcm_vars(v)).eq."wice" |
---|
| 1240 | |
---|
| 1241 | !=============================================================================== |
---|
| 1242 | ! 2.5 Opening of the associated MCS variables |
---|
| 1243 | !=============================================================================== |
---|
| 1244 | ! Observer variables to extract : |
---|
[2443] | 1245 | IF ((index(GCMvarname,"dust").ne.0).or.(index(GCMvarname,"dso").ne.0)) THEN |
---|
| 1246 | ! if the variable name contains "dust" or "dso". Especially for the targeted variables : |
---|
| 1247 | ! dustq,dustN,dsodust,reffdust,opadust, and their equivalents for stormdust & topdust |
---|
[2404] | 1248 | OBSvarname = "dust" |
---|
| 1249 | numbin_name = "numbindust" |
---|
| 1250 | ELSE IF ((trim(GCMvarname).eq."h2o_ice").or.(trim(GCMvarname).eq."wice") & |
---|
[2443] | 1251 | .or.(trim(GCMvarname).eq."reffice").or.(trim(GCMvarname).eq."opawice")) THEN |
---|
[2404] | 1252 | OBSvarname = "wice" |
---|
| 1253 | numbin_name = "numbinwice" |
---|
| 1254 | ELSE ! default case is temp binning, since it contains the most values |
---|
| 1255 | OBSvarname = "temp" |
---|
| 1256 | numbin_name = "numbintemp" |
---|
| 1257 | ENDIF |
---|
| 1258 | |
---|
| 1259 | SELECT CASE (dayornight) |
---|
| 1260 | CASE ("dayside") |
---|
| 1261 | OBSvarname = "d"//OBSvarname |
---|
| 1262 | numbin_name = "d"//numbin_name |
---|
| 1263 | CASE ("nightside") |
---|
| 1264 | OBSvarname = "n"//OBSvarname |
---|
| 1265 | numbin_name = "n"//numbin_name |
---|
| 1266 | END SELECT |
---|
| 1267 | |
---|
| 1268 | !================================================================ |
---|
| 1269 | ! 2.5.1 MCS reference variable (for the missing values) |
---|
| 1270 | !================================================================ |
---|
| 1271 | ! Check that the observation file contains that variable |
---|
| 1272 | status=nf90_inq_varid(obsfid,trim(OBSvarname),OBSvarid) |
---|
| 1273 | error_text="Failed to find variable "//trim(OBSvarname)//" in "//trim(obsfile) |
---|
| 1274 | call status_check(status,error_text) |
---|
| 1275 | |
---|
| 1276 | ! Sanity checks on the variable |
---|
| 1277 | status=nf90_inquire_variable(obsfid,OBSvarid,ndims=nbdim,dimids=OBSvarshape) |
---|
| 1278 | error_text="Failed to obtain information on variable "//trim(OBSvarname) |
---|
| 1279 | call status_check(status,error_text) |
---|
| 1280 | |
---|
| 1281 | ! Check that it is a 4D variable |
---|
| 1282 | if (nbdim.ne.4) then |
---|
| 1283 | write(*,*) "Error, expected a 4D (lon-lat-alt-time) variable for ",trim(OBSvarname) |
---|
| 1284 | stop |
---|
| 1285 | endif |
---|
| 1286 | ! Check that its dimensions are indeed lon,lat,alt,time (in the right order) |
---|
| 1287 | if (OBSvarshape(1).ne.lon_dimid_obs) then |
---|
| 1288 | write(*,*) "Error, expected first dimension of ",trim(OBSvarname)," to be longitude!" |
---|
| 1289 | stop |
---|
| 1290 | endif |
---|
| 1291 | if (OBSvarshape(2).ne.lat_dimid_obs) then |
---|
| 1292 | write(*,*) "Error, expected second dimension of ",trim(OBSvarname)," to be latitude!" |
---|
| 1293 | stop |
---|
| 1294 | endif |
---|
| 1295 | if (OBSvarshape(3).ne.alt_dimid_obs) then |
---|
| 1296 | write(*,*) "Error, expected third dimension of ",trim(OBSvarname)," to be altitude!" |
---|
| 1297 | stop |
---|
| 1298 | endif |
---|
| 1299 | if (OBSvarshape(4).ne.time_dimid_obs) then |
---|
| 1300 | write(*,*) "Error, expected fourth dimension of ",trim(OBSvarname)," to be time!" |
---|
| 1301 | stop |
---|
| 1302 | endif |
---|
| 1303 | |
---|
| 1304 | ! Length allocation for each dimension of the 4D variable |
---|
| 1305 | allocate(OBS_var(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1306 | |
---|
| 1307 | ! Load datasets |
---|
| 1308 | status=NF90_GET_VAR(obsfid,OBSvarid,OBS_var) |
---|
| 1309 | error_text="Failed to load "//trim(OBSvarname)//" from the obsfile" |
---|
| 1310 | call status_check(status,error_text) |
---|
| 1311 | write(*,*) trim(OBSvarname)," loaded from the obsfile as reference variable" |
---|
| 1312 | |
---|
| 1313 | ! Get OBS_var missing_value attribute |
---|
| 1314 | status=nf90_get_att(obsfid,OBSvarid,"_FillValue",OBSmiss_val) |
---|
| 1315 | error_text="Failed to load missing_value attribute" |
---|
| 1316 | call status_check(status,error_text) |
---|
| 1317 | |
---|
| 1318 | !================================================================ |
---|
| 1319 | ! 2.5.2 Number of values in the OBS bin (for the sol binning) |
---|
| 1320 | !================================================================ |
---|
| 1321 | ! Check that the observation file contains that variable |
---|
| 1322 | status=nf90_inq_varid(obsfid,trim(numbin_name),numbin_id) |
---|
| 1323 | |
---|
| 1324 | ! Checks have already been done in section 2.3 |
---|
| 1325 | |
---|
| 1326 | ! Length allocation for each dimension of the 4D variable |
---|
| 1327 | allocate(numbin(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1328 | |
---|
| 1329 | ! Load datasets |
---|
| 1330 | status=NF90_GET_VAR(obsfid,numbin_id,numbin) |
---|
| 1331 | |
---|
| 1332 | !=============================================================================== |
---|
| 1333 | ! 2.6 Definition of GCM variables in outfile |
---|
| 1334 | !=============================================================================== |
---|
| 1335 | ! Switch to netcdf define mode |
---|
| 1336 | status=nf90_redef(outfid) |
---|
| 1337 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 1338 | call status_check(status,error_text) |
---|
| 1339 | |
---|
| 1340 | ! Definition of the variable |
---|
| 1341 | SELECT CASE (dayornight) |
---|
| 1342 | CASE ("dayside") |
---|
| 1343 | outvarname = "d"//GCMvarname |
---|
| 1344 | CASE ("nightside") |
---|
| 1345 | outvarname = "n"//GCMvarname |
---|
| 1346 | END SELECT |
---|
| 1347 | status=NF90_DEF_VAR(outfid,trim(outvarname),nf90_float,OBSvarshape,outvarid) |
---|
| 1348 | error_text="Error: could not define the variable "//trim(outvarname)//" in the outfile" |
---|
| 1349 | call status_check(status,error_text) |
---|
| 1350 | |
---|
| 1351 | ! Write the attributes |
---|
| 1352 | SELECT CASE (dayornight) |
---|
| 1353 | CASE ("dayside") |
---|
| 1354 | long_name = trim(long_name)//" - day side" |
---|
| 1355 | status=nf90_put_att(outfid,outvarid,"long_name",long_name) |
---|
| 1356 | CASE ("nightside") |
---|
| 1357 | long_name = trim(long_name)//" - night side" |
---|
| 1358 | status=nf90_put_att(outfid,outvarid,"long_name",long_name) |
---|
| 1359 | END SELECT |
---|
| 1360 | status=nf90_put_att(outfid,outvarid,"units",units) |
---|
| 1361 | status=nf90_put_att(outfid,outvarid,"_FillValue",OBSmiss_val) |
---|
| 1362 | comment = "Reference numbin: "//trim(numbin_name) |
---|
| 1363 | status=nf90_put_att(outfid,outvarid,"comment",comment) |
---|
| 1364 | |
---|
| 1365 | write(*,*)trim(outvarname)," has been created in the outfile" |
---|
| 1366 | write(*,'(" with missing_value attribute : ",1pe12.5)')OBSmiss_val |
---|
| 1367 | write(*,*)"" |
---|
| 1368 | |
---|
| 1369 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 1370 | status=nf90_enddef(outfid) |
---|
| 1371 | error_text="Error: could not close the define mode of the outfile" |
---|
| 1372 | call status_check(status,error_text) |
---|
| 1373 | |
---|
| 1374 | allocate(outvar(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1375 | |
---|
| 1376 | |
---|
| 1377 | !=================================================================================== |
---|
| 1378 | ! 3. EXTRACTION OF THE VARIABLE |
---|
| 1379 | !=================================================================================== |
---|
| 1380 | !=============================================================================== |
---|
| 1381 | ! 3.1 Do some checks and preparations before the extraction |
---|
| 1382 | !=============================================================================== |
---|
| 1383 | Ls: do l=1,OBSLslen ! loop on all observed Solar longitudes |
---|
| 1384 | Ls_val=OBSLs(l) ! Ls_val=center of the output bin |
---|
| 1385 | if ((Ls_val.lt.0.).or.(Ls_val.gt.360.)) then |
---|
| 1386 | write(*,*) "Unexpected value for OBSLs: ",Ls_val |
---|
| 1387 | stop |
---|
| 1388 | endif |
---|
| 1389 | |
---|
| 1390 | ! Convert the Ls bin into a sol interval on which the binning is done : |
---|
| 1391 | !---------------------------------------------------------------------- |
---|
| 1392 | !-> get the index of the maximum value of GCM sol (m_maxsol) that is lower than Ls bin's superior bound (maxsol) |
---|
| 1393 | call ls2sol(Ls_val+OBSdeltaLs/2.,maxsol) |
---|
| 1394 | m_maxsol=1 |
---|
| 1395 | do while (((GCMtime(m_maxsol)+0.5).lt.maxsol).AND.(m_maxsol.le.(GCMtimelen-1))) |
---|
| 1396 | !! The +0.5 is there to take into account the whole planet (lon=[-180°;180°]) and not just the lon=0° from the GCM |
---|
| 1397 | m_maxsol=m_maxsol+1 |
---|
| 1398 | enddo |
---|
| 1399 | !-> get the index of the minimum value of GCM sol (m_minsol) that is greater than Ls bin's inferior bound (minsol) |
---|
| 1400 | call ls2sol(Ls_val-OBSdeltaLs/2.,minsol) |
---|
| 1401 | m_minsol=1 |
---|
| 1402 | do while (((GCMtime(m_minsol)-0.5).le.minsol).AND.(m_minsol.le.(GCMtimelen-1))) |
---|
| 1403 | !! Same comment for the -0.5 |
---|
| 1404 | m_minsol=m_minsol+1 |
---|
| 1405 | enddo |
---|
| 1406 | if (m_minsol.gt.m_maxsol) then |
---|
| 1407 | write(*,*) "No value in gcmfile between sol=",minsol," and sol=",maxsol," (Ls=",Ls_val,"°)" |
---|
| 1408 | ! Write a missing_value to output |
---|
| 1409 | outvar(:,:,:,l) = OBSmiss_val |
---|
| 1410 | CYCLE Ls ! go directly to the next Ls |
---|
| 1411 | endif |
---|
| 1412 | ! Get all the integer values of GCM sols that fit in this interval |
---|
| 1413 | solcount=floor(GCMtime(m_maxsol))-ceiling(GCMtime(m_minsol))+1 |
---|
| 1414 | ! sols that are not fully in the interval are not counted |
---|
| 1415 | allocate(int_sol_list(solcount)) |
---|
| 1416 | ! write(*,*) "GCMminsol=", GCMtime(m_minsol) |
---|
| 1417 | ! write(*,*)"GCMmaxsol=", GCMtime(m_maxsol) |
---|
| 1418 | do m=1,solcount |
---|
| 1419 | int_sol_list(m)=ceiling(GCMtime(m_minsol)) + m-1 |
---|
| 1420 | enddo |
---|
| 1421 | ! write(*,*)"int_sol_list=",int_sol_list |
---|
| 1422 | |
---|
| 1423 | |
---|
| 1424 | latitude: do j=1,OBSlatlen ! loop on all observed latitudes |
---|
| 1425 | lat_val=OBSlat(j) |
---|
| 1426 | if ((lat_val.lt.-90.).or.(lat_val.gt.90.)) then |
---|
| 1427 | write(*,*) "Unexpected value for OBSlat: ",lat_val |
---|
| 1428 | stop |
---|
| 1429 | endif |
---|
| 1430 | |
---|
| 1431 | longitude: do i=1,OBSlonlen ! loop on all observed longitudes |
---|
| 1432 | lon_val=OBSlon(i) |
---|
| 1433 | if ((lon_val.lt.-360.).or.(lon_val.gt.360.)) then |
---|
| 1434 | write(*,*) "Unexpected value for lon_val: ",lon_val |
---|
| 1435 | stop |
---|
| 1436 | endif |
---|
| 1437 | ! We want lon_val in [-180:180] for the subroutine extraction |
---|
| 1438 | if (lon_val.lt.-180.) lon_val=lon_val+360. |
---|
| 1439 | if (lon_val.gt.180.) lon_val=lon_val-360. |
---|
| 1440 | |
---|
| 1441 | LT_val=OBSLT(i,j,l) ! find the Observer average LT value at bin(lon_val, lat_val, Ls_val) |
---|
| 1442 | |
---|
| 1443 | if ((LT_val.lt.0.).or.(LT_val.gt.24.)) then |
---|
| 1444 | if (LT_val.eq.LTmiss_val) then |
---|
| 1445 | ! write(*,*) "Missing value in obsfile for LT_val" |
---|
| 1446 | ! Write a missing_value to output |
---|
| 1447 | outvar(i,j,:,l) = OBSmiss_val |
---|
| 1448 | CYCLE longitude ! go directly to the next longitude |
---|
| 1449 | else |
---|
| 1450 | write(*,*) "Unexpected value for LT_val: ",LT_val |
---|
| 1451 | stop |
---|
| 1452 | endif |
---|
| 1453 | endif |
---|
| 1454 | |
---|
| 1455 | altitude: do k=1,OBSaltlen ! loop on all observed altitudes |
---|
| 1456 | alt_val=OBSalt(k) |
---|
| 1457 | if (OBS_var(i,j,k,l).eq.OBSmiss_val) then |
---|
| 1458 | ! write(*,*) "Missing value in obsfile for ",OBSvarname |
---|
| 1459 | ! Write a missing_value to output |
---|
| 1460 | outvar(i,j,k,l) = OBSmiss_val |
---|
| 1461 | CYCLE altitude ! go directly to the next altitude |
---|
| 1462 | endif |
---|
| 1463 | |
---|
| 1464 | !=============================================================================== |
---|
| 1465 | ! 3.2 Compute GCM sol date corresponding to Observer Ls (via m_(min/max)sol) |
---|
| 1466 | ! and LT (via OBSLT(min/max)) |
---|
| 1467 | !=============================================================================== |
---|
| 1468 | LTcount=floor(numbin(i,j,k,l)) ! find the Observer number of temp values |
---|
| 1469 | ! at bin(lon_val,lat_val,alt_val,Ls_val) |
---|
| 1470 | if (LTcount.eq.0.) then |
---|
| 1471 | ! Write a missing_value to output |
---|
| 1472 | outvar(i,j,k,l) = OBSmiss_val |
---|
| 1473 | CYCLE altitude ! go directly to the next altitude |
---|
| 1474 | endif |
---|
| 1475 | if (LTcount.lt.0.) then |
---|
| 1476 | write(*,*) "Unexpected value for LTcount: ",LTcount |
---|
| 1477 | stop |
---|
| 1478 | endif |
---|
| 1479 | |
---|
| 1480 | ! Generate the sol list for the interpolation |
---|
| 1481 | allocate(sol_list(solcount*LTcount)) |
---|
| 1482 | call gen_sol_list(solcount,int_sol_list,LTcount,LT_val,OBSLTmax(i,j,l),& |
---|
| 1483 | OBSLTmin(i,j,l),lon_val,LTmod,dayornight,& |
---|
| 1484 | sol_list) |
---|
| 1485 | |
---|
| 1486 | solerrcount=0 |
---|
| 1487 | solbinned_value=0 |
---|
| 1488 | sol_bin: do m=1,solcount*LTcount ! loop on all GCM sols of the bin |
---|
| 1489 | sol=sol_list(m) |
---|
| 1490 | ! write(*,*)"sol=",sol |
---|
| 1491 | !=============================================================================== |
---|
| 1492 | ! 3.3 Do the interpolation and binning for the given location |
---|
| 1493 | !=============================================================================== |
---|
| 1494 | call extraction(lon_val,lat_val,alt_val,sol,& |
---|
| 1495 | GCMlonlen,GCMlatlen,GCMaltlen,GCMtimelen,& |
---|
| 1496 | GCMlon,GCMlat,GCMalt,GCMtime,& |
---|
| 1497 | GCM_var,GCMmiss_val,GCMalttype,GCMvarname,extr_value) |
---|
| 1498 | |
---|
| 1499 | if (extr_value.eq.GCMmiss_val) then |
---|
| 1500 | ! write(*,*) "Missing value in gcmfile at lon=",lon_val,"; lat=",lat_val,"; alt=",alt_val,"; sol=",sol |
---|
| 1501 | solerrcount=solerrcount+1 |
---|
| 1502 | CYCLE sol_bin ! go directly to the next GCM sol of the bin |
---|
| 1503 | endif |
---|
| 1504 | solbinned_value=solbinned_value+extr_value |
---|
| 1505 | |
---|
| 1506 | enddo sol_bin ! end loop on all GCM sols of the bin |
---|
| 1507 | if ((solcount*LTcount-solerrcount).ne.0) then |
---|
| 1508 | solbinned_value=solbinned_value/(solcount*LTcount-solerrcount) |
---|
| 1509 | else |
---|
| 1510 | ! write(*,*)"No GCM value in this sol bin" |
---|
| 1511 | solbinned_value=OBSmiss_val |
---|
| 1512 | endif |
---|
| 1513 | ! Write value to output |
---|
| 1514 | outvar(i,j,k,l)=solbinned_value |
---|
| 1515 | |
---|
| 1516 | errcount=errcount+solerrcount |
---|
| 1517 | |
---|
| 1518 | deallocate(sol_list) |
---|
| 1519 | |
---|
| 1520 | enddo altitude ! end loop on observed altitudes |
---|
| 1521 | enddo longitude ! end loop on observed longitudes |
---|
| 1522 | enddo latitude ! end loop on observed latitudes |
---|
| 1523 | |
---|
| 1524 | deallocate(int_sol_list) |
---|
| 1525 | |
---|
| 1526 | enddo Ls ! end loop on observed Solar longitudes |
---|
| 1527 | |
---|
| 1528 | ! write(*,*)"Nb of GCM missing values :",errcount |
---|
| 1529 | |
---|
| 1530 | !=============================================================================== |
---|
| 1531 | ! 3.4 Write the data in the netcdf output file |
---|
| 1532 | !=============================================================================== |
---|
| 1533 | status = nf90_put_var(outfid, outvarid, outvar) |
---|
| 1534 | error_text="Error: could not write "//trim(outvarname)//" data in the outfile" |
---|
| 1535 | call status_check(status,error_text) |
---|
| 1536 | |
---|
| 1537 | ! Deallocations before going to the next GCM variable |
---|
| 1538 | deallocate(GCM_var) |
---|
| 1539 | deallocate(OBS_var) |
---|
| 1540 | deallocate(numbin) |
---|
| 1541 | deallocate(outvar) |
---|
| 1542 | ENDDO VAR ! end loop on variables |
---|
[2558] | 1543 | |
---|
| 1544 | |
---|
| 1545 | !=================================================================================== |
---|
| 1546 | ! 4. CDOD RATIO (tau_GCM/tau_MCS) |
---|
| 1547 | !=================================================================================== |
---|
| 1548 | !=============================================================================== |
---|
| 1549 | ! 4.1 Preliminary checks |
---|
| 1550 | !=============================================================================== |
---|
| 1551 | if (OBSalttype.eq.'p') then |
---|
| 1552 | tau_ratio_ok = .true. |
---|
| 1553 | else |
---|
| 1554 | tau_ratio_ok = .false. |
---|
| 1555 | endif |
---|
[2404] | 1556 | |
---|
[2558] | 1557 | IF (dayornight.EQ."dayside") THEN |
---|
| 1558 | out_dztauname = "d" |
---|
| 1559 | OBS_dztauname = "ddust" |
---|
| 1560 | OBS_tempname = "dtemp" |
---|
| 1561 | ELSE ! i.e. dayornight="nightside" |
---|
| 1562 | out_dztauname = "n" |
---|
| 1563 | OBS_dztauname = "ndust" |
---|
| 1564 | OBS_tempname = "ntemp" |
---|
| 1565 | ENDIF |
---|
| 1566 | |
---|
| 1567 | ! Check that the output file contains dust opacity |
---|
| 1568 | ! NB: a dust opacity in outfile requires the dust |
---|
| 1569 | ! opacity variable to be present in the obsfile, |
---|
| 1570 | ! so we don't have to check it out |
---|
| 1571 | if (tau_ratio_ok) then |
---|
| 1572 | out_dztauname = out_dztauname(1:1)//"dust" |
---|
| 1573 | status=nf90_inq_varid(outfid,trim(out_dztauname),out_dztauid) |
---|
| 1574 | if (status.ne.nf90_noerr) then |
---|
[2979] | 1575 | ! if no "d/ndust" in outfile, look for "d/nopa_dust" |
---|
| 1576 | out_dztauname = out_dztauname(1:1)//"opa_dust" |
---|
[2558] | 1577 | status=nf90_inq_varid(outfid,trim(out_dztauname),out_dztauid) |
---|
| 1578 | if (status.ne.nf90_noerr) then |
---|
[2979] | 1579 | ! if no "d/ndust" nor d/nopa_dust" in outfile, look for "d/nopadust" (old variable name from aeroptical, before r2817) |
---|
| 1580 | out_dztauname = out_dztauname(1:1)//"opadust" |
---|
| 1581 | status=nf90_inq_varid(outfid,trim(out_dztauname),out_dztauid) |
---|
| 1582 | if (status.ne.nf90_noerr) then |
---|
| 1583 | tau_ratio_ok = .false. |
---|
| 1584 | write(*,*)"No dust opacity in the outfile, we skip the computation of the CDOD ratio" |
---|
| 1585 | endif |
---|
[2558] | 1586 | endif |
---|
| 1587 | endif |
---|
| 1588 | endif |
---|
| 1589 | |
---|
| 1590 | ! Check that the obsfile contains temperature (for rho computation) |
---|
| 1591 | if (tau_ratio_ok) then |
---|
| 1592 | status=nf90_inq_varid(obsfid,trim(OBS_tempname),OBS_tempid) |
---|
| 1593 | if (status.ne.nf90_noerr) then |
---|
| 1594 | tau_ratio_ok = .false. |
---|
[2979] | 1595 | write(*,*)"No temperature in the obsfile, we skip the computation of the CDOD ratio" |
---|
[2558] | 1596 | endif |
---|
| 1597 | endif |
---|
| 1598 | |
---|
| 1599 | |
---|
| 1600 | !=============================================================================== |
---|
| 1601 | ! 4.2 tau_ratio computation |
---|
| 1602 | !=============================================================================== |
---|
| 1603 | if (tau_ratio_ok) then |
---|
| 1604 | write(*,*)"Computing the Column Dust Optical Depths ratio (tau_ratio = tau_GCM/tau_MCS)..." |
---|
| 1605 | write(*,*)"(it is recommended to use it to normalize the GCM dust profiles when comparing with MCS dust profiles)" |
---|
| 1606 | |
---|
| 1607 | ! Creation of the OBS_plev table |
---|
| 1608 | IF (dayornight.EQ."dayside") THEN ! only do it once |
---|
| 1609 | allocate(OBS_plev(OBSaltlen+1)) |
---|
| 1610 | DO k=2,OBSaltlen |
---|
| 1611 | ! mid-altitude pressure between 2 OBS pressure layers (OBSalt) |
---|
| 1612 | ! NB: OBS_plev(k) > OBSalt(k) > OBS_plev(k+1) |
---|
| 1613 | OBS_plev(k) = sqrt(OBSalt(k-1)*OBSalt(k)) |
---|
| 1614 | ENDDO |
---|
| 1615 | OBS_plev(1)=2000. |
---|
| 1616 | OBS_plev(OBSaltlen+1)=0 |
---|
| 1617 | write(*,*)"OBS_plev = ",OBS_plev |
---|
| 1618 | ENDIF |
---|
| 1619 | |
---|
| 1620 | ! Arrays allocation |
---|
| 1621 | allocate(out_dztau(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1622 | allocate(OBS_dztau(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1623 | allocate(OBS_temp(OBSlonlen,OBSlatlen,OBSaltlen,OBSLslen)) |
---|
| 1624 | |
---|
| 1625 | allocate(tau_ratio(OBSlonlen,OBSlatlen,OBSLslen)) |
---|
| 1626 | tau_ratio(:,:,:) = 0 |
---|
| 1627 | |
---|
| 1628 | ! Load GCM opacity |
---|
| 1629 | status=NF90_GET_VAR(outfid,out_dztauid,out_dztau) |
---|
| 1630 | error_text="Failed to load outfile dust opacity" |
---|
| 1631 | call status_check(status,error_text) |
---|
| 1632 | |
---|
| 1633 | ! Load OBS opacity |
---|
| 1634 | status=nf90_inq_varid(obsfid,trim(OBS_dztauname),OBS_dztauid) |
---|
| 1635 | status=NF90_GET_VAR(OBSfid,OBS_dztauid,OBS_dztau) |
---|
| 1636 | error_text="Failed to load obsfile dust opacity" |
---|
| 1637 | call status_check(status,error_text) |
---|
| 1638 | |
---|
| 1639 | ! Load OBS temperature |
---|
| 1640 | status=NF90_GET_VAR(OBSfid,OBS_tempid,OBS_temp) |
---|
| 1641 | error_text="Failed to load obsfile temperature" |
---|
| 1642 | call status_check(status,error_text) |
---|
| 1643 | |
---|
| 1644 | do l=1,OBSLslen |
---|
| 1645 | do j=1,OBSlatlen |
---|
| 1646 | do i=1,OBSlonlen |
---|
| 1647 | ! initialize both columns to 0 |
---|
| 1648 | out_tau = 0 |
---|
| 1649 | OBS_tau = 0 |
---|
| 1650 | do k=1,OBSlatlen |
---|
| 1651 | ! Compute optical depths where both GCM and OBS have non-missing values |
---|
| 1652 | ! (especially above the highest (in altitude) of the lowest (in altitude) valid values in the GCM and Observer profile) |
---|
[3119] | 1653 | IF ((out_dztau(i,j,k,l).ne.OBSmiss_val).and.& |
---|
| 1654 | (OBS_dztau(i,j,k,l).ne.OBSmiss_val).and.& |
---|
| 1655 | (OBS_temp(i,j,k,l).ne.OBSmiss_val)) then |
---|
[2558] | 1656 | OBS_rho = OBSalt(k) / (r_atm*OBS_temp(i,j,k,l)) |
---|
| 1657 | |
---|
| 1658 | OBS_tau = OBS_tau + OBS_dztau(i,j,k,l)*1.e-3 / OBS_rho /g *(OBS_plev(k)-OBS_plev(k+1)) |
---|
| 1659 | |
---|
| 1660 | out_tau = out_tau + out_dztau(i,j,k,l)*1.e-3 / OBS_rho /g *(OBS_plev(k)-OBS_plev(k+1)) |
---|
| 1661 | |
---|
| 1662 | ENDIF |
---|
| 1663 | enddo !alt |
---|
| 1664 | |
---|
| 1665 | ! Compute the Column Dust Optical Depth ratio |
---|
| 1666 | if ((out_tau.eq.0).or.(OBS_tau.eq.0)) then |
---|
| 1667 | tau_ratio(i,j,l) = LTmiss_val |
---|
| 1668 | else |
---|
| 1669 | tau_ratio(i,j,l) = out_tau/OBS_tau |
---|
| 1670 | endif |
---|
| 1671 | |
---|
| 1672 | enddo !lon |
---|
| 1673 | enddo !lat |
---|
| 1674 | enddo !Ls |
---|
| 1675 | |
---|
| 1676 | |
---|
| 1677 | !=============================================================================== |
---|
| 1678 | ! 4.3 Write tau_ratio in the netcdf output file |
---|
| 1679 | !=============================================================================== |
---|
| 1680 | ! Switch to netcdf define mode |
---|
| 1681 | status=nf90_redef(outfid) |
---|
| 1682 | error_text="Error: could not switch to define mode in the outfile" |
---|
| 1683 | call status_check(status,error_text) |
---|
| 1684 | |
---|
| 1685 | ! Definition of the variable |
---|
| 1686 | SELECT CASE (dayornight) |
---|
| 1687 | CASE ("dayside") |
---|
| 1688 | outvarname = "dtau_ratio" |
---|
| 1689 | CASE ("nightside") |
---|
| 1690 | outvarname = "ntau_ratio" |
---|
| 1691 | END SELECT |
---|
| 1692 | status=NF90_DEF_VAR(outfid,trim(outvarname),nf90_float,LTshape,outvarid) |
---|
| 1693 | error_text="Error: could not define the variable "//trim(outvarname)//" in the outfile" |
---|
| 1694 | call status_check(status,error_text) |
---|
| 1695 | |
---|
| 1696 | ! Write the attributes |
---|
| 1697 | long_name = "Integrated optical depths ratio tau_GCM/tau_MCS" |
---|
| 1698 | SELECT CASE (dayornight) |
---|
| 1699 | CASE ("dayside") |
---|
| 1700 | long_name = trim(long_name)//" - day side" |
---|
| 1701 | status=nf90_put_att(outfid,outvarid,"long_name",long_name) |
---|
| 1702 | CASE ("nightside") |
---|
| 1703 | long_name = trim(long_name)//" - night side" |
---|
| 1704 | status=nf90_put_att(outfid,outvarid,"long_name",long_name) |
---|
| 1705 | END SELECT |
---|
| 1706 | status=nf90_put_att(outfid,outvarid,"units","/") |
---|
| 1707 | status=nf90_put_att(outfid,outvarid,"_FillValue",LTmiss_val) |
---|
| 1708 | comment = "Reference numbin: dust+temp" |
---|
| 1709 | status=nf90_put_att(outfid,outvarid,"comment",comment) |
---|
| 1710 | |
---|
| 1711 | write(*,*)trim(outvarname)," has been created in the outfile" |
---|
| 1712 | write(*,'(" with missing_value attribute : ",1pe12.5)')LTmiss_val |
---|
| 1713 | write(*,*)"" |
---|
| 1714 | |
---|
| 1715 | ! End the netcdf define mode (and thus enter the "data writing" mode) |
---|
| 1716 | status=nf90_enddef(outfid) |
---|
| 1717 | error_text="Error: could not close the define mode of the outfile" |
---|
| 1718 | call status_check(status,error_text) |
---|
| 1719 | |
---|
| 1720 | ! Write data |
---|
| 1721 | status = nf90_put_var(outfid, outvarid, tau_ratio) |
---|
| 1722 | error_text="Error: could not write "//trim(outvarname)//" data in the outfile" |
---|
| 1723 | call status_check(status,error_text) |
---|
| 1724 | |
---|
| 1725 | ! Deallocations |
---|
| 1726 | deallocate(out_dztau) |
---|
| 1727 | deallocate(OBS_dztau) |
---|
| 1728 | deallocate(OBS_temp) |
---|
| 1729 | deallocate(tau_ratio) |
---|
| 1730 | |
---|
| 1731 | endif !tau_ratio_ok |
---|
| 1732 | |
---|
| 1733 | |
---|
[2404] | 1734 | !=================================================================================== |
---|
[2558] | 1735 | ! 5. END OF THE DAY/NIGHT LOOP |
---|
[2404] | 1736 | !=================================================================================== |
---|
| 1737 | IF (dayornight.EQ."dayside") THEN |
---|
| 1738 | ! this is the end of the first loop (on daytime values) |
---|
| 1739 | ! and we still have to loop on nighttime values |
---|
| 1740 | dayornight="nightside" |
---|
| 1741 | deallocate(OBSLT) |
---|
| 1742 | deallocate(OBSLTmax) |
---|
| 1743 | deallocate(OBSLTmin) |
---|
| 1744 | write(*,*)"" |
---|
| 1745 | write(*,*) "Beginning the 2nd loop, on nighttime values"; write(*,*)"" |
---|
| 1746 | CYCLE DAY_OR_NIGHT |
---|
| 1747 | ELSE ! i.e. dayornight="nightside" |
---|
| 1748 | ! this is the end of the second loop (on nighttime values) |
---|
| 1749 | ! and thus the end of the program |
---|
| 1750 | EXIT DAY_OR_NIGHT |
---|
| 1751 | ENDIF |
---|
| 1752 | ENDDO DAY_OR_NIGHT ! end of the day/night loop that begins in section 1.3 |
---|
| 1753 | |
---|
| 1754 | !=================================================================================== |
---|
[2558] | 1755 | ! 6. CLOSE THE FILES |
---|
[2404] | 1756 | !=================================================================================== |
---|
| 1757 | status = nf90_close(gcmfid) |
---|
| 1758 | error_text="Error: could not close file "//trim(gcmfile) |
---|
| 1759 | call status_check(status,error_text) |
---|
| 1760 | status = nf90_close(obsfid) |
---|
| 1761 | error_text="Error: could not close file "//trim(obsfile) |
---|
| 1762 | call status_check(status,error_text) |
---|
| 1763 | status = nf90_close(outfid) |
---|
| 1764 | error_text="Error: could not close file "//trim(outfile) |
---|
| 1765 | call status_check(status,error_text) |
---|
| 1766 | |
---|
| 1767 | write(*,*)"";write(*,*)"-> Program simu_MCS completed!" |
---|
| 1768 | |
---|
| 1769 | end program simu_MCS |
---|
| 1770 | |
---|
| 1771 | |
---|
| 1772 | |
---|
| 1773 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 1774 | !=================================================================================================== |
---|
| 1775 | ! Subroutines |
---|
| 1776 | !=================================================================================================== |
---|
| 1777 | |
---|
| 1778 | subroutine extraction(lon,lat,alt,sol,& |
---|
| 1779 | lonlen,latlen,altlen,timelen,& |
---|
| 1780 | longitude,latitude,altitude,time,& |
---|
| 1781 | field,missing_value,alttype,varname,value) |
---|
| 1782 | |
---|
| 1783 | implicit none |
---|
| 1784 | !================================================================ |
---|
| 1785 | ! Arguments: |
---|
| 1786 | !================================================================ |
---|
| 1787 | real,intent(in) :: lon ! sought longitude (deg, in [-180:180]) |
---|
| 1788 | real,intent(in) :: lat ! sought latitude (deg, in [-90:90]) |
---|
| 1789 | real,intent(in) :: alt ! sought altitude (m or Pa) |
---|
| 1790 | real,intent(in) :: sol ! sought date (sols) |
---|
| 1791 | integer,intent(in) :: lonlen |
---|
| 1792 | integer,intent(in) :: latlen |
---|
| 1793 | integer,intent(in) :: altlen |
---|
| 1794 | integer,intent(in) :: timelen |
---|
| 1795 | real,intent(in) :: longitude(lonlen) |
---|
| 1796 | real,intent(in) :: latitude(latlen) |
---|
| 1797 | real,intent(in) :: altitude(altlen) |
---|
| 1798 | real,intent(in) :: time(timelen) |
---|
| 1799 | real,intent(in) :: field(lonlen,latlen,altlen,timelen) |
---|
| 1800 | real,intent(in) :: missing_value ! default value for "no data" |
---|
| 1801 | character,intent(in) :: alttype ! altitude coord. type:'z' (altitude, m) |
---|
| 1802 | ! 'p' (pressure, Pa) |
---|
| 1803 | character(len=*),intent(in) :: varname ! variable name (in GCM file) |
---|
| 1804 | real,intent(out) :: value |
---|
| 1805 | |
---|
| 1806 | !================================================================ |
---|
| 1807 | ! Local variables: |
---|
| 1808 | !================================================================ |
---|
| 1809 | real,save :: prev_lon=-99 ! previous value of 'lon' routine was called with |
---|
| 1810 | real,save :: prev_lat=-99 ! previous value of 'lat' routine was called with |
---|
| 1811 | real,save :: prev_alt=-9.e20 ! ! previous value of 'alt' |
---|
| 1812 | real,save :: prev_sol=-99 ! previous value of 'sol' routine was called with |
---|
| 1813 | |
---|
| 1814 | ! encompasing indexes: |
---|
| 1815 | integer,save :: ilon_inf=-1,ilon_sup=-1 ! along longitude |
---|
| 1816 | integer,save :: ilat_inf=-1,ilat_sup=-1 ! along latitude |
---|
| 1817 | integer,save :: ialt_inf=-1,ialt_sup=-1 ! along altitude |
---|
| 1818 | integer,save :: itim_inf=-1,itim_sup=-1 ! along time |
---|
| 1819 | |
---|
| 1820 | ! intermediate interpolated values |
---|
| 1821 | real :: t_interp(2,2,2) ! after time interpolation |
---|
| 1822 | real :: zt_interp(2,2) ! after altitude interpolation |
---|
| 1823 | real :: yzt_interp(2) ! after latitude interpolation |
---|
| 1824 | real :: coeff ! interpolation coefficient |
---|
| 1825 | |
---|
| 1826 | integer :: i |
---|
| 1827 | |
---|
| 1828 | ! By default, set value to missing_value |
---|
| 1829 | value=missing_value |
---|
| 1830 | |
---|
| 1831 | !================================================================ |
---|
| 1832 | ! 1. Find encompassing indexes |
---|
| 1833 | !================================================================ |
---|
| 1834 | if (lon.ne.prev_lon) then |
---|
| 1835 | do i=1,lonlen-1 |
---|
| 1836 | if (longitude(i).le.lon) then |
---|
| 1837 | ilon_inf=i |
---|
| 1838 | else |
---|
| 1839 | exit |
---|
| 1840 | endif |
---|
| 1841 | enddo |
---|
| 1842 | ilon_sup=ilon_inf+1 |
---|
| 1843 | endif ! of if (lon.ne.prev_lon) |
---|
| 1844 | !write(*,*) 'ilon_inf=',ilon_inf," longitude(ilon_inf)=",longitude(ilon_inf) |
---|
| 1845 | |
---|
| 1846 | if (lat.ne.prev_lat) then |
---|
[2540] | 1847 | if (latitude(1).gt.latitude(2)) then |
---|
| 1848 | ! decreasing latitudes, from 90N to (-)90S (LMDZ-like) |
---|
| 1849 | do i=1,latlen-1 |
---|
| 1850 | if (latitude(i).ge.lat) then |
---|
| 1851 | ilat_inf=i |
---|
| 1852 | else |
---|
| 1853 | exit |
---|
| 1854 | endif |
---|
| 1855 | enddo |
---|
| 1856 | else |
---|
| 1857 | ! increasing latitudes, from (-)90S to 90N (DYNAMICO-like) |
---|
| 1858 | do i=1,latlen-1 |
---|
| 1859 | if (latitude(i).le.lat) then |
---|
| 1860 | ilat_inf=i |
---|
| 1861 | else |
---|
| 1862 | exit |
---|
| 1863 | endif |
---|
| 1864 | enddo |
---|
| 1865 | endif |
---|
[2404] | 1866 | ilat_sup=ilat_inf+1 |
---|
| 1867 | endif ! of if (lat.ne.prev_lat) |
---|
| 1868 | !write(*,*) 'ilat_inf=',ilat_inf," latitude(ilat_inf)=",latitude(ilat_inf) |
---|
| 1869 | |
---|
| 1870 | if (alt.ne.prev_alt) then |
---|
| 1871 | if (alttype.eq.'p') then ! pressures are ordered from max to min |
---|
| 1872 | !handle special case for alt not in the altitude(1:altlen) interval |
---|
| 1873 | if ((alt.gt.altitude(1)).or.(alt.lt.altitude(altlen))) then |
---|
| 1874 | ialt_inf=-1 |
---|
| 1875 | ialt_sup=-1 |
---|
| 1876 | ! return to main program (with value=missing_value) |
---|
| 1877 | ! write(*,*)"Problem in extraction : GCM alt p" |
---|
| 1878 | return |
---|
| 1879 | else ! general case |
---|
| 1880 | do i=1,altlen-1 |
---|
| 1881 | if (altitude(i).ge.alt) then |
---|
| 1882 | ialt_inf=i |
---|
| 1883 | else |
---|
| 1884 | exit |
---|
| 1885 | endif |
---|
| 1886 | enddo |
---|
| 1887 | ialt_sup=ialt_inf+1 |
---|
| 1888 | endif ! of if ((alt.gt.altitude(1)).or.(alt.lt.altitude(altlen))) |
---|
| 1889 | else ! altitudes (m) are ordered from min to max |
---|
| 1890 | !handle special case for alt not in the altitude(1:altlen) interval |
---|
| 1891 | if ((alt.lt.altitude(1)).or.(alt.gt.altitude(altlen))) then |
---|
| 1892 | ialt_inf=-1 |
---|
| 1893 | ialt_sup=-1 |
---|
| 1894 | ! return to main program (with value=missing_value) |
---|
| 1895 | ! write(*,*)"Problem in extraction : GCM alt z" |
---|
| 1896 | return |
---|
| 1897 | else ! general case |
---|
| 1898 | do i=1,altlen-1 |
---|
| 1899 | if (altitude(i).le.alt) then |
---|
| 1900 | ialt_inf=i |
---|
| 1901 | else |
---|
| 1902 | exit |
---|
| 1903 | endif |
---|
| 1904 | enddo |
---|
| 1905 | ialt_sup=ialt_inf+1 |
---|
| 1906 | endif ! of if ((alt.lt.altitude(1)).or.(alt.gt.altitude(altlen))) |
---|
| 1907 | endif ! of if (alttype.eq.'p') |
---|
| 1908 | endif ! of if (alt.ne.prev_alt) |
---|
| 1909 | !write(*,*) 'ialt_inf=',ialt_inf," altitude(ialt_inf)=",altitude(ialt_inf) |
---|
| 1910 | |
---|
| 1911 | if (sol.ne.prev_sol) then |
---|
| 1912 | !handle special case for sol not in the time(1):time(timenlen) interval |
---|
| 1913 | if ((sol.lt.time(1)).or.(sol.gt.time(timelen))) then |
---|
| 1914 | itim_inf=-1 |
---|
| 1915 | itim_sup=-1 |
---|
| 1916 | ! return to main program (with value=missing_value) |
---|
| 1917 | ! write(*,*)"Problem in extraction : GCM sol" |
---|
| 1918 | return |
---|
| 1919 | else ! general case |
---|
| 1920 | do i=1,timelen-1 |
---|
| 1921 | if (time(i).le.sol) then |
---|
| 1922 | itim_inf=i |
---|
| 1923 | else |
---|
| 1924 | exit |
---|
| 1925 | endif |
---|
| 1926 | enddo |
---|
| 1927 | itim_sup=itim_inf+1 |
---|
| 1928 | endif ! of if ((sol.lt.time(1)).or.(sol.gt.time(timenlen))) |
---|
| 1929 | endif ! of if (sol.ne.prev_sol) |
---|
| 1930 | !write(*,*) 'itim_inf=',itim_inf," time(itim_inf)=",time(itim_inf) |
---|
| 1931 | !write(*,*) 'itim_sup=',itim_sup," time(itim_sup)=",time(itim_sup) |
---|
| 1932 | |
---|
| 1933 | !================================================================ |
---|
| 1934 | ! 2. Interpolate |
---|
| 1935 | !================================================================ |
---|
| 1936 | ! check that there are no "missing_value" in the field() elements we need |
---|
| 1937 | ! otherwise return to main program (with value=missing_value) |
---|
| 1938 | if (field(ilon_inf,ilat_inf,ialt_inf,itim_inf).eq.missing_value) then |
---|
| 1939 | ! write(*,*)"Error 1 in extraction" |
---|
| 1940 | return |
---|
| 1941 | endif |
---|
| 1942 | if (field(ilon_inf,ilat_inf,ialt_inf,itim_sup).eq.missing_value) then |
---|
| 1943 | ! write(*,*)"Error 2 in extraction" |
---|
| 1944 | return |
---|
| 1945 | endif |
---|
| 1946 | if (field(ilon_inf,ilat_inf,ialt_sup,itim_inf).eq.missing_value) then |
---|
| 1947 | ! write(*,*)"Error 3 in extraction" |
---|
| 1948 | return |
---|
| 1949 | endif |
---|
| 1950 | if (field(ilon_inf,ilat_inf,ialt_sup,itim_sup).eq.missing_value) then |
---|
| 1951 | ! write(*,*)"Error 4 in extraction" |
---|
| 1952 | return |
---|
| 1953 | endif |
---|
| 1954 | if (field(ilon_inf,ilat_sup,ialt_inf,itim_inf).eq.missing_value) then |
---|
| 1955 | ! write(*,*)"Error 5 in extraction" |
---|
| 1956 | return |
---|
| 1957 | endif |
---|
| 1958 | if (field(ilon_inf,ilat_sup,ialt_inf,itim_sup).eq.missing_value) then |
---|
| 1959 | ! write(*,*)"Error 6 in extraction" |
---|
| 1960 | return |
---|
| 1961 | endif |
---|
| 1962 | if (field(ilon_inf,ilat_sup,ialt_sup,itim_inf).eq.missing_value) then |
---|
| 1963 | ! write(*,*)"Error 7 in extraction" |
---|
| 1964 | return |
---|
| 1965 | endif |
---|
| 1966 | if (field(ilon_inf,ilat_sup,ialt_sup,itim_sup).eq.missing_value) then |
---|
| 1967 | ! write(*,*)"Error 8 in extraction" |
---|
| 1968 | return |
---|
| 1969 | endif |
---|
| 1970 | if (field(ilon_sup,ilat_inf,ialt_inf,itim_inf).eq.missing_value) then |
---|
| 1971 | ! write(*,*)"Error 9 in extraction" |
---|
| 1972 | return |
---|
| 1973 | endif |
---|
| 1974 | if (field(ilon_sup,ilat_inf,ialt_inf,itim_sup).eq.missing_value) then |
---|
| 1975 | ! write(*,*)"Error 10 in extraction" |
---|
| 1976 | return |
---|
| 1977 | endif |
---|
| 1978 | if (field(ilon_sup,ilat_inf,ialt_sup,itim_inf).eq.missing_value) then |
---|
| 1979 | ! write(*,*)"Error 11 in extraction" |
---|
| 1980 | return |
---|
| 1981 | endif |
---|
| 1982 | if (field(ilon_sup,ilat_inf,ialt_sup,itim_sup).eq.missing_value) then |
---|
| 1983 | ! write(*,*)"Error 12 in extraction" |
---|
| 1984 | return |
---|
| 1985 | endif |
---|
| 1986 | if (field(ilon_sup,ilat_sup,ialt_inf,itim_inf).eq.missing_value) then |
---|
| 1987 | ! write(*,*)"Error 13 in extraction" |
---|
| 1988 | return |
---|
| 1989 | endif |
---|
| 1990 | if (field(ilon_sup,ilat_sup,ialt_inf,itim_sup).eq.missing_value) then |
---|
| 1991 | ! write(*,*)"Error 14 in extraction" |
---|
| 1992 | return |
---|
| 1993 | endif |
---|
| 1994 | if (field(ilon_sup,ilat_sup,ialt_sup,itim_inf).eq.missing_value) then |
---|
| 1995 | ! write(*,*)"Error 15 in extraction" |
---|
| 1996 | return |
---|
| 1997 | endif |
---|
| 1998 | if (field(ilon_sup,ilat_sup,ialt_sup,itim_sup).eq.missing_value) then |
---|
| 1999 | ! write(*,*)"Error 16 in extraction" |
---|
| 2000 | return |
---|
| 2001 | endif |
---|
| 2002 | |
---|
| 2003 | !================================================================ |
---|
| 2004 | ! 2.1 Linear interpolation in time |
---|
| 2005 | !================================================================ |
---|
| 2006 | coeff=(sol-time(itim_inf))/(time(itim_sup)-time(itim_inf)) |
---|
| 2007 | t_interp(1,1,1)=field(ilon_inf,ilat_inf,ialt_inf,itim_inf)+ & |
---|
| 2008 | coeff*(field(ilon_inf,ilat_inf,ialt_inf,itim_sup)- & |
---|
| 2009 | field(ilon_inf,ilat_inf,ialt_inf,itim_inf)) |
---|
| 2010 | t_interp(1,1,2)=field(ilon_inf,ilat_inf,ialt_sup,itim_inf)+ & |
---|
| 2011 | coeff*(field(ilon_inf,ilat_inf,ialt_sup,itim_sup)- & |
---|
| 2012 | field(ilon_inf,ilat_inf,ialt_sup,itim_inf)) |
---|
| 2013 | t_interp(1,2,1)=field(ilon_inf,ilat_sup,ialt_inf,itim_inf)+ & |
---|
| 2014 | coeff*(field(ilon_inf,ilat_sup,ialt_inf,itim_sup)- & |
---|
| 2015 | field(ilon_inf,ilat_sup,ialt_inf,itim_inf)) |
---|
| 2016 | t_interp(1,2,2)=field(ilon_inf,ilat_sup,ialt_sup,itim_inf)+ & |
---|
| 2017 | coeff*(field(ilon_inf,ilat_sup,ialt_sup,itim_sup)- & |
---|
| 2018 | field(ilon_inf,ilat_sup,ialt_sup,itim_inf)) |
---|
| 2019 | t_interp(2,1,1)=field(ilon_sup,ilat_inf,ialt_inf,itim_inf)+ & |
---|
| 2020 | coeff*(field(ilon_sup,ilat_inf,ialt_inf,itim_sup)- & |
---|
| 2021 | field(ilon_sup,ilat_inf,ialt_inf,itim_inf)) |
---|
| 2022 | t_interp(2,1,2)=field(ilon_sup,ilat_inf,ialt_sup,itim_inf)+ & |
---|
| 2023 | coeff*(field(ilon_sup,ilat_inf,ialt_sup,itim_sup)- & |
---|
| 2024 | field(ilon_sup,ilat_inf,ialt_sup,itim_inf)) |
---|
| 2025 | t_interp(2,2,1)=field(ilon_sup,ilat_sup,ialt_inf,itim_inf)+ & |
---|
| 2026 | coeff*(field(ilon_sup,ilat_sup,ialt_inf,itim_sup)- & |
---|
| 2027 | field(ilon_sup,ilat_sup,ialt_inf,itim_inf)) |
---|
| 2028 | t_interp(2,2,2)=field(ilon_sup,ilat_sup,ialt_sup,itim_inf)+ & |
---|
| 2029 | coeff*(field(ilon_sup,ilat_sup,ialt_sup,itim_sup)- & |
---|
| 2030 | field(ilon_sup,ilat_sup,ialt_sup,itim_inf)) |
---|
| 2031 | |
---|
| 2032 | !================================================================ |
---|
| 2033 | ! 2.2 Vertical interpolation |
---|
| 2034 | !================================================================ |
---|
| 2035 | if (((varname=='rho').or.(varname=='pressure')).and.(alttype=='z')) then |
---|
| 2036 | ! do the interpolation on the log of the quantity |
---|
| 2037 | coeff=(alt-altitude(ialt_inf))/(altitude(ialt_sup)-altitude(ialt_inf)) |
---|
| 2038 | zt_interp(1,1)=log(t_interp(1,1,1))+coeff* & |
---|
| 2039 | (log(t_interp(1,1,2))-log(t_interp(1,1,1))) |
---|
| 2040 | zt_interp(1,2)=log(t_interp(1,2,1))+coeff* & |
---|
| 2041 | (log(t_interp(1,2,2))-log(t_interp(1,2,1))) |
---|
| 2042 | zt_interp(2,1)=log(t_interp(2,1,1))+coeff* & |
---|
| 2043 | (log(t_interp(2,1,2))-log(t_interp(2,1,1))) |
---|
| 2044 | zt_interp(2,2)=log(t_interp(2,2,1))+coeff* & |
---|
| 2045 | (log(t_interp(2,2,2))-log(t_interp(2,2,1))) |
---|
| 2046 | zt_interp(1:2,1:2)=exp(zt_interp(1:2,1:2)) |
---|
| 2047 | else ! general case |
---|
| 2048 | coeff=(alt-altitude(ialt_inf))/(altitude(ialt_sup)-altitude(ialt_inf)) |
---|
| 2049 | zt_interp(1,1)=t_interp(1,1,1)+coeff*(t_interp(1,1,2)-t_interp(1,1,1)) |
---|
| 2050 | zt_interp(1,2)=t_interp(1,2,1)+coeff*(t_interp(1,2,2)-t_interp(1,2,1)) |
---|
| 2051 | zt_interp(2,1)=t_interp(2,1,1)+coeff*(t_interp(2,1,2)-t_interp(2,1,1)) |
---|
| 2052 | zt_interp(2,2)=t_interp(2,2,1)+coeff*(t_interp(2,2,2)-t_interp(2,2,1)) |
---|
| 2053 | endif |
---|
| 2054 | |
---|
| 2055 | !================================================================ |
---|
| 2056 | ! 2.3 Latitudinal interpolation |
---|
| 2057 | !================================================================ |
---|
| 2058 | coeff=(lat-latitude(ilat_inf))/(latitude(ilat_sup)-latitude(ilat_inf)) |
---|
| 2059 | yzt_interp(1)=zt_interp(1,1)+coeff*(zt_interp(1,2)-zt_interp(1,1)) |
---|
| 2060 | yzt_interp(2)=zt_interp(2,1)+coeff*(zt_interp(2,2)-zt_interp(2,1)) |
---|
| 2061 | |
---|
| 2062 | !================================================================ |
---|
| 2063 | ! 2.4 longitudinal interpolation |
---|
| 2064 | !================================================================ |
---|
| 2065 | coeff=(lon-longitude(ilon_inf))/(longitude(ilon_sup)-longitude(ilon_inf)) |
---|
| 2066 | value=yzt_interp(1)+coeff*(yzt_interp(2)-yzt_interp(1)) |
---|
| 2067 | |
---|
| 2068 | end subroutine extraction |
---|
| 2069 | |
---|
| 2070 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2071 | |
---|
| 2072 | subroutine inidim(outfid,lonlen,latlen,altlen,timelen,lon,lat,alt,time,units_alt,& |
---|
| 2073 | londimid,latdimid,altdimid,timedimid) |
---|
| 2074 | !================================================================ |
---|
| 2075 | ! Purpose: |
---|
| 2076 | ! Initialize a data file (NetCDF format) |
---|
| 2077 | !================================================================ |
---|
| 2078 | ! Remarks: |
---|
| 2079 | ! The NetCDF file remains open |
---|
| 2080 | !================================================================ |
---|
| 2081 | use netcdf ! NetCDF definitions |
---|
| 2082 | implicit none |
---|
| 2083 | !================================================================ |
---|
| 2084 | ! Arguments: |
---|
| 2085 | !================================================================ |
---|
| 2086 | integer, intent(in):: outfid |
---|
| 2087 | ! outfid: [netcdf] file ID |
---|
| 2088 | integer, intent(in):: lonlen |
---|
| 2089 | ! lonlen: longitude length (# of longitude values) |
---|
| 2090 | integer, intent(in):: latlen |
---|
| 2091 | ! latlen: latitude length (# of latitude values) |
---|
| 2092 | integer, intent(in):: altlen |
---|
| 2093 | ! altlen: altitude length (# of altitude values) |
---|
| 2094 | integer, intent(in):: timelen |
---|
| 2095 | ! timelen: time length (# of time values) |
---|
| 2096 | real, intent(in):: lon(lonlen) |
---|
| 2097 | ! lon(): longitude |
---|
| 2098 | real, intent(in):: lat(latlen) |
---|
| 2099 | ! lat(): latitude |
---|
| 2100 | real, intent(in):: alt(altlen) |
---|
| 2101 | ! alt(): altitude |
---|
| 2102 | real, intent(in):: time(timelen) |
---|
| 2103 | ! time(): time (Ls) |
---|
| 2104 | character(len=1), intent(in) :: units_alt |
---|
| 2105 | ! units_alt: altitude coord. type:'z' (altitude, m) 'p' (pressure, Pa) |
---|
| 2106 | integer,intent(inout):: londimid |
---|
| 2107 | ! londimid: [netcdf] lon() (i.e.: longitude) ID in MCS and output files (they are the same) |
---|
| 2108 | integer,intent(inout) :: latdimid |
---|
| 2109 | ! latdimid: [netcdf] lat() ID in MCS and output files (they are the same) |
---|
| 2110 | integer,intent(inout):: altdimid |
---|
| 2111 | ! altdimid: [netcdf] alt() ID in MCS and output files (they are the same) |
---|
| 2112 | integer,intent(inout):: timedimid |
---|
| 2113 | ! timedimid: [netcdf] time() ID in MCS and output files (they are the same) |
---|
| 2114 | |
---|
| 2115 | !================================================================ |
---|
| 2116 | ! Local variables: |
---|
| 2117 | !================================================================ |
---|
| 2118 | integer :: lonvarid,latvarid,altvarid,timevarid,status |
---|
| 2119 | ! *varid: [netcdf] ID of a variable |
---|
| 2120 | ! status: [netcdf] return error code (from called subroutines) |
---|
| 2121 | character(len=256) :: error_text |
---|
| 2122 | integer :: d ! loop index on dimensions ID |
---|
| 2123 | |
---|
| 2124 | !=============================================================== |
---|
| 2125 | ! 1. Define/write "dimensions" in the same order than MCS file |
---|
| 2126 | ! and get their IDs |
---|
| 2127 | !================================================================ |
---|
| 2128 | do d=1,4 |
---|
| 2129 | if (altdimid.eq.d) then |
---|
| 2130 | status=nf90_def_dim(outfid, "altitude", altlen, altdimid) |
---|
| 2131 | error_text="Error: could not define the altitude dimension in the outfile" |
---|
| 2132 | call status_check(status,error_text) |
---|
| 2133 | else if (timedimid.eq.d) then |
---|
| 2134 | status=nf90_def_dim(outfid, "time", timelen, timedimid) |
---|
| 2135 | error_text="Error: could not define the time dimension in the outfile" |
---|
| 2136 | call status_check(status,error_text) |
---|
| 2137 | else if (latdimid.eq.d) then |
---|
| 2138 | status=nf90_def_dim(outfid, "latitude", latlen, latdimid) |
---|
| 2139 | error_text="Error: could not define the latitude dimension in the outfile" |
---|
| 2140 | call status_check(status,error_text) |
---|
| 2141 | else if (londimid.eq.d) then |
---|
| 2142 | status=nf90_def_dim(outfid, "longitude", lonlen, londimid) |
---|
| 2143 | error_text="Error: could not define the longitude dimension in the outfile" |
---|
| 2144 | call status_check(status,error_text) |
---|
| 2145 | endif |
---|
| 2146 | enddo |
---|
| 2147 | |
---|
| 2148 | !================================================================ |
---|
| 2149 | ! 2. Define "variables" and their attributes |
---|
| 2150 | !================================================================ |
---|
| 2151 | !================================================================ |
---|
| 2152 | ! 2.1 Write "longitude" (data and attributes) |
---|
| 2153 | !================================================================ |
---|
| 2154 | |
---|
| 2155 | ! Insert the definition of the variable |
---|
| 2156 | status=nf90_def_var(outfid,"longitude",nf90_float,(/londimid/),lonvarid) |
---|
| 2157 | error_text="Error: could not define the longitude variable in the outfile" |
---|
| 2158 | call status_check(status,error_text) |
---|
| 2159 | |
---|
| 2160 | ! Write the attributes |
---|
| 2161 | status=nf90_put_att(outfid,lonvarid,"long_name","longitude") |
---|
| 2162 | error_text="Error: could not put the long_name attribute for the longitude variable in the outfile" |
---|
| 2163 | call status_check(status,error_text) |
---|
| 2164 | |
---|
| 2165 | status=nf90_put_att(outfid,lonvarid,"units","degrees_east") |
---|
| 2166 | error_text="Error: could not put the units attribute for the longitude variable in the outfile" |
---|
| 2167 | call status_check(status,error_text) |
---|
| 2168 | |
---|
| 2169 | !================================================================ |
---|
| 2170 | ! 2.2 "latitude" |
---|
| 2171 | !================================================================ |
---|
| 2172 | |
---|
| 2173 | ! Insert the definition of the variable |
---|
| 2174 | status=nf90_def_var(outfid,"latitude",nf90_float,(/latdimid/),latvarid) |
---|
| 2175 | error_text="Error: could not define the latitude variable in the outfile" |
---|
| 2176 | call status_check(status,error_text) |
---|
| 2177 | |
---|
| 2178 | ! Write the attributes |
---|
| 2179 | status=nf90_put_att(outfid,latvarid,"long_name","latitude") |
---|
| 2180 | error_text="Error: could not put the long_name attribute for the latitude variable in the outfile" |
---|
| 2181 | call status_check(status,error_text) |
---|
| 2182 | |
---|
| 2183 | status=nf90_put_att(outfid,latvarid,"units","degrees_north") |
---|
| 2184 | error_text="Error: could not put the units attribute for the latitude variable in the outfile" |
---|
| 2185 | call status_check(status,error_text) |
---|
| 2186 | |
---|
| 2187 | !================================================================ |
---|
| 2188 | ! 2.3 Write "altitude" (data and attributes) |
---|
| 2189 | !================================================================ |
---|
| 2190 | |
---|
| 2191 | ! Insert the definition of the variable |
---|
| 2192 | status=nf90_def_var(outfid,"altitude",nf90_float,(/altdimid/),altvarid) |
---|
| 2193 | error_text="Error: could not define the altitude variable in the outfile" |
---|
| 2194 | call status_check(status,error_text) |
---|
| 2195 | |
---|
| 2196 | ! Write the attributes |
---|
| 2197 | if (units_alt.eq.'p') then ! pressure coordinate |
---|
| 2198 | status=nf90_put_att(outfid,altvarid,"long_name","pressure") |
---|
| 2199 | error_text="Error: could not put the long_name attribute for the altitude variable in the outfile" |
---|
| 2200 | call status_check(status,error_text) |
---|
| 2201 | |
---|
| 2202 | status=nf90_put_att(outfid,altvarid,'units',"Pa") |
---|
| 2203 | error_text="Error: could not put the units attribute for the altitude variable in the outfile" |
---|
| 2204 | call status_check(status,error_text) |
---|
| 2205 | |
---|
| 2206 | status=nf90_put_att(outfid,altvarid,'positive',"down") |
---|
| 2207 | error_text="Error: could not put the positive attribute for the altitude variable in the outfile" |
---|
| 2208 | call status_check(status,error_text) |
---|
| 2209 | |
---|
| 2210 | status=nf90_put_att(outfid,altvarid,'comment',& |
---|
| 2211 | "The vertical variable is in fact pressure, not altitude. We just call it altitude for easy reading in Ferret and Grads.") |
---|
| 2212 | error_text="Error: could not put the comment attribute for the altitude variable in the outfile" |
---|
| 2213 | call status_check(status,error_text) |
---|
| 2214 | |
---|
| 2215 | else if (units_alt.eq.'z') then ! altitude coordinate |
---|
| 2216 | status=nf90_put_att(outfid,altvarid,"long_name","altitude") |
---|
| 2217 | error_text="Error: could not put the long_name attribute for the altitude variable in the outfile" |
---|
| 2218 | call status_check(status,error_text) |
---|
| 2219 | |
---|
| 2220 | status=nf90_put_att(outfid,altvarid,'units',"m") |
---|
| 2221 | error_text="Error: could not put the units attribute for the altitude variable in the outfile" |
---|
| 2222 | call status_check(status,error_text) |
---|
| 2223 | |
---|
| 2224 | status=nf90_put_att(outfid,altvarid,'positive',"up") |
---|
| 2225 | error_text="Error: could not put the positive attribute for the altitude variable in the outfile" |
---|
| 2226 | call status_check(status,error_text) |
---|
| 2227 | |
---|
| 2228 | else |
---|
| 2229 | write(*,*)"I do not understand this unit type ",units_alt," for outfile altitude!" |
---|
| 2230 | stop |
---|
| 2231 | end if |
---|
| 2232 | |
---|
| 2233 | !================================================================ |
---|
| 2234 | ! 2.4 "time" |
---|
| 2235 | !================================================================ |
---|
| 2236 | |
---|
| 2237 | ! Insert the definition of the variable |
---|
| 2238 | status=nf90_def_var(outfid,"time",nf90_float,(/timedimid/),timevarid) |
---|
| 2239 | error_text="Error: could not define the time variable in the outfile" |
---|
| 2240 | call status_check(status,error_text) |
---|
| 2241 | |
---|
| 2242 | ! Write the attributes |
---|
| 2243 | status=nf90_put_att(outfid,timevarid,"long_name","Solar longitude") |
---|
| 2244 | error_text="Error: could not put the long_name attribute for the time variable in the outfile" |
---|
| 2245 | call status_check(status,error_text) |
---|
| 2246 | |
---|
| 2247 | status=nf90_put_att(outfid,timevarid,"units","days since 0000-01-1 00:00:00") |
---|
| 2248 | error_text="Error: could not put the units attribute for the time variable in the outfile" |
---|
| 2249 | call status_check(status,error_text) |
---|
| 2250 | |
---|
| 2251 | status=nf90_put_att(outfid,timevarid,"comment",& |
---|
| 2252 | "Units is in fact degrees, but set to a dummy value of days for compatibility with Ferret and Grads.") |
---|
| 2253 | error_text="Error: could not put the comment attribute for the time variable in the outfile" |
---|
| 2254 | call status_check(status,error_text) |
---|
| 2255 | |
---|
| 2256 | !================================================================ |
---|
| 2257 | ! 2.5 End netcdf define mode |
---|
| 2258 | !================================================================ |
---|
| 2259 | status=nf90_enddef(outfid) |
---|
| 2260 | error_text="Error: could not end the define mode of the outfile" |
---|
| 2261 | call status_check(status,error_text) |
---|
| 2262 | |
---|
| 2263 | !================================================================ |
---|
| 2264 | ! 3. Write "variables" (data of the dimension variables) |
---|
| 2265 | !================================================================ |
---|
| 2266 | ! "time" |
---|
| 2267 | status=nf90_put_var(outfid,timevarid,time) |
---|
| 2268 | error_text="Error: could not write the time variable's data in the outfile" |
---|
| 2269 | call status_check(status,error_text) |
---|
| 2270 | |
---|
| 2271 | ! "latitude" |
---|
| 2272 | status=nf90_put_var(outfid,latvarid,lat) |
---|
| 2273 | error_text="Error: could not write the latitude variable's data in the outfile" |
---|
| 2274 | call status_check(status,error_text) |
---|
| 2275 | |
---|
| 2276 | ! "longitude" |
---|
| 2277 | status=nf90_put_var(outfid,lonvarid,lon) |
---|
| 2278 | error_text="Error: could not write the longitude variable's data in the outfile" |
---|
| 2279 | call status_check(status,error_text) |
---|
| 2280 | |
---|
| 2281 | ! "altitude" |
---|
| 2282 | status=nf90_put_var(outfid,altvarid,alt) |
---|
| 2283 | error_text="Error: could not write the altitude variable's data in the outfile" |
---|
| 2284 | call status_check(status,error_text) |
---|
| 2285 | |
---|
| 2286 | end subroutine inidim |
---|
| 2287 | |
---|
| 2288 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2289 | |
---|
| 2290 | subroutine ls2sol(ls,sol) |
---|
| 2291 | |
---|
| 2292 | implicit none |
---|
| 2293 | !================================================================ |
---|
| 2294 | ! Arguments: |
---|
| 2295 | !================================================================ |
---|
| 2296 | real,intent(in) :: ls |
---|
| 2297 | real,intent(out) :: sol |
---|
| 2298 | |
---|
| 2299 | !================================================================ |
---|
| 2300 | ! Local: |
---|
| 2301 | !================================================================ |
---|
| 2302 | double precision xref,zx0,zteta,zz |
---|
| 2303 | !xref: mean anomaly, zteta: true anomaly, zx0: eccentric anomaly |
---|
| 2304 | double precision year_day |
---|
| 2305 | double precision peri_day,timeperi,e_elips |
---|
| 2306 | double precision pi,degrad |
---|
| 2307 | parameter (year_day=668.6d0) ! number of sols in a martian year |
---|
| 2308 | parameter (peri_day=485.35d0) ! date (in sols) of perihelion |
---|
| 2309 | !timeperi: 2*pi*( 1 - Ls(perihelion)/ 360 ); Ls(perihelion)=250.99 |
---|
| 2310 | parameter (timeperi=1.90258341759902d0) |
---|
| 2311 | parameter (e_elips=0.0934d0) ! eccentricity of orbit |
---|
| 2312 | parameter (pi=3.14159265358979d0) |
---|
| 2313 | parameter (degrad=57.2957795130823d0) |
---|
| 2314 | |
---|
| 2315 | if (abs(ls).lt.1.0e-5) then |
---|
| 2316 | if (ls.ge.0.0) then |
---|
| 2317 | sol = 0.0 |
---|
| 2318 | else |
---|
| 2319 | sol = real(year_day) |
---|
| 2320 | end if |
---|
| 2321 | return |
---|
| 2322 | end if |
---|
| 2323 | |
---|
| 2324 | zteta = ls/degrad + timeperi |
---|
| 2325 | zx0 = 2.0*datan(dtan(0.5*zteta)/dsqrt((1.+e_elips)/(1.-e_elips))) |
---|
| 2326 | xref = zx0-e_elips*dsin(zx0) |
---|
| 2327 | zz = xref/(2.*pi) |
---|
| 2328 | sol = real(zz*year_day + peri_day) |
---|
| 2329 | if (sol.lt.0.0) sol = sol + real(year_day) |
---|
| 2330 | if (sol.ge.year_day) sol = sol - real(year_day) |
---|
| 2331 | |
---|
| 2332 | |
---|
| 2333 | end subroutine ls2sol |
---|
| 2334 | |
---|
| 2335 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2336 | |
---|
| 2337 | subroutine gen_sol_list(solcount,int_sol_list,LTcount,LTave,LTmax,LTmin,lon,f_LT,dayornight,& |
---|
| 2338 | sol_list) |
---|
| 2339 | !================================================================ |
---|
| 2340 | ! Purpose: |
---|
| 2341 | ! Generate a list that is the combination of two lists : |
---|
| 2342 | ! - int_sol_list, which is a list of integer values of sol |
---|
| 2343 | ! - LT_list, which contains a number LTcount of LT values in the |
---|
| 2344 | ! interval [LTmin;LTmax] which are evenly distributed around |
---|
| 2345 | ! LTave (so that their average is LTave) |
---|
| 2346 | !================================================================ |
---|
| 2347 | |
---|
| 2348 | implicit none |
---|
| 2349 | !================================================================ |
---|
| 2350 | ! Arguments: |
---|
| 2351 | !================================================================ |
---|
| 2352 | integer, intent(in) :: solcount ! nb of integer values of sol |
---|
| 2353 | real, intent(in) :: int_sol_list(solcount) ! list of the integer values of sol |
---|
| 2354 | integer, intent(in) :: LTcount ! nb of LT per sol to be interpolated |
---|
| 2355 | real, intent(in) :: LTave ! average of LT |
---|
| 2356 | real, intent(in) :: LTmax, LTmin ! bounds of the LT interval for the interpolation |
---|
| 2357 | real, intent(in) :: lon ! longitude value for the transformation into a sol value at lon=0° |
---|
| 2358 | external f_LT ! function that changes LT interval for night LT |
---|
| 2359 | real f_LT |
---|
| 2360 | character (len=10), intent(in) :: dayornight ! to know if we have day or night values |
---|
| 2361 | |
---|
| 2362 | real, intent(out) :: sol_list(solcount*LTcount) ! all the sol values at lon=0° to interpolate |
---|
| 2363 | |
---|
| 2364 | !================================================================ |
---|
| 2365 | ! Local variables: |
---|
| 2366 | !================================================================ |
---|
| 2367 | integer :: N |
---|
| 2368 | integer :: a,b,c ! loop iteration indexes |
---|
| 2369 | real :: LT_list(LTcount) |
---|
| 2370 | |
---|
| 2371 | N = floor(LTcount/2.) |
---|
| 2372 | |
---|
| 2373 | !================================================================ |
---|
| 2374 | ! 1. Filling of LT_list |
---|
| 2375 | !================================================================ |
---|
| 2376 | SELECT CASE (dayornight) |
---|
| 2377 | CASE ("dayside") |
---|
| 2378 | if (abs(LTave-LTmax).le.abs(LTave-LTmin)) then ! LTave is closer to LTmax than to LTmin |
---|
| 2379 | do a=1,N |
---|
| 2380 | LT_list(a)=LTave+a*abs(LTave-LTmax)/(N+1) |
---|
| 2381 | LT_list(N+a)=LTave-a*abs(LTave-LTmax)/(N+1) |
---|
| 2382 | enddo |
---|
| 2383 | else ! LTave is closer to LTmin than to LTmax |
---|
| 2384 | do a=1,N |
---|
| 2385 | LT_list(a)=LTave+a*abs(LTave-LTmin)/(N+1) |
---|
| 2386 | LT_list(N+a)=LTave-a*abs(LTave-LTmin)/(N+1) |
---|
| 2387 | enddo |
---|
| 2388 | endif |
---|
| 2389 | CASE ("nightside") |
---|
| 2390 | if (abs(f_LT(LTave)-f_LT(LTmax)).le.abs(f_LT(LTave)-f_LT(LTmin))) then ! LTave is closer to LTmax than to LTmin |
---|
| 2391 | do a=1,N |
---|
| 2392 | LT_list(a)=LTave+a*abs(f_LT(LTave)-f_LT(LTmax))/(N+1) |
---|
| 2393 | LT_list(a)=mod(LT_list(a),24.) ! reput the LT in a [0;24[ interval if needed |
---|
| 2394 | LT_list(N+a)=LTave-a*abs(f_LT(LTave)-f_LT(LTmax))/(N+1) |
---|
| 2395 | LT_list(N+a)=mod(LT_list(N+a),24.) |
---|
| 2396 | enddo |
---|
| 2397 | else ! LTave is closer to LTmin than to LTmax |
---|
| 2398 | do a=1,N |
---|
| 2399 | LT_list(a)=LTave+a*abs(f_LT(LTave)-f_LT(LTmin))/(N+1) |
---|
| 2400 | LT_list(a)=mod(LT_list(a),24.) |
---|
| 2401 | LT_list(N+a)=LTave-a*abs(f_LT(LTave)-f_LT(LTmin))/(N+1) |
---|
| 2402 | LT_list(N+a)=mod(LT_list(N+a),24.) |
---|
| 2403 | enddo |
---|
| 2404 | endif |
---|
| 2405 | END SELECT |
---|
| 2406 | |
---|
| 2407 | if (mod(LTcount,2).eq.1) then ! if LTcount is an odd number |
---|
| 2408 | LT_list(LTcount)=LTave ! add LTave to the list |
---|
| 2409 | endif |
---|
| 2410 | |
---|
| 2411 | !================================================================ |
---|
| 2412 | ! 2. Combination of int_sol_list & LT_list into sol_list |
---|
| 2413 | !================================================================ |
---|
| 2414 | c=1 |
---|
| 2415 | do a=1,solcount |
---|
| 2416 | do b=1,LTcount |
---|
| 2417 | sol_list(c)=int_sol_list(a)+(LT_list(b)-lon/15.)/24. |
---|
| 2418 | c=c+1 |
---|
| 2419 | enddo |
---|
| 2420 | enddo |
---|
| 2421 | |
---|
| 2422 | end subroutine gen_sol_list |
---|
| 2423 | |
---|
| 2424 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2425 | |
---|
| 2426 | subroutine status_check(status,error_text) |
---|
| 2427 | |
---|
| 2428 | use netcdf |
---|
| 2429 | implicit none |
---|
| 2430 | !================================================================ |
---|
| 2431 | ! Arguments: |
---|
| 2432 | !================================================================ |
---|
| 2433 | integer,intent(in) :: status |
---|
| 2434 | character(len=256),intent(in) :: error_text |
---|
| 2435 | |
---|
| 2436 | if (status.ne.nf90_noerr) then |
---|
| 2437 | write(*,*)trim(error_text) |
---|
| 2438 | write(*,*)trim(nf90_strerror(status)) |
---|
| 2439 | stop |
---|
| 2440 | endif |
---|
| 2441 | |
---|
| 2442 | end subroutine status_check |
---|
| 2443 | |
---|
| 2444 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2445 | |
---|
| 2446 | real function LTmod(LT) |
---|
| 2447 | !================================================================ |
---|
| 2448 | ! Purpose: |
---|
| 2449 | ! For night local times management, replace hours |
---|
| 2450 | ! from a [0;24[ interval to a [-12;12[ interval in which |
---|
| 2451 | ! midnight = 0 (in order to ensure continuity when comparing |
---|
| 2452 | ! night local times) |
---|
| 2453 | !================================================================ |
---|
| 2454 | implicit none |
---|
| 2455 | real,intent(in) :: LT |
---|
| 2456 | |
---|
| 2457 | LTmod = 2*mod(LT,12.)-LT |
---|
| 2458 | return |
---|
| 2459 | end function LTmod |
---|