1 | program aeroptical |
---|
2 | ! program for computation of aerosol opacities |
---|
3 | ! author : Antoine Bierjon, April-May 2020 |
---|
4 | ! Rebranding in November 2022 |
---|
5 | ! |
---|
6 | !=============================================================================== |
---|
7 | ! PREFACE |
---|
8 | !=============================================================================== |
---|
9 | ! This program can be used to compute the opacity of any aerosol. |
---|
10 | ! |
---|
11 | ! It takes as input a GCM output file (diagfi,stats,concat) that |
---|
12 | ! contains: |
---|
13 | ! - the Mass Mixing Ratios of the aerosols to be computed |
---|
14 | ! - their effective radius (replaceable by a single value given by the user) |
---|
15 | ! - the atmospheric density |
---|
16 | ! The user also inputs the wavelength to calibrate the opacities, |
---|
17 | ! the type of opacity (extinction or absorption) and the directory |
---|
18 | ! containing the ASCII files of the aerosols' optical properties. |
---|
19 | ! |
---|
20 | ! Finally, for each aerosol to be computed, the user inputs : |
---|
21 | ! <aerosol_name>,<aerosol_mmr_varname>,<aerosol_reff_varname/value>,<aerosol_rho_value>,<optpropfile_name> |
---|
22 | ! |
---|
23 | ! The program outputs a netcdf file containing the opacities [1/km] of |
---|
24 | ! the aerosols calibrated at the prescribed wavelength. |
---|
25 | ! |
---|
26 | ! NB: this program requires to compile the module aeropt_mod.F90 at the same time |
---|
27 | !=============================================================================== |
---|
28 | |
---|
29 | |
---|
30 | use netcdf |
---|
31 | use aeropt_mod |
---|
32 | |
---|
33 | !=============================================================================== |
---|
34 | ! Variable declarations |
---|
35 | !=============================================================================== |
---|
36 | |
---|
37 | implicit none ! for no implicitly typed variables |
---|
38 | |
---|
39 | ! User inputs |
---|
40 | character(len=256) :: datadir ! directory containing the ASCII optical properties files |
---|
41 | real :: wvl_val ! reference wavelength of the output opacities [m] |
---|
42 | character(len=3) :: opatype ! opacity type : extinction (ext) or absorption (abs) |
---|
43 | character(len=3) :: compute_col ! does the user want to compute the column-integrated optical depth [yes/no] |
---|
44 | character(len=50),dimension(20) :: name_aer,mmrname_aer,reffname_aer ! aerosols' name, MMR varname, reff varname (or value) |
---|
45 | integer,dimension(20) :: reffval_ok ! to know if the user gave a variable name or a value for reff_aer |
---|
46 | real,dimension(20) :: reffval_aer ! value given by the user for reff_aer [m] |
---|
47 | real,dimension(20) :: rho_aer ! Density of the aerosols [kg.m-3] |
---|
48 | character(len=128),dimension(20) :: optpropfile_aer ! name of the ASCII optical properties file |
---|
49 | integer :: read_ok ! to know if the user input is well read |
---|
50 | integer :: naerkind ! nb of aerosols to consider |
---|
51 | |
---|
52 | ! GCM file |
---|
53 | character(len=128) :: gcmfile ! name of the netcdf GCM input file |
---|
54 | integer :: gcmfid ! [netcdf] GCM input file ID number |
---|
55 | integer :: ierr ! [netcdf] subroutine returned error code |
---|
56 | character(len=256) :: error_text ! text to display in case of error |
---|
57 | integer :: lonlen,latlen,altlen,timelen ! nb of grid points along longitude,latitude,altitude,time |
---|
58 | integer :: GCM_layers ! number of GCM layers |
---|
59 | integer :: layerdimout,interlayerdimout ! "GCM_layers" and "GCM_layers+1" IDs |
---|
60 | |
---|
61 | logical,dimension(:),allocatable :: aerok ! to know if the needed fields are in the GCM file |
---|
62 | integer,dimension(:),allocatable :: mmrvarid ! stores a MMR variable ID number |
---|
63 | integer,dimension(:),allocatable :: reffvarid ! stores a reff variable ID number |
---|
64 | integer :: tmpvarid ! temporarily stores a variable ID number |
---|
65 | real,dimension(:,:,:,:),allocatable :: mmr ! aerosols mass mixing ratios [kg/kg] |
---|
66 | real,dimension(:,:,:,:),allocatable :: reff ! aerosols effective radii [m] |
---|
67 | real,dimension(:,:,:,:),allocatable :: rho ! atmospheric density [kg.m-3] |
---|
68 | integer :: iaer ! Loop index for aerosols |
---|
69 | integer :: ilon,ilat,ialt,it ! Loop indices for coordinates |
---|
70 | |
---|
71 | ! Opacity computation |
---|
72 | character(len=128) :: aeroptlogfile ! name of the aeropt_mod log file |
---|
73 | real :: Qext_val ! Qext after interpolations |
---|
74 | real :: omeg_val ! omeg after interpolations |
---|
75 | |
---|
76 | ! Output file |
---|
77 | character(len=128) :: outfile ! name of the netcdf output file |
---|
78 | integer :: outfid ! [netcdf] file ID number |
---|
79 | integer :: latdimout,londimout,altdimout,timedimout |
---|
80 | ! latdimout: stores latitude dimension ID number |
---|
81 | ! londimout: stores longitude dimension ID number |
---|
82 | ! altdimout: stores altitude dimension ID number |
---|
83 | ! timedimout: stores time dimension ID number |
---|
84 | integer :: tmpvaridout ! temporarily stores a variable ID number |
---|
85 | integer :: varshape(4) ! stores a variable shape (of dimensions' IDs) |
---|
86 | character(len=16) :: tmpvarname ! temporarily stores a variable name |
---|
87 | character(len=128) :: tmplongname ! temporarily stores a variable long_name attribute |
---|
88 | real,dimension(:,:,:,:),allocatable :: opa_aer ! Aerosols opacities dtau/dz [1/km] |
---|
89 | real,dimension(:),allocatable :: missval ! Value to put in outfile when the reff is out of bounds |
---|
90 | ! or when there is a mising value in input file |
---|
91 | integer :: tauvarshape(3) ! stores the column variable shape (of dimensions' IDs) |
---|
92 | real,dimension(:,:,:,:),allocatable :: delta_z ! Layers' height for column computation [km] |
---|
93 | real,dimension(:,:,:),allocatable :: tau_aer ! Aerosols column-integrated optical depth [/] |
---|
94 | |
---|
95 | |
---|
96 | !=============================================================================== |
---|
97 | ! 0. USER INPUTS FOR NEW AEROPTICAL.DEF |
---|
98 | |
---|
99 | !1) Name of the GCM input file |
---|
100 | !2) Directory where the optical properties files are |
---|
101 | !3) The wavelength at which the output opacities are calibrated (in meters) |
---|
102 | !4) Opacity type : extinction (ext) or absorption (abs) |
---|
103 | !5) Does the user want to compute the column-integrated optical depth? [yes/no] |
---|
104 | !6) aerosol_name,aerosol_mmr(variable name),aerosol_reff(variable name or single value in m),aerosol_rho(value in kg/m3),optpropfile_name |
---|
105 | !7) aerosol_name,aerosol_mmr(variable name),aerosol_reff(variable name or single value in m),aerosol_rho(value in kg/m3),optpropfile_name |
---|
106 | !... |
---|
107 | !N) aerosol_name,aerosol_mmr(variable name),aerosol_reff(variable name or single value in m),aerosol_rho(value in kg/m3),optpropfile_name |
---|
108 | |
---|
109 | !=============================================================================== |
---|
110 | write(*,*) "Welcome in the aerosol opacities' computation program !" |
---|
111 | write(*,*) |
---|
112 | |
---|
113 | ! 1) Ask the GCM file name |
---|
114 | WRITE(*,*) "Enter an input file name (GCM diagfi/stats/concat) :" |
---|
115 | READ(*,*) gcmfile |
---|
116 | |
---|
117 | ! 2) Ask the directory containing the optical properties files |
---|
118 | write(*,*)"" |
---|
119 | WRITE(*,*) "In which directory should we look for the optical properties' files ?" |
---|
120 | READ(*,'(a)') datadir |
---|
121 | |
---|
122 | ! 3) Ask the wavelength to the user |
---|
123 | write(*,*)"" |
---|
124 | WRITE(*,*) "Value of the reference wavelength for the opacities' computation (in meters) ?" |
---|
125 | READ(*,*) wvl_val |
---|
126 | |
---|
127 | ! 4) Ask the type of opacity that has to be computed |
---|
128 | write(*,*)"" |
---|
129 | WRITE(*,*) "Do you want to compute extinction or absorption opacities ? (ext/abs)" |
---|
130 | READ(*,*) opatype |
---|
131 | if ((trim(opatype).ne."ext").and.(trim(opatype).ne."abs")) then |
---|
132 | write(*,*) "opatype = "//opatype |
---|
133 | write(*,*) "Error: the opacity type must be ext or abs" |
---|
134 | stop |
---|
135 | endif |
---|
136 | |
---|
137 | ! 5) Computation of the column-integrated optical depth? |
---|
138 | write(*,*)"" |
---|
139 | WRITE(*,*) "Do you want to compute the column-integrated optical depth ? (yes/no)" |
---|
140 | READ(*,*) compute_col |
---|
141 | if ((trim(compute_col).ne."yes").and.(trim(compute_col).ne."no")) then |
---|
142 | write(*,*) "compute_col = "//compute_col |
---|
143 | write(*,*) "Error: please answer yes or no" |
---|
144 | stop |
---|
145 | endif |
---|
146 | |
---|
147 | |
---|
148 | ! 6-N) Ask for aerosols to compute to the user |
---|
149 | write(*,*)"" |
---|
150 | write(*,*) "For which aerosols do you want to compute the opacity?" |
---|
151 | write(*,*) "List of aerosols (up to 20), separated by <Enter>s, with each line containing :" |
---|
152 | write(*,*) "<aerosol_name>,<aerosol_mmr>,<aerosol_reff>,<aerosol_rho>,<optpropfile_name>" |
---|
153 | write(*,*) "(aerosol_mmr is a variable name, aerosol_reff can be variable name or value in meter ; aerosol_rho is a value in kg/m3)" |
---|
154 | write(*,*)"" |
---|
155 | write(*,*) "An empty line , i.e: just <Enter>, implies end of list." |
---|
156 | |
---|
157 | reffval_ok(:) = 0 ! initialize reffval_ok to 0 for all aerosols |
---|
158 | |
---|
159 | naerkind=0 |
---|
160 | read(*,*,iostat=read_ok) name_aer(1),mmrname_aer(1),reffname_aer(1),rho_aer(1),optpropfile_aer(1) |
---|
161 | |
---|
162 | do while ((read_ok.eq.0).and.(name_aer(naerkind+1)/="")) |
---|
163 | naerkind=naerkind+1 |
---|
164 | |
---|
165 | ! Check if reffname_aer is a variable name or a value (in meters) |
---|
166 | read(reffname_aer(naerkind),*,iostat=reffval_ok(naerkind)) reffval_aer(naerkind) |
---|
167 | if (reffval_ok(naerkind).eq.0) then |
---|
168 | write(*,*) "This value for reff of ",trim(name_aer(naerkind))," will be broadcasted to the whole grid :",reffval_aer(naerkind) |
---|
169 | endif |
---|
170 | |
---|
171 | ! Read next line |
---|
172 | read(*,*,iostat=read_ok) name_aer(naerkind+1),mmrname_aer(naerkind+1),reffname_aer(naerkind+1),rho_aer(naerkind+1),optpropfile_aer(naerkind+1) |
---|
173 | enddo |
---|
174 | |
---|
175 | if(naerkind==0) then |
---|
176 | write(*,*) "no aerosol... better stop now." |
---|
177 | stop |
---|
178 | endif |
---|
179 | |
---|
180 | write(*,*) "naerkind=",naerkind |
---|
181 | write(*,*) "name_aer=",name_aer(1:naerkind) |
---|
182 | write(*,*) "mmrname_aer=",mmrname_aer(1:naerkind) |
---|
183 | write(*,*) "reffname_aer=",reffname_aer(1:naerkind) |
---|
184 | write(*,*) "rho_aer=",rho_aer(1:naerkind) |
---|
185 | write(*,*) "optpropfile_aer=",optpropfile_aer(1:naerkind) |
---|
186 | |
---|
187 | |
---|
188 | !=============================================================================== |
---|
189 | ! 1. GCM FILE & OUTPUT FILE INITIALIZATIONS |
---|
190 | !=============================================================================== |
---|
191 | |
---|
192 | !========================================================================== |
---|
193 | ! 1.1 GCM file opening & dimensions - Output file initializations |
---|
194 | !========================================================================== |
---|
195 | |
---|
196 | !========================================================================== |
---|
197 | ! --> 1.1.1 Open the netcdf GCM file given by the user |
---|
198 | |
---|
199 | ierr=nf90_open(gcmfile,nf90_nowrite,gcmfid) ! nowrite mode=the program can only read the opened file |
---|
200 | error_text="Error: could not open file "//trim(gcmfile) |
---|
201 | call status_check(ierr,error_text) |
---|
202 | |
---|
203 | !========================================================================== |
---|
204 | ! --> 1.1.2 Creation of the outfile |
---|
205 | if (opatype.eq."ext") then |
---|
206 | outfile=gcmfile(1:index(gcmfile, ".nc")-1)//"_OPAext.nc" |
---|
207 | else ! opatype.eq."abs" |
---|
208 | outfile=gcmfile(1:index(gcmfile, ".nc")-1)//"_OPAabs.nc" |
---|
209 | endif |
---|
210 | |
---|
211 | |
---|
212 | ierr=NF90_CREATE(outfile,IOR(nf90_clobber,nf90_64bit_offset),outfid) |
---|
213 | ! NB: clobber mode=overwrite any dataset with the same file name ! |
---|
214 | ! 64bit_offset enables the creation of large netcdf files with variables>2GB |
---|
215 | error_text="Error: could not create file "//trim(outfile) |
---|
216 | call status_check(ierr,error_text) |
---|
217 | write(*,*)"";WRITE(*,*)"Output file is: ",trim(outfile);write(*,*)"" |
---|
218 | |
---|
219 | !========================================================================== |
---|
220 | ! --> 1.1.3 Get the dimensions and create them in the output file |
---|
221 | |
---|
222 | ! Initialize output file's lat,lon,alt,time and controle dimensions |
---|
223 | call inidims(gcmfile,gcmfid,outfile,outfid,& |
---|
224 | lonlen,latlen,altlen,timelen,& |
---|
225 | latdimout,londimout,altdimout,timedimout,& |
---|
226 | GCM_layers,layerdimout,interlayerdimout) |
---|
227 | |
---|
228 | ! Length allocation for each dimension of the 4D variable delta_z |
---|
229 | allocate(delta_z(lonlen,latlen,altlen,timelen)) |
---|
230 | |
---|
231 | ! Initialize output file's aps,bps,ap,bp and phisinit variables |
---|
232 | ! + compute delta_z for the column integration |
---|
233 | call init2(gcmfid,lonlen,latlen,altlen,timelen,GCM_layers,& |
---|
234 | outfid,londimout,latdimout,altdimout,timedimout,& |
---|
235 | layerdimout,interlayerdimout,& |
---|
236 | compute_col,delta_z) |
---|
237 | |
---|
238 | !========================================================================== |
---|
239 | ! 1.2 Check GCM aerosols |
---|
240 | !========================================================================== |
---|
241 | |
---|
242 | ! Variables length allocation |
---|
243 | allocate(mmrvarid(naerkind)) |
---|
244 | allocate(reffvarid(naerkind)) |
---|
245 | |
---|
246 | ! Initialize missing value |
---|
247 | allocate(missval(naerkind)) |
---|
248 | missval(:)=1.e+20 |
---|
249 | |
---|
250 | ! Initialize aerok to .true. for all aerosols |
---|
251 | allocate(aerok(naerkind)) |
---|
252 | aerok(:)=.true. |
---|
253 | |
---|
254 | ! Loop on all aerosols |
---|
255 | DO iaer = 1, naerkind |
---|
256 | ! MASS MIXING RATIO |
---|
257 | ! Check that the GCM file contains that variable |
---|
258 | ierr=nf90_inq_varid(gcmfid,mmrname_aer(iaer),mmrvarid(iaer)) |
---|
259 | error_text="Failed to find variable "//trim(mmrname_aer(iaer))//& |
---|
260 | " in "//trim(gcmfile)//& |
---|
261 | " We'll skip the "//name_aer(iaer)//" aerosol." |
---|
262 | if (ierr.ne.nf90_noerr) then |
---|
263 | write(*,*)trim(error_text) |
---|
264 | aerok(iaer)=.false. |
---|
265 | endif |
---|
266 | |
---|
267 | ! EFFECTIVE RADIUS |
---|
268 | if (aerok(iaer)) then |
---|
269 | IF (reffval_ok(iaer).ne.0) THEN ! if the user gave a variable's name and not a value |
---|
270 | ! Check that the GCM file contains that variable |
---|
271 | ierr=nf90_inq_varid(gcmfid,reffname_aer(iaer),reffvarid(iaer)) |
---|
272 | error_text="Failed to find variable "//trim(reffname_aer(iaer))//& |
---|
273 | " in "//trim(gcmfile)//& |
---|
274 | " We'll skip the "//trim(name_aer(iaer))//" aerosol." |
---|
275 | if (ierr.ne.nf90_noerr) then |
---|
276 | write(*,*)trim(error_text) |
---|
277 | aerok(iaer)=.false. |
---|
278 | endif |
---|
279 | ENDIF |
---|
280 | endif |
---|
281 | ENDDO !iaer |
---|
282 | |
---|
283 | ! Check if there is still at least one aerosol to compute |
---|
284 | IF (.NOT.ANY(aerok)) THEN |
---|
285 | write(*,*) "No aerosol among those listed was found in the file. Better stop now..." |
---|
286 | stop |
---|
287 | ENDIF |
---|
288 | |
---|
289 | !========================================================================== |
---|
290 | ! 1.3 Get GCM atmospheric density |
---|
291 | !========================================================================== |
---|
292 | |
---|
293 | ! Check that the GCM file contains that variable |
---|
294 | ierr=nf90_inq_varid(gcmfid,"rho",tmpvarid) |
---|
295 | error_text="Failed to find variable rho in "//trim(gcmfile)& |
---|
296 | //" We need it to compute the opacity [1/km]." |
---|
297 | call status_check(ierr,error_text) |
---|
298 | ! Length allocation for each dimension of the 4D variable |
---|
299 | allocate(rho(lonlen,latlen,altlen,timelen)) |
---|
300 | ! Load dataset |
---|
301 | ierr=nf90_get_var(gcmfid,tmpvarid,rho) |
---|
302 | error_text="Failed to load atmospheric density" |
---|
303 | call status_check(ierr,error_text) |
---|
304 | write(*,*) "Atmospheric density rho loaded from "//trim(gcmfile) |
---|
305 | |
---|
306 | |
---|
307 | !=============================================================================== |
---|
308 | ! 2. OUTPUT FILE VARIABLES |
---|
309 | !=============================================================================== |
---|
310 | !========================================================================== |
---|
311 | ! 2.1 Output variable's initializations and datasets loading |
---|
312 | !========================================================================== |
---|
313 | ! Define the ID shape of the output variables |
---|
314 | varshape(1)=londimout |
---|
315 | varshape(2)=latdimout |
---|
316 | varshape(3)=altdimout |
---|
317 | varshape(4)=timedimout |
---|
318 | |
---|
319 | tauvarshape(1)=londimout |
---|
320 | tauvarshape(2)=latdimout |
---|
321 | tauvarshape(3)=timedimout |
---|
322 | |
---|
323 | |
---|
324 | ! Loop on all aerosols |
---|
325 | DO iaer = 1, naerkind |
---|
326 | if (aerok(iaer)) then ! check if this aerosol can be computed |
---|
327 | write(*,*)"" |
---|
328 | write(*,*)"Computing ",trim(name_aer(iaer))," opacity..." |
---|
329 | |
---|
330 | ! Length allocation of the input variables |
---|
331 | ALLOCATE(mmr(lonlen,latlen,altlen,timelen)) |
---|
332 | ALLOCATE(reff(lonlen,latlen,altlen,timelen)) |
---|
333 | |
---|
334 | ! Datasets loading |
---|
335 | ! ->MMR |
---|
336 | ierr=NF90_GET_VAR(gcmfid,mmrvarid(iaer),mmr(:,:,:,:)) |
---|
337 | error_text="Failed to load mass mixing ratio" |
---|
338 | call status_check(ierr,error_text) |
---|
339 | write(*,*) "Mass mixing ratio loaded from "//trim(gcmfile) |
---|
340 | ! Get missing value |
---|
341 | ierr=nf90_get_att(gcmfid,mmrvarid(iaer),"missing_value",missval(iaer)) |
---|
342 | if (ierr.ne.nf90_noerr) then |
---|
343 | missval(iaer) = 1.e+20 |
---|
344 | endif |
---|
345 | |
---|
346 | ! ->REFF |
---|
347 | IF (reffval_ok(iaer).ne.0) THEN ! if the user gave a variable's name and not a value |
---|
348 | ! Load dataset |
---|
349 | ierr=NF90_GET_VAR(gcmfid,reffvarid(iaer),reff(:,:,:,:)) |
---|
350 | error_text="Failed to load effective radius" |
---|
351 | call status_check(ierr,error_text) |
---|
352 | write(*,*) "Effective radius loaded from "//trim(gcmfile) |
---|
353 | ELSE ! if reffval_ok(iaer)=0 |
---|
354 | reff(:,:,:,:) = reffval_aer(iaer) |
---|
355 | write(*,*) "Effective radius loaded from the user input" |
---|
356 | ENDIF |
---|
357 | |
---|
358 | |
---|
359 | ! Length allocation of the output variables |
---|
360 | ALLOCATE(opa_aer(lonlen,latlen,altlen,timelen)) |
---|
361 | |
---|
362 | if (trim(compute_col).eq."yes") then |
---|
363 | ALLOCATE(tau_aer(lonlen,latlen,timelen)) |
---|
364 | tau_aer(:,:,:) = 0. |
---|
365 | endif ! trim(compute_col).eq."yes" |
---|
366 | |
---|
367 | |
---|
368 | !========================================================================== |
---|
369 | ! 2.2 Reading optical properties |
---|
370 | !========================================================================== |
---|
371 | |
---|
372 | write(*,*)"" |
---|
373 | ! create the aeropt_mod log file for this aerosol |
---|
374 | write(aeroptlogfile,*) "aeropt_log_",trim(name_aer(iaer)),".txt" |
---|
375 | aeroptlogfileID = 100+iaer ! update the log file unit |
---|
376 | CALL create_logfile(aeroptlogfile) |
---|
377 | |
---|
378 | write(*,*)"Reading the optical properties..." |
---|
379 | ! Read the properties tables |
---|
380 | CALL read_optpropfile(datadir,optpropfile_aer(iaer)) |
---|
381 | |
---|
382 | !========================================================================== |
---|
383 | ! 2.3 Creation of the output aerosol opacity variable |
---|
384 | !========================================================================== |
---|
385 | ! Switch to netcdf define mode |
---|
386 | ierr=nf90_redef(outfid) |
---|
387 | error_text="ERROR: Couldn't start netcdf define mode" |
---|
388 | call status_check(ierr,error_text) |
---|
389 | write(*,*)"" |
---|
390 | |
---|
391 | ! Insert the definition of the variable |
---|
392 | tmpvarname="opa_"//trim(name_aer(iaer)) |
---|
393 | ierr=NF90_DEF_VAR(outfid,trim(tmpvarname),nf90_float,varshape,tmpvaridout) |
---|
394 | error_text="ERROR: Couldn't create "//trim(tmpvarname)//" in the outfile" |
---|
395 | call status_check(ierr,error_text) |
---|
396 | write(*,*) trim(tmpvarname)//" has been created in the outfile" |
---|
397 | |
---|
398 | ! Write the attributes |
---|
399 | if (opatype.eq."ext") then |
---|
400 | tmplongname=trim(name_aer(iaer))//" extinction opacity at reference wavelength" |
---|
401 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",tmplongname) |
---|
402 | else ! opatype.eq."abs" |
---|
403 | tmplongname=trim(name_aer(iaer))//" absorption opacity at reference wavelength" |
---|
404 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",tmplongname) |
---|
405 | endif |
---|
406 | |
---|
407 | ierr=nf90_put_att(outfid,tmpvaridout,"units","opacity/km") |
---|
408 | ierr=nf90_put_att(outfid,tmpvaridout,"refwavelength",wvl_val) |
---|
409 | ierr=nf90_put_att(outfid,tmpvaridout,"missing_value",missval(iaer)) |
---|
410 | write(*,*)"with missing value = ",missval(iaer) |
---|
411 | |
---|
412 | ! End netcdf define mode |
---|
413 | ierr=nf90_enddef(outfid) |
---|
414 | error_text="ERROR: Couldn't end netcdf define mode" |
---|
415 | call status_check(ierr,error_text) |
---|
416 | |
---|
417 | !========================================================================== |
---|
418 | ! 2.4 Computation of the opacities |
---|
419 | !========================================================================== |
---|
420 | IF (reffval_ok(iaer).eq.0) THEN |
---|
421 | ! if the user gave a reff value, |
---|
422 | ! do the interpolation only once |
---|
423 | CALL interp_wvl_reff(wvl_val,reffval_aer(iaer),missval(iaer),& |
---|
424 | Qext_val,omeg_val) |
---|
425 | write(*,*) "Qext(reff_val)=",Qext_val," ; omeg(reff_val)=",omeg_val |
---|
426 | ENDIF |
---|
427 | |
---|
428 | do ilon=1,lonlen |
---|
429 | do ilat=1,latlen |
---|
430 | do it=1,timelen |
---|
431 | do ialt=1,altlen |
---|
432 | IF (reffval_ok(iaer).ne.0) THEN ! if the user gave a variable's name and not a value |
---|
433 | CALL interp_wvl_reff(wvl_val,reff(ilon,ilat,ialt,it),missval(iaer),& |
---|
434 | Qext_val,omeg_val) |
---|
435 | ENDIF |
---|
436 | |
---|
437 | if ((Qext_val.eq.missval(iaer))) then |
---|
438 | ! if the user's effective radius is |
---|
439 | ! out of the optpropfile boundaries |
---|
440 | opa_aer(ilon,ilat,ialt,it)=missval(iaer) |
---|
441 | write(aeroptlogfileID,*) "(mmr = ",mmr(ilon,ilat,ialt,it),"kg/kg)" |
---|
442 | write(aeroptlogfileID,*)"" |
---|
443 | |
---|
444 | else if ((mmr(ilon,ilat,ialt,it).eq.missval(iaer))) then |
---|
445 | ! if there is a missing value in input file |
---|
446 | opa_aer(ilon,ilat,ialt,it)=missval(iaer) |
---|
447 | |
---|
448 | else |
---|
449 | ! COMPUTATION OF THE EXTINCTION OPACITY [1/km] |
---|
450 | opa_aer(ilon,ilat,ialt,it)= 750.*Qext_val*mmr(ilon,ilat,ialt,it)*rho(ilon,ilat,ialt,it)& |
---|
451 | / ( rho_aer(iaer) * reff(ilon,ilat,ialt,it) ) |
---|
452 | |
---|
453 | if (opatype.eq."abs") then |
---|
454 | ! COMPUTATION OF THE ABSORPTION OPACITY [1/km] |
---|
455 | opa_aer(ilon,ilat,ialt,it)= (1 - omeg_val) * opa_aer(ilon,ilat,ialt,it) |
---|
456 | endif ! opatype.eq."abs" |
---|
457 | |
---|
458 | endif |
---|
459 | |
---|
460 | if (trim(compute_col).eq."yes") then |
---|
461 | ! COMPUTATION OF THE COLUMN-INTEGRATED OPTICAL DEPTH [/] |
---|
462 | if (opa_aer(ilon,ilat,ialt,it).ne.missval(iaer)) then |
---|
463 | tau_aer(ilon,ilat,it)= tau_aer(ilon,ilat,it)+opa_aer(ilon,ilat,ialt,it)*delta_z(ilon,ilat,ialt,it) |
---|
464 | endif |
---|
465 | endif ! trim(compute_col).eq."yes" |
---|
466 | |
---|
467 | enddo ! ialt |
---|
468 | |
---|
469 | if (trim(compute_col).eq."yes") then |
---|
470 | ! COMPUTATION OF THE COLUMN-INTEGRATED OPTICAL DEPTH [/] |
---|
471 | if (tau_aer(ilon,ilat,it).eq.0) then |
---|
472 | tau_aer(ilon,ilat,it)= missval(iaer) |
---|
473 | endif |
---|
474 | endif ! trim(compute_col).eq."yes" |
---|
475 | |
---|
476 | enddo ! it |
---|
477 | enddo ! ilat |
---|
478 | enddo ! ilon |
---|
479 | |
---|
480 | !========================================================================== |
---|
481 | ! 3.3 Writing opacity in the output file |
---|
482 | !========================================================================== |
---|
483 | ierr = NF90_PUT_VAR(outfid, tmpvaridout, opa_aer(:,:,:,:)) |
---|
484 | error_text="Error: could not write "//trim(tmpvarname)//" data in the outfile" |
---|
485 | call status_check(ierr,error_text) |
---|
486 | |
---|
487 | !========================================================================== |
---|
488 | ! 3.4 Writing column-integrated optical depth in the output file |
---|
489 | !========================================================================== |
---|
490 | if (trim(compute_col).eq."yes") then |
---|
491 | ! Switch to netcdf define mode |
---|
492 | ierr=nf90_redef(outfid) |
---|
493 | error_text="ERROR: Couldn't start netcdf define mode" |
---|
494 | call status_check(ierr,error_text) |
---|
495 | write(*,*)"" |
---|
496 | |
---|
497 | ! Insert the definition of the variable |
---|
498 | tmpvarname="tau_"//trim(name_aer(iaer)) |
---|
499 | ierr=NF90_DEF_VAR(outfid,trim(tmpvarname),nf90_float,tauvarshape,tmpvaridout) |
---|
500 | error_text="ERROR: Couldn't create "//trim(tmpvarname)//" in the outfile" |
---|
501 | call status_check(ierr,error_text) |
---|
502 | write(*,*) trim(tmpvarname)//" has been created in the outfile" |
---|
503 | |
---|
504 | ! Write the attributes |
---|
505 | if (opatype.eq."ext") then |
---|
506 | tmplongname=trim(name_aer(iaer))//" extinction column-integrated optical depth at ref. wavelength" |
---|
507 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",tmplongname) |
---|
508 | else ! opatype.eq."abs" |
---|
509 | tmplongname=trim(name_aer(iaer))//" absorption column-integrated optical depth at ref. wavelength" |
---|
510 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",tmplongname) |
---|
511 | endif |
---|
512 | |
---|
513 | ierr=nf90_put_att(outfid,tmpvaridout,"units","/") |
---|
514 | ierr=nf90_put_att(outfid,tmpvaridout,"refwavelength",wvl_val) |
---|
515 | ierr=nf90_put_att(outfid,tmpvaridout,"missing_value",missval(iaer)) |
---|
516 | write(*,*)"with missing value = ",missval(iaer) |
---|
517 | |
---|
518 | ! End netcdf define mode |
---|
519 | ierr=nf90_enddef(outfid) |
---|
520 | error_text="ERROR: Couldn't end netcdf define mode" |
---|
521 | call status_check(ierr,error_text) |
---|
522 | |
---|
523 | |
---|
524 | ! Writing variable in output file |
---|
525 | ierr = NF90_PUT_VAR(outfid, tmpvaridout, tau_aer(:,:,:)) |
---|
526 | error_text="Error: could not write "//trim(tmpvarname)//" data in the outfile" |
---|
527 | call status_check(ierr,error_text) |
---|
528 | endif ! trim(compute_col).eq."yes" |
---|
529 | |
---|
530 | |
---|
531 | !========================================================================== |
---|
532 | ! 3.5 Prepare next loop |
---|
533 | !========================================================================== |
---|
534 | DEALLOCATE(mmr) |
---|
535 | DEALLOCATE(reff) |
---|
536 | DEALLOCATE(opa_aer) |
---|
537 | IF (allocated(tau_aer)) DEALLOCATE(tau_aer) |
---|
538 | |
---|
539 | CALL end_aeropt_mod |
---|
540 | |
---|
541 | endif ! if aerok(iaer) |
---|
542 | |
---|
543 | ENDDO ! iaer |
---|
544 | |
---|
545 | !=============================================================================== |
---|
546 | ! 4. Close the files and end the program |
---|
547 | !=============================================================================== |
---|
548 | ierr = nf90_close(gcmfid) |
---|
549 | error_text="Error: could not close file "//trim(gcmfile) |
---|
550 | call status_check(ierr,error_text) |
---|
551 | |
---|
552 | ierr = nf90_close(outfid) |
---|
553 | error_text="Error: could not close file "//trim(outfile) |
---|
554 | call status_check(ierr,error_text) |
---|
555 | |
---|
556 | write(*,*)"";write(*,*)"Program aeroptical completed!" |
---|
557 | |
---|
558 | end program aeroptical |
---|
559 | |
---|
560 | |
---|
561 | |
---|
562 | |
---|
563 | |
---|
564 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
565 | !! SUBROUTINES |
---|
566 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
567 | |
---|
568 | subroutine status_check(ierr,error_text) |
---|
569 | |
---|
570 | use netcdf |
---|
571 | implicit none |
---|
572 | !================================================================ |
---|
573 | ! Arguments: |
---|
574 | !================================================================ |
---|
575 | integer,intent(in) :: ierr ! NetCDF status number |
---|
576 | character(len=256),intent(in) :: error_text |
---|
577 | |
---|
578 | if (ierr.ne.nf90_noerr) then |
---|
579 | write(*,*)trim(error_text) |
---|
580 | write(*,*)trim(nf90_strerror(ierr)) |
---|
581 | stop |
---|
582 | endif |
---|
583 | |
---|
584 | end subroutine status_check |
---|
585 | |
---|
586 | |
---|
587 | !******************************************************************************* |
---|
588 | |
---|
589 | subroutine inidims(gcmfile,gcmfid,outfile,outfid,lonlen,latlen,altlen,timelen,& |
---|
590 | latdimout,londimout,altdimout,timedimout,& |
---|
591 | GCM_layers,layerdimout,interlayerdimout) |
---|
592 | |
---|
593 | !============================================================================== |
---|
594 | ! Purpose: |
---|
595 | ! Read the dimensions of the input file |
---|
596 | ! and write them in the output file |
---|
597 | !============================================================================== |
---|
598 | ! Remarks: |
---|
599 | ! The NetCDF files must be open |
---|
600 | !============================================================================== |
---|
601 | use netcdf |
---|
602 | |
---|
603 | implicit none |
---|
604 | |
---|
605 | !============================================================================== |
---|
606 | ! Arguments: |
---|
607 | !============================================================================== |
---|
608 | character(len=128),intent(in) :: gcmfile ! name of the netcdf GCM input file |
---|
609 | integer,intent(in) :: gcmfid ! [netcdf] GCM input file ID number |
---|
610 | character(len=128),intent(in) :: outfile ! name of the netcdf output file |
---|
611 | integer,intent(in) :: outfid ! [netcdf] file ID number |
---|
612 | |
---|
613 | integer,intent(out) :: lonlen,latlen,altlen,timelen |
---|
614 | ! nb of grid points along longitude,latitude,altitude,time |
---|
615 | integer,intent(out) :: latdimout,londimout,altdimout,timedimout |
---|
616 | ! latdimout: stores latitude dimension ID number |
---|
617 | ! londimout: stores longitude dimension ID number |
---|
618 | ! altdimout: stores altitude dimension ID number |
---|
619 | ! timedimout: stores time dimension ID number |
---|
620 | integer,intent(out) :: GCM_layers ! number of GCM layers |
---|
621 | integer,intent(out) :: layerdimout ! "GCM_layers" ID |
---|
622 | integer,intent(out) :: interlayerdimout ! "GCM_layers+1" ID |
---|
623 | |
---|
624 | !============================================================================== |
---|
625 | ! Local variables: |
---|
626 | !============================================================================== |
---|
627 | integer :: ierr ! [netcdf] subroutine returned error code |
---|
628 | character(len=256) :: error_text ! text to display in case of error |
---|
629 | |
---|
630 | integer :: tmpdimid,tmpvarid ! temporary store a dimension/variable ID number |
---|
631 | character (len=64) :: long_name,units,positive |
---|
632 | ! long_name(): [netcdf] long_name attribute |
---|
633 | ! units(): [netcdf] units attribute |
---|
634 | ! positive(): [netcdf] positive direction attribute (for altitude) |
---|
635 | real, dimension(:), allocatable:: lat,lon,alt,time,ctl |
---|
636 | ! lat(): array, stores latitude coordinates |
---|
637 | ! lon(): array, stores longitude coordinates |
---|
638 | ! alt(): array, stores altitude coordinates |
---|
639 | ! time(): array, stores time coordinates |
---|
640 | ! ctl(): array, stores controle variable |
---|
641 | integer :: ctllen ! nb of elements in the controle array |
---|
642 | integer :: tmpdimidout,tmpvaridout ! temporary stores a dimension/variable ID number |
---|
643 | |
---|
644 | !============================================================================== |
---|
645 | ! LONGITUDE |
---|
646 | !============================================================================== |
---|
647 | ! Get the dimension in GCM file |
---|
648 | ierr=nf90_inq_dimid(gcmfid,"longitude",tmpdimid) |
---|
649 | error_text="Error: Dimension <longitude> is missing in file "//trim(gcmfile) |
---|
650 | call status_check(ierr,error_text) |
---|
651 | |
---|
652 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=lonlen) |
---|
653 | allocate(lon(lonlen)) |
---|
654 | |
---|
655 | ! Create the dimension in output file |
---|
656 | ierr=NF90_DEF_DIM(outfid,"longitude",lonlen,londimout) |
---|
657 | error_text="Error: could not define the longitude dimension in the outfile" |
---|
658 | call status_check(ierr,error_text) |
---|
659 | |
---|
660 | ! Get the field in GCM file |
---|
661 | ierr=nf90_inq_varid(gcmfid,"longitude",tmpvarid) |
---|
662 | error_text="Error: Field <longitude> is missing in file "//trim(gcmfile) |
---|
663 | call status_check(ierr,error_text) |
---|
664 | |
---|
665 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,lon) |
---|
666 | error_text="Failed to load longitude" |
---|
667 | call status_check(ierr,error_text) |
---|
668 | |
---|
669 | ! Create the field in the output file |
---|
670 | ierr=NF90_DEF_VAR(outfid,"longitude",nf90_float,(/londimout/),tmpvaridout) |
---|
671 | error_text="Error: could not define the longitude variable in the outfile" |
---|
672 | call status_check(ierr,error_text) |
---|
673 | |
---|
674 | ! Get the field attributes in the GCM file |
---|
675 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
676 | if (ierr.ne.nf90_noerr) then |
---|
677 | ! if no attribute "long_name", try "title" |
---|
678 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
679 | endif |
---|
680 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
681 | |
---|
682 | ! Put the field attributes in the output file |
---|
683 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
684 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
685 | |
---|
686 | ! Write the field values in the output file |
---|
687 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
688 | error_text="Error: could not end the define mode of the outfile" |
---|
689 | call status_check(ierr,error_text) |
---|
690 | |
---|
691 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,lon) |
---|
692 | error_text="Error: could not write the longitude field in the outfile" |
---|
693 | call status_check(ierr,error_text) |
---|
694 | |
---|
695 | !============================================================================== |
---|
696 | ! LATITUDE |
---|
697 | !============================================================================== |
---|
698 | ! Switch to netcdf define mode |
---|
699 | ierr=nf90_redef(outfid) |
---|
700 | error_text="Error: could not switch to define mode in the outfile" |
---|
701 | call status_check(ierr,error_text) |
---|
702 | |
---|
703 | ! Get the dimension in GCM file |
---|
704 | ierr=nf90_inq_dimid(gcmfid,"latitude",tmpdimid) |
---|
705 | error_text="Error: Dimension <latitude> is missing in file "//trim(gcmfile) |
---|
706 | call status_check(ierr,error_text) |
---|
707 | |
---|
708 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=latlen) |
---|
709 | allocate(lat(latlen)) |
---|
710 | |
---|
711 | ! Create the dimension in output file |
---|
712 | ierr=NF90_DEF_DIM(outfid,"latitude",latlen,latdimout) |
---|
713 | error_text="Error: could not define the latitude dimension in the outfile" |
---|
714 | call status_check(ierr,error_text) |
---|
715 | |
---|
716 | ! Get the field in GCM file |
---|
717 | ierr=nf90_inq_varid(gcmfid,"latitude",tmpvarid) |
---|
718 | error_text="Error: Field <latitude> is missing in file "//trim(gcmfile) |
---|
719 | call status_check(ierr,error_text) |
---|
720 | |
---|
721 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,lat) |
---|
722 | error_text="Failed to load latitude" |
---|
723 | call status_check(ierr,error_text) |
---|
724 | |
---|
725 | ! Create the field in the output file |
---|
726 | ierr=NF90_DEF_VAR(outfid,"latitude",nf90_float,(/latdimout/),tmpvaridout) |
---|
727 | error_text="Error: could not define the latitude variable in the outfile" |
---|
728 | call status_check(ierr,error_text) |
---|
729 | |
---|
730 | ! Get the field attributes in the GCM file |
---|
731 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
732 | if (ierr.ne.nf90_noerr) then |
---|
733 | ! if no attribute "long_name", try "title" |
---|
734 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
735 | endif |
---|
736 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
737 | |
---|
738 | ! Put the field attributes in the output file |
---|
739 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
740 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
741 | |
---|
742 | ! Write the field values in the output file |
---|
743 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
744 | error_text="Error: could not end the define mode of the outfile" |
---|
745 | call status_check(ierr,error_text) |
---|
746 | |
---|
747 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,lat) |
---|
748 | error_text="Error: could not write the latitude field in the outfile" |
---|
749 | call status_check(ierr,error_text) |
---|
750 | |
---|
751 | !============================================================================== |
---|
752 | ! ALTITUDE |
---|
753 | !============================================================================== |
---|
754 | ! Switch to netcdf define mode |
---|
755 | ierr=nf90_redef(outfid) |
---|
756 | error_text="Error: could not switch to define mode in the outfile" |
---|
757 | call status_check(ierr,error_text) |
---|
758 | |
---|
759 | ! Get the dimension in GCM file |
---|
760 | ierr=nf90_inq_dimid(gcmfid,"altitude",tmpdimid) |
---|
761 | error_text="Error: Dimension <altitude> is missing in file "//trim(gcmfile) |
---|
762 | call status_check(ierr,error_text) |
---|
763 | |
---|
764 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=altlen) |
---|
765 | allocate(alt(altlen)) |
---|
766 | |
---|
767 | ! Create the dimension in output file |
---|
768 | ierr=NF90_DEF_DIM(outfid,"altitude",altlen,altdimout) |
---|
769 | error_text="Error: could not define the altitude dimension in the outfile" |
---|
770 | call status_check(ierr,error_text) |
---|
771 | |
---|
772 | ! Get the field in GCM file |
---|
773 | ierr=nf90_inq_varid(gcmfid,"altitude",tmpvarid) |
---|
774 | error_text="Error: Field <altitude> is missing in file "//trim(gcmfile) |
---|
775 | call status_check(ierr,error_text) |
---|
776 | |
---|
777 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,alt) |
---|
778 | error_text="Failed to load altitude" |
---|
779 | call status_check(ierr,error_text) |
---|
780 | |
---|
781 | ! Create the field in the output file |
---|
782 | ierr=NF90_DEF_VAR(outfid,"altitude",nf90_float,(/altdimout/),tmpvaridout) |
---|
783 | error_text="Error: could not define the altitude variable in the outfile" |
---|
784 | call status_check(ierr,error_text) |
---|
785 | |
---|
786 | ! Get the field attributes in the GCM file |
---|
787 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
788 | if (ierr.ne.nf90_noerr) then |
---|
789 | ! if no attribute "long_name", try "title" |
---|
790 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
791 | endif |
---|
792 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
793 | ierr=nf90_get_att(gcmfid,tmpvarid,"positive",positive) |
---|
794 | |
---|
795 | ! Put the field attributes in the output file |
---|
796 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
797 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
798 | ierr=nf90_put_att(outfid,tmpvaridout,"positive",positive) |
---|
799 | |
---|
800 | ! Write the field values in the output file |
---|
801 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
802 | error_text="Error: could not end the define mode of the outfile" |
---|
803 | call status_check(ierr,error_text) |
---|
804 | |
---|
805 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,alt) |
---|
806 | error_text="Error: could not write the altitude field in the outfile" |
---|
807 | call status_check(ierr,error_text) |
---|
808 | |
---|
809 | !============================================================================== |
---|
810 | ! TIME |
---|
811 | !============================================================================== |
---|
812 | ! Switch to netcdf define mode |
---|
813 | ierr=nf90_redef(outfid) |
---|
814 | error_text="Error: could not switch to define mode in the outfile" |
---|
815 | call status_check(ierr,error_text) |
---|
816 | |
---|
817 | ! Get the dimension in GCM file |
---|
818 | ierr=nf90_inq_dimid(gcmfid,"Time",tmpdimid) |
---|
819 | error_text="Error: Dimension <Time> is missing in file "//trim(gcmfile) |
---|
820 | call status_check(ierr,error_text) |
---|
821 | |
---|
822 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=timelen) |
---|
823 | allocate(time(timelen)) |
---|
824 | |
---|
825 | ! Create the dimension in output file |
---|
826 | ierr=NF90_DEF_DIM(outfid,"Time",timelen,timedimout) |
---|
827 | error_text="Error: could not define the time dimension in the outfile" |
---|
828 | call status_check(ierr,error_text) |
---|
829 | |
---|
830 | ! Get the field in GCM file |
---|
831 | ierr=nf90_inq_varid(gcmfid,"Time",tmpvarid) |
---|
832 | error_text="Error: Field <Time> is missing in file "//trim(gcmfile) |
---|
833 | call status_check(ierr,error_text) |
---|
834 | |
---|
835 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,time) |
---|
836 | error_text="Failed to load Time" |
---|
837 | call status_check(ierr,error_text) |
---|
838 | |
---|
839 | ! Create the field in the output file |
---|
840 | ierr=NF90_DEF_VAR(outfid,"Time",nf90_float,(/timedimout/),tmpvaridout) |
---|
841 | error_text="Error: could not define the Time variable in the outfile" |
---|
842 | call status_check(ierr,error_text) |
---|
843 | |
---|
844 | ! Get the field attributes in the GCM file |
---|
845 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
846 | if (ierr.ne.nf90_noerr) then |
---|
847 | ! if no attribute "long_name", try "title" |
---|
848 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
849 | endif |
---|
850 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
851 | |
---|
852 | ! Put the field attributes in the output file |
---|
853 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
854 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
855 | |
---|
856 | ! Write the field values in the output file |
---|
857 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
858 | error_text="Error: could not end the define mode of the outfile" |
---|
859 | call status_check(ierr,error_text) |
---|
860 | |
---|
861 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,time) |
---|
862 | error_text="Error: could not write the Time field in the outfile" |
---|
863 | call status_check(ierr,error_text) |
---|
864 | |
---|
865 | !============================================================================== |
---|
866 | ! CONTROLE |
---|
867 | !============================================================================== |
---|
868 | ! Switch to netcdf define mode |
---|
869 | ierr=nf90_redef(outfid) |
---|
870 | error_text="Error: could not switch to define mode in the outfile" |
---|
871 | call status_check(ierr,error_text) |
---|
872 | |
---|
873 | ! Get the dimension in GCM file |
---|
874 | ierr=nf90_inq_dimid(gcmfid,"index",tmpdimid) |
---|
875 | error_text="Dimension <index> is missing in file "//trim(gcmfile)& |
---|
876 | //". We'll skip that one." |
---|
877 | if (ierr.ne.nf90_noerr) then |
---|
878 | write(*,*)trim(error_text) |
---|
879 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
880 | error_text="Error: could not end the define mode of the outfile" |
---|
881 | call status_check(ierr,error_text) |
---|
882 | else |
---|
883 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=ctllen) |
---|
884 | allocate(ctl(ctllen)) |
---|
885 | |
---|
886 | ! Create the dimension in output file |
---|
887 | ierr=NF90_DEF_DIM(outfid,"index",ctllen,tmpdimidout) |
---|
888 | error_text="Error: could not define the index dimension in the outfile" |
---|
889 | call status_check(ierr,error_text) |
---|
890 | |
---|
891 | ! Get the field in GCM file |
---|
892 | ierr=nf90_inq_varid(gcmfid,"controle",tmpvarid) |
---|
893 | error_text="Error: Field <controle> is missing in file "//trim(gcmfile) |
---|
894 | call status_check(ierr,error_text) |
---|
895 | |
---|
896 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,ctl) |
---|
897 | error_text="Failed to load ctl" |
---|
898 | call status_check(ierr,error_text) |
---|
899 | |
---|
900 | ! Create the field in the output file |
---|
901 | ierr=NF90_DEF_VAR(outfid,"controle",nf90_float,(/tmpdimidout/),tmpvaridout) |
---|
902 | error_text="Error: could not define the controle variable in the outfile" |
---|
903 | call status_check(ierr,error_text) |
---|
904 | |
---|
905 | ! Get the field attributes in the GCM file |
---|
906 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
907 | if (ierr.ne.nf90_noerr) then |
---|
908 | ! if no attribute "long_name", try "title" |
---|
909 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
910 | endif |
---|
911 | |
---|
912 | ! Put the field attributes in the output file |
---|
913 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
914 | |
---|
915 | ! Write the field values in the output file |
---|
916 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
917 | error_text="Error: could not end the define mode of the outfile" |
---|
918 | call status_check(ierr,error_text) |
---|
919 | |
---|
920 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,ctl) |
---|
921 | error_text="Error: could not write the controle field in the outfile" |
---|
922 | call status_check(ierr,error_text) |
---|
923 | endif |
---|
924 | |
---|
925 | |
---|
926 | !============================================================================== |
---|
927 | ! Load size of aps() or sigma() (in case it is not altlen) |
---|
928 | !============================================================================== |
---|
929 | ! Switch to netcdf define mode |
---|
930 | ierr=nf90_redef(outfid) |
---|
931 | error_text="Error: could not switch to define mode in the outfile" |
---|
932 | call status_check(ierr,error_text) |
---|
933 | |
---|
934 | ! Default is that GCM_layers=altlen |
---|
935 | ! but for outputs of zrecast, it may be a different value |
---|
936 | ierr=nf90_inq_dimid(gcmfid,"GCM_layers",tmpdimid) |
---|
937 | if (ierr.ne.nf90_noerr) then |
---|
938 | ! didn't find a GCM_layers dimension; |
---|
939 | ! we look for interlayer dimension |
---|
940 | ierr=nf90_inq_dimid(gcmfid,"interlayer",tmpdimid) |
---|
941 | if (ierr.eq.nf90_noerr) then |
---|
942 | ! load value of GCM_layers |
---|
943 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=GCM_layers) |
---|
944 | GCM_layers=GCM_layers-1 |
---|
945 | else |
---|
946 | ! didn't find a GCM_layers or interlayer dimension; therefore we have: |
---|
947 | GCM_layers=altlen |
---|
948 | endif |
---|
949 | else |
---|
950 | ! load value of GCM_layers |
---|
951 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=GCM_layers) |
---|
952 | endif |
---|
953 | |
---|
954 | ! Create the dimensions in output file |
---|
955 | ierr = NF90_DEF_DIM(outfid,"GCM_layers",GCM_layers,layerdimout) |
---|
956 | error_text="Error: could not define the GCM_layers dimension in the outfile" |
---|
957 | call status_check(ierr,error_text) |
---|
958 | ierr = NF90_DEF_DIM(outfid,"GCM_interlayers",GCM_layers+1,interlayerdimout) |
---|
959 | error_text="Error: could not define the GCM_interlayers dimension in the outfile" |
---|
960 | call status_check(ierr,error_text) |
---|
961 | |
---|
962 | ! End netcdf define mode |
---|
963 | ierr=nf90_enddef(outfid) |
---|
964 | error_text="Error: could not end the define mode of the outfile" |
---|
965 | call status_check(ierr,error_text) |
---|
966 | |
---|
967 | end subroutine inidims |
---|
968 | |
---|
969 | |
---|
970 | !******************************************************************************* |
---|
971 | |
---|
972 | subroutine init2(gcmfid,lonlen,latlen,altlen,timelen,GCM_layers, & |
---|
973 | outfid,londimout,latdimout,altdimout,timedimout, & |
---|
974 | layerdimout,interlayerdimout,& |
---|
975 | compute_col,delta_z) |
---|
976 | !============================================================================== |
---|
977 | ! Purpose: |
---|
978 | ! Copy ap() , bp(), aps(), bps(), aire() and phisinit() |
---|
979 | ! from input file to output file |
---|
980 | ! + compute layers' height if needed |
---|
981 | !============================================================================== |
---|
982 | ! Remarks: |
---|
983 | ! The NetCDF files must be open |
---|
984 | !============================================================================== |
---|
985 | |
---|
986 | use netcdf |
---|
987 | |
---|
988 | implicit none |
---|
989 | |
---|
990 | !============================================================================== |
---|
991 | ! Arguments: |
---|
992 | !============================================================================== |
---|
993 | integer, intent(in) :: gcmfid ! NetCDF output file ID |
---|
994 | integer, intent(in) :: lonlen ! # of grid points along longitude |
---|
995 | integer, intent(in) :: latlen ! # of grid points along latitude |
---|
996 | integer, intent(in) :: altlen ! # of grid points along altitude |
---|
997 | integer, intent(in) :: timelen ! # of grid points along time |
---|
998 | integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers |
---|
999 | integer, intent(in) :: outfid ! NetCDF output file ID |
---|
1000 | integer, intent(in) :: londimout ! longitude dimension ID |
---|
1001 | integer, intent(in) :: latdimout ! latitude dimension ID |
---|
1002 | integer, intent(in) :: altdimout ! altitude dimension ID |
---|
1003 | integer, intent(in) :: timedimout ! time dimension ID |
---|
1004 | integer, intent(in) :: layerdimout ! GCM_layers dimension ID |
---|
1005 | integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID |
---|
1006 | character(len=3), intent(in) :: compute_col ! does the user want to compute the column-integrated optical depth [yes/no] |
---|
1007 | |
---|
1008 | real,dimension(lonlen,latlen,altlen,timelen), intent(out) :: delta_z ! altitude difference between interlayers [km] |
---|
1009 | |
---|
1010 | !============================================================================== |
---|
1011 | ! Local variables: |
---|
1012 | !============================================================================== |
---|
1013 | real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates |
---|
1014 | real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates |
---|
1015 | real,dimension(:),allocatable :: sigma ! sigma levels |
---|
1016 | real,dimension(:,:),allocatable :: aire ! mesh areas |
---|
1017 | real,dimension(:,:),allocatable :: phisinit ! Ground geopotential |
---|
1018 | |
---|
1019 | integer :: ierr ! [netcdf] subroutine returned error code |
---|
1020 | character(len=256) :: error_text ! text to display in case of error |
---|
1021 | integer :: tmpvarid ! temporary variable ID |
---|
1022 | |
---|
1023 | logical :: area ! is "aire" available ? |
---|
1024 | logical :: phis ! is "phisinit" available ? |
---|
1025 | logical :: hybrid ! are "aps" and "bps" available ? |
---|
1026 | logical :: apbp ! are "ap" and "bp" available ? + are they usable for delta_z computation ? |
---|
1027 | |
---|
1028 | |
---|
1029 | ! for delta_z computation |
---|
1030 | real,dimension(:,:,:,:),allocatable :: zlev ! altitude of the interlayers [km] |
---|
1031 | logical :: is_zlev ! is zzlev available ? |
---|
1032 | integer,dimension(4) :: zlevdimids ! dimensions of zzlev in GCM file |
---|
1033 | integer :: zlevvertlen ! length of the vertical dimension of zlev |
---|
1034 | |
---|
1035 | real,dimension(altlen) :: alt ! altitude coordinates |
---|
1036 | character (len=64) :: altlong_name,altunits ! altitude long_name and units attributes |
---|
1037 | |
---|
1038 | real,dimension(:,:,:,:),allocatable :: plev ! pressure of the interlayers [km] |
---|
1039 | real,dimension(:,:,:),allocatable :: ps ! surface pressure [Pa] |
---|
1040 | real,dimension(:,:,:,:),allocatable :: rho ! atmospheric density [kg/m3] |
---|
1041 | real,parameter :: g0 =3.7257964 ! Mars gravity [m.s-2] |
---|
1042 | |
---|
1043 | real,dimension(:,:,:,:),allocatable :: zlay ! altitude of the layers [km] |
---|
1044 | |
---|
1045 | integer :: ilon,ilat,it,ialt ! loop indices |
---|
1046 | real :: missval ! GCM variables missing value attribute |
---|
1047 | |
---|
1048 | !============================================================================== |
---|
1049 | ! 1. Read data from input file |
---|
1050 | !============================================================================== |
---|
1051 | |
---|
1052 | ! hybrid coordinate aps |
---|
1053 | !write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
1054 | allocate(aps(GCM_layers),stat=ierr) |
---|
1055 | if (ierr.ne.0) then |
---|
1056 | write(*,*) "init2: failed to allocate aps!" |
---|
1057 | stop |
---|
1058 | endif |
---|
1059 | ierr=nf90_inq_varid(gcmfid,"aps",tmpvarid) |
---|
1060 | if (ierr.ne.nf90_noerr) then |
---|
1061 | write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord." |
---|
1062 | hybrid=.false. |
---|
1063 | else |
---|
1064 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,aps) |
---|
1065 | hybrid=.true. |
---|
1066 | if (ierr.ne.nf90_noerr) then |
---|
1067 | write(*,*) "init2 Error: Failed reading aps" |
---|
1068 | stop |
---|
1069 | endif |
---|
1070 | |
---|
1071 | ! hybrid coordinate bps |
---|
1072 | ! write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
1073 | allocate(bps(GCM_layers),stat=ierr) |
---|
1074 | if (ierr.ne.0) then |
---|
1075 | write(*,*) "init2: failed to allocate bps!" |
---|
1076 | stop |
---|
1077 | endif |
---|
1078 | ierr=nf90_inq_varid(gcmfid,"bps",tmpvarid) |
---|
1079 | if (ierr.ne.nf90_noerr) then |
---|
1080 | write(*,*) "init2 Error: Failed to get bps ID." |
---|
1081 | stop |
---|
1082 | endif |
---|
1083 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,bps) |
---|
1084 | if (ierr.ne.nf90_noerr) then |
---|
1085 | write(*,*) "init2 Error: Failed reading bps" |
---|
1086 | stop |
---|
1087 | endif |
---|
1088 | endif |
---|
1089 | |
---|
1090 | ! hybrid coordinate ap |
---|
1091 | allocate(ap(GCM_layers+1),stat=ierr) |
---|
1092 | if (ierr.ne.0) then |
---|
1093 | write(*,*) "init2: failed to allocate ap!" |
---|
1094 | stop |
---|
1095 | else |
---|
1096 | ierr=nf90_inq_varid(gcmfid,"ap",tmpvarid) |
---|
1097 | if (ierr.ne.nf90_noerr) then |
---|
1098 | write(*,*) "Ooops. Failed to get ap ID. OK." |
---|
1099 | apbp=.false. |
---|
1100 | else |
---|
1101 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,ap) |
---|
1102 | apbp=.true. |
---|
1103 | if (ierr.ne.nf90_noerr) then |
---|
1104 | write(*,*) "Error: Failed reading ap" |
---|
1105 | stop |
---|
1106 | endif |
---|
1107 | endif |
---|
1108 | endif |
---|
1109 | |
---|
1110 | ! hybrid coordinate bp |
---|
1111 | allocate(bp(GCM_layers+1),stat=ierr) |
---|
1112 | if (ierr.ne.0) then |
---|
1113 | write(*,*) "init2: failed to allocate bp!" |
---|
1114 | stop |
---|
1115 | else |
---|
1116 | ierr=nf90_inq_varid(gcmfid,"bp",tmpvarid) |
---|
1117 | if (ierr.ne.nf90_noerr) then |
---|
1118 | write(*,*) "Ooops. Failed to get bp ID. OK." |
---|
1119 | apbp=.false. |
---|
1120 | else |
---|
1121 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,bp) |
---|
1122 | apbp=.true. |
---|
1123 | if (ierr.ne.nf90_noerr) then |
---|
1124 | write(*,*) "Error: Failed reading bp" |
---|
1125 | stop |
---|
1126 | endif |
---|
1127 | endif |
---|
1128 | endif |
---|
1129 | |
---|
1130 | ! sigma levels (if any) |
---|
1131 | if (.not.hybrid) then |
---|
1132 | allocate(sigma(GCM_layers),stat=ierr) |
---|
1133 | if (ierr.ne.0) then |
---|
1134 | write(*,*) "init2: failed to allocate sigma" |
---|
1135 | stop |
---|
1136 | endif |
---|
1137 | ierr=nf90_inq_varid(gcmfid,"sigma",tmpvarid) |
---|
1138 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,sigma) |
---|
1139 | if (ierr.ne.nf90_noerr) then |
---|
1140 | write(*,*) "init2 Error: Failed reading sigma" |
---|
1141 | stop |
---|
1142 | endif |
---|
1143 | endif ! of if (.not.hybrid) |
---|
1144 | |
---|
1145 | ! mesh area |
---|
1146 | allocate(aire(lonlen,latlen),stat=ierr) |
---|
1147 | if (ierr.ne.0) then |
---|
1148 | write(*,*) "init2: failed to allocate aire!" |
---|
1149 | stop |
---|
1150 | endif |
---|
1151 | ierr=nf90_inq_varid(gcmfid,"aire",tmpvarid) |
---|
1152 | if (ierr.ne.nf90_noerr) then |
---|
1153 | write(*,*)"init2 warning: Failed to get aire ID." |
---|
1154 | area = .false. |
---|
1155 | else |
---|
1156 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,aire) |
---|
1157 | if (ierr.ne.nf90_noerr) then |
---|
1158 | write(*,*) "init2 Error: Failed reading aire" |
---|
1159 | stop |
---|
1160 | endif |
---|
1161 | area = .true. |
---|
1162 | endif |
---|
1163 | |
---|
1164 | ! ground geopotential phisinit |
---|
1165 | allocate(phisinit(lonlen,latlen),stat=ierr) |
---|
1166 | if (ierr.ne.0) then |
---|
1167 | write(*,*) "init2: failed to allocate phisinit!" |
---|
1168 | stop |
---|
1169 | endif |
---|
1170 | ierr=nf90_inq_varid(gcmfid,"phisinit",tmpvarid) |
---|
1171 | if (ierr.ne.nf90_noerr) then |
---|
1172 | write(*,*)"init2 warning: Failed to get phisinit ID." |
---|
1173 | phis = .false. |
---|
1174 | else |
---|
1175 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,phisinit) |
---|
1176 | if (ierr.ne.nf90_noerr) then |
---|
1177 | write(*,*) "init2 Error: Failed reading phisinit" |
---|
1178 | stop |
---|
1179 | endif |
---|
1180 | phis = .true. |
---|
1181 | endif |
---|
1182 | |
---|
1183 | !============================================================================== |
---|
1184 | ! 2. Write |
---|
1185 | !============================================================================== |
---|
1186 | |
---|
1187 | !========================================================================= |
---|
1188 | ! 2.1 Hybrid coordinates ap(), bp(), aps() and bps() |
---|
1189 | !========================================================================= |
---|
1190 | if(hybrid) then |
---|
1191 | ! define aps |
---|
1192 | ! Switch to netcdf define mode |
---|
1193 | ierr=nf90_redef(outfid) |
---|
1194 | ! Insert the definition of the variable |
---|
1195 | ierr=NF90_DEF_VAR(outfid,"aps",nf90_float,(/layerdimout/),tmpvarid) |
---|
1196 | if (ierr.ne.nf90_noerr) then |
---|
1197 | write(*,*) "init2 Error: Failed to define the variable aps" |
---|
1198 | stop |
---|
1199 | endif |
---|
1200 | ! Write the attributes |
---|
1201 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid pressure at midlayers") |
---|
1202 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1203 | ! End netcdf define mode |
---|
1204 | ierr=nf90_enddef(outfid) |
---|
1205 | |
---|
1206 | ! write aps |
---|
1207 | ierr=NF90_PUT_VAR(outfid,tmpvarid,aps) |
---|
1208 | if (ierr.ne.nf90_noerr) then |
---|
1209 | write(*,*) "init2 Error: Failed to write aps" |
---|
1210 | stop |
---|
1211 | endif |
---|
1212 | |
---|
1213 | ! define bps |
---|
1214 | ! Switch to netcdf define mode |
---|
1215 | ierr=nf90_redef(outfid) |
---|
1216 | ! Insert the definition of the variable |
---|
1217 | ierr=NF90_DEF_VAR(outfid,"bps",nf90_float,(/layerdimout/),tmpvarid) |
---|
1218 | if (ierr.ne.nf90_noerr) then |
---|
1219 | write(*,*) "init2 Error: Failed to define the variable bps" |
---|
1220 | stop |
---|
1221 | endif |
---|
1222 | ! Write the attributes |
---|
1223 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at midlayers") |
---|
1224 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1225 | ! End netcdf define mode |
---|
1226 | ierr=nf90_enddef(outfid) |
---|
1227 | |
---|
1228 | ! write bps |
---|
1229 | ierr=NF90_PUT_VAR(outfid,tmpvarid,bps) |
---|
1230 | if (ierr.ne.nf90_noerr) then |
---|
1231 | write(*,*) "init2 Error: Failed to write bps" |
---|
1232 | stop |
---|
1233 | endif |
---|
1234 | |
---|
1235 | if (apbp) then |
---|
1236 | ! define ap |
---|
1237 | |
---|
1238 | ! Switch to netcdf define mode |
---|
1239 | ierr=nf90_redef(outfid) |
---|
1240 | ! Insert the definition of the variable |
---|
1241 | ierr=NF90_DEF_VAR(outfid,"ap",nf90_float,(/interlayerdimout/),tmpvarid) |
---|
1242 | if (ierr.ne.nf90_noerr) then |
---|
1243 | write(*,*) "init2 Error: Failed to define the variable ap" |
---|
1244 | stop |
---|
1245 | endif |
---|
1246 | ! Write the attributes |
---|
1247 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at interlayers") |
---|
1248 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1249 | ! End netcdf define mode |
---|
1250 | ierr=nf90_enddef(outfid) |
---|
1251 | |
---|
1252 | ! write ap |
---|
1253 | ierr=NF90_PUT_VAR(outfid,tmpvarid,ap) |
---|
1254 | if (ierr.ne.nf90_noerr) then |
---|
1255 | write(*,*) "Error: Failed to write ap" |
---|
1256 | stop |
---|
1257 | endif |
---|
1258 | |
---|
1259 | ! define bp |
---|
1260 | |
---|
1261 | ! Switch to netcdf define mode |
---|
1262 | ierr=nf90_redef(outfid) |
---|
1263 | ! Insert the definition of the variable |
---|
1264 | ierr=NF90_DEF_VAR(outfid,"bp",nf90_float,(/interlayerdimout/),tmpvarid) |
---|
1265 | if (ierr.ne.nf90_noerr) then |
---|
1266 | write(*,*) "init2 Error: Failed to define the variable bp" |
---|
1267 | stop |
---|
1268 | endif |
---|
1269 | ! Write the attributes |
---|
1270 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at interlayers") |
---|
1271 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1272 | ! End netcdf define mode |
---|
1273 | ierr=nf90_enddef(outfid) |
---|
1274 | |
---|
1275 | ! write bp |
---|
1276 | ierr=NF90_PUT_VAR(outfid,tmpvarid,bp) |
---|
1277 | if (ierr.ne.nf90_noerr) then |
---|
1278 | write(*,*) "Error: Failed to write bp" |
---|
1279 | stop |
---|
1280 | endif |
---|
1281 | endif ! of if (apbp) |
---|
1282 | |
---|
1283 | else |
---|
1284 | |
---|
1285 | ! Switch to netcdf define mode |
---|
1286 | ierr=nf90_redef(outfid) |
---|
1287 | ! Insert the definition of the variable |
---|
1288 | ierr=NF90_DEF_VAR(outfid,"sigma",nf90_float,(/layerdimout/),tmpvarid) |
---|
1289 | if (ierr.ne.nf90_noerr) then |
---|
1290 | write(*,*) "init2 Error: Failed to define the variable sigma" |
---|
1291 | stop |
---|
1292 | endif |
---|
1293 | ! Write the attributes |
---|
1294 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","sigma at midlayers") |
---|
1295 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1296 | ! End netcdf define mode |
---|
1297 | ierr=nf90_enddef(outfid) |
---|
1298 | |
---|
1299 | ! write sigma |
---|
1300 | ierr=NF90_PUT_VAR(outfid,tmpvarid,sigma) |
---|
1301 | if (ierr.ne.nf90_noerr) then |
---|
1302 | write(*,*) "init2 Error: Failed to write sigma" |
---|
1303 | stop |
---|
1304 | endif |
---|
1305 | endif ! of if (hybrid) |
---|
1306 | |
---|
1307 | !========================================================================= |
---|
1308 | ! 2.2 aire() and phisinit() |
---|
1309 | !========================================================================= |
---|
1310 | |
---|
1311 | if (area) then |
---|
1312 | |
---|
1313 | ! Switch to netcdf define mode |
---|
1314 | ierr=nf90_redef(outfid) |
---|
1315 | ! Insert the definition of the variable |
---|
1316 | ierr=NF90_DEF_VAR(outfid,"aire",nf90_float,(/londimout,latdimout/),tmpvarid) |
---|
1317 | if (ierr.ne.nf90_noerr) then |
---|
1318 | write(*,*) "init2 Error: Failed to define the variable aire" |
---|
1319 | stop |
---|
1320 | endif |
---|
1321 | ! Write the attributes |
---|
1322 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","Mesh area") |
---|
1323 | ierr=nf90_put_att(outfid,tmpvarid,"units","m2") |
---|
1324 | ! End netcdf define mode |
---|
1325 | ierr=nf90_enddef(outfid) |
---|
1326 | |
---|
1327 | ! write aire |
---|
1328 | ierr=NF90_PUT_VAR(outfid,tmpvarid,aire) |
---|
1329 | if (ierr.ne.nf90_noerr) then |
---|
1330 | write(*,*) "init2 Error: Failed to write aire" |
---|
1331 | stop |
---|
1332 | endif |
---|
1333 | endif ! of if (area) |
---|
1334 | |
---|
1335 | if (phis) then |
---|
1336 | |
---|
1337 | ! Switch to netcdf define mode |
---|
1338 | ierr=nf90_redef(outfid) |
---|
1339 | ! Insert the definition of the variable |
---|
1340 | ierr=NF90_DEF_VAR(outfid,"phisinit",nf90_float,(/londimout,latdimout/),tmpvarid) |
---|
1341 | if (ierr.ne.nf90_noerr) then |
---|
1342 | write(*,*) "init2 Error: Failed to define the variable phisinit" |
---|
1343 | stop |
---|
1344 | endif |
---|
1345 | ! Write the attributes |
---|
1346 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","Ground level geopotential") |
---|
1347 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1348 | ! End netcdf define mode |
---|
1349 | ierr=nf90_enddef(outfid) |
---|
1350 | |
---|
1351 | ! write phisinit |
---|
1352 | ierr=NF90_PUT_VAR(outfid,tmpvarid,phisinit) |
---|
1353 | if (ierr.ne.nf90_noerr) then |
---|
1354 | write(*,*) "init2 Error: Failed to write phisinit" |
---|
1355 | stop |
---|
1356 | endif |
---|
1357 | |
---|
1358 | endif ! of if (phis) |
---|
1359 | |
---|
1360 | |
---|
1361 | !============================================================================== |
---|
1362 | ! 3. Compute delta_z for column integration |
---|
1363 | !============================================================================== |
---|
1364 | if (trim(compute_col).eq."yes") then |
---|
1365 | write(*,*) "Computing the layers' height delta_z..." |
---|
1366 | |
---|
1367 | !========================================================================= |
---|
1368 | ! 3.1 CASE 1: the GCM file contains zzlev variable |
---|
1369 | !========================================================================= |
---|
1370 | ! Does the field exist in GCM file? |
---|
1371 | is_zlev=.false. |
---|
1372 | ierr=nf90_inq_varid(gcmfid,"zzlev",tmpvarid) |
---|
1373 | |
---|
1374 | IF (ierr.eq.nf90_noerr) THEN |
---|
1375 | ! zzlev is present in the GCM file |
---|
1376 | is_zlev=.true. |
---|
1377 | write(*,*) "zzlev was found in GCM file" |
---|
1378 | |
---|
1379 | ! Get the field's dimensions |
---|
1380 | ierr=nf90_inquire_variable(gcmfid,tmpvarid,dimids=zlevdimids) |
---|
1381 | |
---|
1382 | ! Get the length of the vertical dimension (assuming vertical dim is #3 like other variables) |
---|
1383 | ierr=nf90_inquire_dimension(gcmfid,zlevdimids(3),len=zlevvertlen) |
---|
1384 | |
---|
1385 | ! Get missing value |
---|
1386 | ierr=nf90_get_att(gcmfid,tmpvarid,"missing_value",missval) |
---|
1387 | if (ierr.ne.nf90_noerr) then |
---|
1388 | missval = 1.e+20 |
---|
1389 | endif |
---|
1390 | |
---|
1391 | if (zlevvertlen.eq.altlen+1) then |
---|
1392 | ! we have every interlayer altitudes, even the top one |
---|
1393 | allocate(zlev(lonlen,latlen,altlen+1,timelen)) |
---|
1394 | |
---|
1395 | ! Get zlev variable |
---|
1396 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,zlev) |
---|
1397 | error_text="Failed to load zzlev" |
---|
1398 | call status_check(ierr,error_text) |
---|
1399 | |
---|
1400 | ! Compute delta_z [km] |
---|
1401 | ! NB : the unit for zzlev variable is "m" |
---|
1402 | do ilon=1,lonlen |
---|
1403 | do ilat=1,latlen |
---|
1404 | do it=1,timelen |
---|
1405 | do ialt=1,altlen |
---|
1406 | |
---|
1407 | if ((zlev(ilon,ilat,ialt+1,it).eq.missval).or.(zlev(ilon,ilat,ialt,it).eq.missval)) then |
---|
1408 | delta_z(ilon,ilat,ialt,it) = 0 |
---|
1409 | else |
---|
1410 | delta_z(ilon,ilat,ialt,it) = (zlev(ilon,ilat,ialt+1,it)-zlev(ilon,ilat,ialt,it)) /1000 |
---|
1411 | endif ! missval |
---|
1412 | |
---|
1413 | enddo ! ialt=1,altlen |
---|
1414 | enddo ! it=1,timelen |
---|
1415 | enddo ! ilat=1,latlen |
---|
1416 | enddo ! ilon=1,lonlen |
---|
1417 | |
---|
1418 | else if (zlevvertlen.eq.altlen) then |
---|
1419 | ! we have every interlayer altitudes, except the top one |
---|
1420 | allocate(zlev(lonlen,latlen,altlen,timelen)) |
---|
1421 | |
---|
1422 | ! Get zlev variable |
---|
1423 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,zlev) |
---|
1424 | error_text="Failed to load zzlev" |
---|
1425 | call status_check(ierr,error_text) |
---|
1426 | |
---|
1427 | ! Compute delta_z [km] |
---|
1428 | ! NB : the unit for zrecasted altitude is "m" |
---|
1429 | do ilon=1,lonlen |
---|
1430 | do ilat=1,latlen |
---|
1431 | do it=1,timelen |
---|
1432 | do ialt=1,altlen-1 |
---|
1433 | |
---|
1434 | if ((zlev(ilon,ilat,ialt+1,it).eq.missval).or.(zlev(ilon,ilat,ialt,it).eq.missval)) then |
---|
1435 | delta_z(ilon,ilat,ialt,it) = 0 |
---|
1436 | else |
---|
1437 | delta_z(ilon,ilat,ialt,it) = (zlev(ilon,ilat,ialt+1,it)-zlev(ilon,ilat,ialt,it)) /1000 |
---|
1438 | endif ! missval |
---|
1439 | |
---|
1440 | enddo ! ialt=1,altlen-1 |
---|
1441 | |
---|
1442 | ! top of last layer : we assume the same height than the layer below |
---|
1443 | delta_z(ilon,ilat,altlen,it) = delta_z(ilon,ilat,altlen-1,it) |
---|
1444 | |
---|
1445 | enddo ! it=1,timelen |
---|
1446 | enddo ! ilat=1,latlen |
---|
1447 | enddo ! ilon=1,lonlen |
---|
1448 | |
---|
1449 | |
---|
1450 | write(*,*) "delta_z has been well computed from variable zzlev" |
---|
1451 | write(*,*)"" |
---|
1452 | |
---|
1453 | else |
---|
1454 | ! weird case where the GCM file would have been zrecasted but not the zzlev variable? |
---|
1455 | ! anyway, we can't use zzlev then |
---|
1456 | write(*,*) "zzlev and altitude of GCM file don't have the same dimension." |
---|
1457 | write(*,*) "We will try to use the altitude dimension to compute the layers' height." |
---|
1458 | write(*,*)"" |
---|
1459 | is_zlev=.false. |
---|
1460 | |
---|
1461 | endif ! (zlevvertlen.eq.altlen+1) |
---|
1462 | |
---|
1463 | ENDIF ! (ierr.eq.nf90_noerr) = zzlev present in file |
---|
1464 | |
---|
1465 | !========================================================================= |
---|
1466 | ! 3.2 CASE 2: try to use the GCM altitude dimension |
---|
1467 | !========================================================================= |
---|
1468 | IF (is_zlev.eq..false.) THEN |
---|
1469 | ! Get the field in GCM file |
---|
1470 | ierr=nf90_inq_varid(gcmfid,"altitude",tmpvarid) |
---|
1471 | |
---|
1472 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,alt) |
---|
1473 | error_text="Failed to load altitude" |
---|
1474 | call status_check(ierr,error_text) |
---|
1475 | |
---|
1476 | ! Get the field attributes in the GCM file |
---|
1477 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",altlong_name) |
---|
1478 | if (ierr.ne.nf90_noerr) then |
---|
1479 | ! if no attribute "long_name", try "title" |
---|
1480 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",altlong_name) |
---|
1481 | endif |
---|
1482 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",altunits) |
---|
1483 | |
---|
1484 | IF (trim(altlong_name).eq."pseudo-alt") THEN |
---|
1485 | ! GCM file has not been zrecasted, so we can use ap,bp |
---|
1486 | ! and compute delta_z from the layer mass |
---|
1487 | ! we need ps and rho to use this method |
---|
1488 | |
---|
1489 | ! Check that ap,bp have the same dimension length than altlen+1 |
---|
1490 | if (GCM_layers+1.ne.altlen+1) then |
---|
1491 | apbp=.false. |
---|
1492 | endif |
---|
1493 | |
---|
1494 | ! Load Psurf to reconstruct pressure |
---|
1495 | ! Does the field exist in GCM file? |
---|
1496 | ierr=nf90_inq_varid(gcmfid,"ps",tmpvarid) |
---|
1497 | if (ierr.eq.nf90_noerr) then |
---|
1498 | allocate(ps(lonlen,latlen,timelen)) |
---|
1499 | ! Get ps variable |
---|
1500 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,ps) |
---|
1501 | error_text="Failed to load surface pressure ps" |
---|
1502 | call status_check(ierr,error_text) |
---|
1503 | else |
---|
1504 | apbp=.false. |
---|
1505 | endif |
---|
1506 | |
---|
1507 | ! Load rho to compute delta_z from the layer mass |
---|
1508 | ! Does the field exist in GCM file? |
---|
1509 | ierr=nf90_inq_varid(gcmfid,"rho",tmpvarid) |
---|
1510 | if (ierr.eq.nf90_noerr) then |
---|
1511 | allocate(rho(lonlen,latlen,altlen,timelen)) |
---|
1512 | ! Get temp variable |
---|
1513 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,rho) |
---|
1514 | error_text="Failed to load atmospheric rho" |
---|
1515 | call status_check(ierr,error_text) |
---|
1516 | else |
---|
1517 | apbp=.false. |
---|
1518 | endif |
---|
1519 | |
---|
1520 | if (apbp) then |
---|
1521 | allocate(plev(lonlen,latlen,altlen+1,timelen)) |
---|
1522 | ! Compute delta_z [km] |
---|
1523 | do ilon=1,lonlen |
---|
1524 | do ilat=1,latlen |
---|
1525 | do it=1,timelen |
---|
1526 | ! plev at first level |
---|
1527 | plev(ilon,ilat,1,it) = ap(1)+bp(1)*ps(ilon,ilat,it) |
---|
1528 | |
---|
1529 | do ialt=1,altlen |
---|
1530 | ! plev just above |
---|
1531 | plev(ilon,ilat,ialt+1,it) = ap(ialt+1)+bp(ialt+1)*ps(ilon,ilat,it) |
---|
1532 | ! delta_z [km] from layer mass |
---|
1533 | delta_z(ilon,ilat,ialt,it) = (plev(ilon,ilat,ialt,it)-plev(ilon,ilat,ialt+1,it)) & |
---|
1534 | / (g0*rho(ilon,ilat,ialt,it)) /1000 |
---|
1535 | |
---|
1536 | enddo ! ialt=1,altlen |
---|
1537 | |
---|
1538 | enddo ! it=1,timelen |
---|
1539 | enddo ! ilat=1,latlen |
---|
1540 | enddo ! ilon=1,lonlen |
---|
1541 | |
---|
1542 | write(*,*) "delta_z has been well computed from variables ap,bp,ps,rho" |
---|
1543 | write(*,*)"" |
---|
1544 | |
---|
1545 | else |
---|
1546 | ! no ap,bp in GCM file |
---|
1547 | ! => take the middle of the layer altitudes as the interlayer altitude |
---|
1548 | ! NB : the unit for "pseudo-alt" altitude (not zrecasted) is "km" |
---|
1549 | do ialt=2,altlen-1 |
---|
1550 | delta_z(:,:,ialt,:) = (alt(ialt+1)-alt(ialt-1))/2 |
---|
1551 | enddo |
---|
1552 | delta_z(:,:,1,:) = alt(2)/2 ! bottom at 0km (ok since pseudo-alt can't be negative) |
---|
1553 | delta_z(:,:,altlen,:) = delta_z(:,:,altlen-1,:) ! top of last layer : we assume the same height than the layer below |
---|
1554 | |
---|
1555 | write(*,*) "delta_z has been well computed from dimension pseudo-altitude" |
---|
1556 | write(*,*)"" |
---|
1557 | |
---|
1558 | endif ! apbp |
---|
1559 | |
---|
1560 | ELSE IF (trim(altunits).eq."m") THEN |
---|
1561 | ! GCM file has been zrecasted in altitude above surface/areoid/center of planet |
---|
1562 | ! NB : the unit for zrecasted altitude is "m" |
---|
1563 | |
---|
1564 | ! Take the arithmetic mean of the layer altitudes as the interlayer altitude |
---|
1565 | do ialt=2,altlen-1 |
---|
1566 | delta_z(:,:,ialt,:) = (alt(ialt+1)-alt(ialt-1))/2 /1000 |
---|
1567 | enddo |
---|
1568 | |
---|
1569 | ! Bottom layer : we can't assume a bottom of the first layer at 0km |
---|
1570 | ! since "altitude above areoid" can be negative |
---|
1571 | ! We thus assume zlev(2)-zlev(1) = zlay(2)-zlay(1) |
---|
1572 | delta_z(:,:,1,:) = (alt(2)-alt(1)) /1000 |
---|
1573 | |
---|
1574 | ! Top of last layer : we assume the same height than the layer below |
---|
1575 | delta_z(:,:,altlen,:) = delta_z(:,:,altlen-1,:) |
---|
1576 | |
---|
1577 | write(*,*) "delta_z has been well computed from dimension altitude" |
---|
1578 | write(*,*)"" |
---|
1579 | |
---|
1580 | ELSE IF (trim(altunits).eq."Pa") THEN |
---|
1581 | ! GCM file has been zrecasted in pressure |
---|
1582 | |
---|
1583 | ! Check if "zareoid" exist in GCM file |
---|
1584 | ierr=nf90_inq_varid(gcmfid,"zareoid",tmpvarid) |
---|
1585 | if (ierr.eq.nf90_noerr) then |
---|
1586 | allocate(zlay(lonlen,latlen,altlen,timelen)) |
---|
1587 | ! Get zareoid variable |
---|
1588 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,zlay) |
---|
1589 | error_text="Failed to load zareoid" |
---|
1590 | call status_check(ierr,error_text) |
---|
1591 | |
---|
1592 | ! Get missing value |
---|
1593 | ierr=nf90_get_att(gcmfid,tmpvarid,"missing_value",missval) |
---|
1594 | if (ierr.ne.nf90_noerr) then |
---|
1595 | missval = 1.e+20 |
---|
1596 | endif |
---|
1597 | |
---|
1598 | ! Take the arithmetic mean of the layer altitudes as the interlayer altitude |
---|
1599 | ! NB : the unit for zareoid is "m" |
---|
1600 | ! Compute delta_z [km] |
---|
1601 | do ilon=1,lonlen |
---|
1602 | do ilat=1,latlen |
---|
1603 | do it=1,timelen |
---|
1604 | |
---|
1605 | ! Bottom layer |
---|
1606 | if ((zlay(ilon,ilat,2,it).eq.missval).or.(zlay(ilon,ilat,1,it).eq.missval)) then |
---|
1607 | delta_z(ilon,ilat,1,it) = 0 |
---|
1608 | else |
---|
1609 | ! We can't assume a bottom of the first layer at 0km |
---|
1610 | ! since "altitude above areoid" can be negative |
---|
1611 | ! We thus assume zlev(2)-zlev(1) = zlay(2)-zlay(1) |
---|
1612 | delta_z(:,:,1,:) = (zlay(ilon,ilat,2,it)-zlay(ilon,ilat,1,it)) /1000 |
---|
1613 | endif ! missval |
---|
1614 | |
---|
1615 | ! rest of the column |
---|
1616 | do ialt=2,altlen-1 |
---|
1617 | if ((zlay(ilon,ilat,ialt+1,it).eq.missval).or.(zlay(ilon,ilat,ialt-1,it).eq.missval)) then |
---|
1618 | delta_z(ilon,ilat,ialt,it) = 0 |
---|
1619 | else |
---|
1620 | delta_z(ilon,ilat,ialt,it) = (zlay(ilon,ilat,ialt+1,it)-zlay(ilon,ilat,ialt-1,it))/2 /1000 |
---|
1621 | endif ! missval |
---|
1622 | enddo ! ialt=2,altlen-1 |
---|
1623 | |
---|
1624 | ! Top of last layer : we assume the same height than the layer below |
---|
1625 | delta_z(ilon,ilat,altlen,it) = delta_z(ilon,ilat,altlen-1,it) |
---|
1626 | |
---|
1627 | enddo ! it=1,timelen |
---|
1628 | enddo ! ilat=1,latlen |
---|
1629 | enddo ! ilon=1,lonlen |
---|
1630 | |
---|
1631 | write(*,*) "delta_z has been well computed from variable zareoid" |
---|
1632 | write(*,*)"" |
---|
1633 | |
---|
1634 | else |
---|
1635 | |
---|
1636 | write(*,*) "No variable zareoid present in GCM file" |
---|
1637 | write(*,*) "Computing the layers' height with the layers' pressure coordinate" |
---|
1638 | write(*,*) "is not (yet?) handled by the program." |
---|
1639 | write(*,*) "The column-integrated optical depths can thus not be computed." |
---|
1640 | write(*,*) "Better stop now..." |
---|
1641 | write(*,*) "You may retry with putting 'no' for the column computation." |
---|
1642 | stop |
---|
1643 | |
---|
1644 | endif ! (ierr.eq.nf90_noerr) = zareoid present in file |
---|
1645 | |
---|
1646 | |
---|
1647 | ELSE |
---|
1648 | ! abort the program |
---|
1649 | write(*,*) "delta_z can not be computed from the GCM file's variables." |
---|
1650 | write(*,*) "The column-integrated optical depths can thus not be computed." |
---|
1651 | write(*,*) "Better stop now..." |
---|
1652 | write(*,*) "You may retry with putting 'no' for the column computation." |
---|
1653 | stop |
---|
1654 | |
---|
1655 | ENDIF ! trim(altlong_name).eq."pseudo-alt" |
---|
1656 | |
---|
1657 | ENDIF ! is_zlev.eq.false |
---|
1658 | |
---|
1659 | !========================================================================= |
---|
1660 | ! 3.3 Write delta_z in output file |
---|
1661 | !========================================================================= |
---|
1662 | |
---|
1663 | ! Switch to netcdf define mode |
---|
1664 | ierr=nf90_redef(outfid) |
---|
1665 | ! Insert the definition of the variable |
---|
1666 | ierr=NF90_DEF_VAR(outfid,"delta_z",nf90_float,(/londimout,latdimout,altdimout,timedimout/),tmpvarid) |
---|
1667 | if (ierr.ne.nf90_noerr) then |
---|
1668 | write(*,*) "init2 Error: Failed to define the variable delta_z" |
---|
1669 | stop |
---|
1670 | endif |
---|
1671 | ! Write the attributes |
---|
1672 | ierr=nf90_put_att(outfid,tmpvarid,"long_name","Layer height used for column optical depth computation") |
---|
1673 | ierr=nf90_put_att(outfid,tmpvarid,"units","km") |
---|
1674 | ! End netcdf define mode |
---|
1675 | ierr=nf90_enddef(outfid) |
---|
1676 | |
---|
1677 | ! write delta_z |
---|
1678 | ierr=NF90_PUT_VAR(outfid,tmpvarid,delta_z) |
---|
1679 | if (ierr.ne.nf90_noerr) then |
---|
1680 | write(*,*) "init2 Error: Failed to write delta_z" |
---|
1681 | stop |
---|
1682 | endif |
---|
1683 | |
---|
1684 | endif ! trim(compute_col).eq."yes" |
---|
1685 | |
---|
1686 | !============================================================================== |
---|
1687 | ! 4. Cleanup |
---|
1688 | !============================================================================== |
---|
1689 | |
---|
1690 | if (allocated(aps)) deallocate(aps) |
---|
1691 | if (allocated(bps)) deallocate(bps) |
---|
1692 | if (allocated(ap)) deallocate(ap) |
---|
1693 | if (allocated(bp)) deallocate(bp) |
---|
1694 | if (allocated(sigma)) deallocate(sigma) |
---|
1695 | if (allocated(phisinit)) deallocate(phisinit) |
---|
1696 | if (allocated(aire)) deallocate(aire) |
---|
1697 | |
---|
1698 | end subroutine init2 |
---|
1699 | |
---|
1700 | |
---|
1701 | !******************************************************************************* |
---|