1 | program aeroptical |
---|
2 | ! program for computation of aerosol opacities |
---|
3 | ! author : Antoine Bierjon, April 2020 |
---|
4 | ! |
---|
5 | !=============================================================================== |
---|
6 | ! PREFACE |
---|
7 | !=============================================================================== |
---|
8 | ! This program takes in input a GCM output file (diagfi,stats,concat) that |
---|
9 | ! contains: |
---|
10 | ! - the Mass Mixing Ratios of dust (dustq) and water ice (h2o_ice) |
---|
11 | ! - their effective radius (reffdust, reffice(*)) |
---|
12 | ! - the atmospheric density (rho) |
---|
13 | ! as well as the 2 boundaries of a wavelength interval and the directory |
---|
14 | ! containing the ASCII files of the optical properties of the aerosols. |
---|
15 | ! |
---|
16 | ! It outputs a netcdf file containing the opacities [1/km] of the dust and |
---|
17 | ! water ice aerosols as well as the combined opacity of the 2 aerosols, |
---|
18 | ! averaged within the prescribed wavelength interval. |
---|
19 | ! |
---|
20 | ! |
---|
21 | ! (*) due to a high uncertainty of reffice in the gcm, the value is asked |
---|
22 | ! directly to the user (if the user returns 0, then the program reads the GCM |
---|
23 | ! file to get reffice) |
---|
24 | ! |
---|
25 | ! NOTA BENE: |
---|
26 | ! 1) if one wanted to add another aerosol to compute, one should look for |
---|
27 | ! the comments + NEW AER? that are disseminated all along the code to indicate |
---|
28 | ! the parts of the code that must be modified. |
---|
29 | ! 2) another enhancement of this program could be the possibility to read its |
---|
30 | ! own product files and recalibrate the opacities at another wavelength |
---|
31 | !=============================================================================== |
---|
32 | |
---|
33 | |
---|
34 | use netcdf |
---|
35 | |
---|
36 | !=============================================================================== |
---|
37 | ! Variable declarations |
---|
38 | !=============================================================================== |
---|
39 | |
---|
40 | implicit none ! for no implicitly typed variables |
---|
41 | |
---|
42 | |
---|
43 | ! GCM file |
---|
44 | character(len=128) :: gcmfile ! name of the netcdf GCM input file |
---|
45 | integer :: gcmfid ! [netcdf] GCM input file ID number |
---|
46 | integer :: ierr ! [netcdf] subroutine returned error code |
---|
47 | character(len=256) :: error_text ! text to display in case of error |
---|
48 | integer :: tmpvarid ! temporary store a variable ID number |
---|
49 | integer :: lonlen,latlen,altlen,timelen ! nb of grid points along longitude,latitude,altitude,time |
---|
50 | integer :: GCM_layers ! number of GCM layers |
---|
51 | integer :: layerdimout,interlayerdimout ! "GCM_layers" and "GCM_layers+1" IDs |
---|
52 | |
---|
53 | logical,dimension(:),allocatable :: aerok ! to know if the needed fields are in the GCM file |
---|
54 | real,dimension(:,:,:,:,:),allocatable :: mmr ! aerosols mass mixing ratios [kg/kg] |
---|
55 | real,dimension(:,:,:,:,:),allocatable :: reff ! aerosols effective radii [m] |
---|
56 | real :: reffwice_val ! value given by the user for reffwice (0 if read in the GCM file) [m] |
---|
57 | real,dimension(:,:,:,:),allocatable :: rho ! atmospheric density [kg.m-3] |
---|
58 | integer :: naerkind ! nb of aerosols to consider |
---|
59 | integer :: iaer ! aerosol kind index (1=dust, 2=water ice) |
---|
60 | integer :: ilon,ilat,ialt,it ! Loop indices for coordinates |
---|
61 | |
---|
62 | |
---|
63 | ! Output file |
---|
64 | character(len=128) :: outfile ! name of the netcdf output file |
---|
65 | integer :: outfid ! [netcdf] file ID number |
---|
66 | integer :: latdimout,londimout,altdimout,timedimout |
---|
67 | ! latdimout: stores latitude dimension ID number |
---|
68 | ! londimout: stores longitude dimension ID number |
---|
69 | ! altdimout: stores altitude dimension ID number |
---|
70 | ! timedimout: stores time dimension ID number |
---|
71 | integer :: tmpvaridout ! temporary stores a variable ID number |
---|
72 | |
---|
73 | real :: wvl_val ! reference wavelength of the output opacities (given by the user) [m] |
---|
74 | integer :: varshape(4) ! stores a variable shape (of dimensions' IDs) |
---|
75 | character(len=16) :: tmpvarname ! temporary stores a variable name |
---|
76 | real,dimension(:,:,:,:,:),allocatable :: opa_aer ! Opacity of the aerosols [1/km] |
---|
77 | real :: missval ! Value to put in outfile when the reff is out of bounds |
---|
78 | PARAMETER(missval=1E+20) |
---|
79 | |
---|
80 | |
---|
81 | ! Optical properties (read in external ASCII files) |
---|
82 | character(len=256) :: datadir ! directory containing the ASCII files |
---|
83 | character(len=128) :: optpropfile ! name of the ASCII optical properties file |
---|
84 | logical :: file_ok ! to know if the file can be opened |
---|
85 | integer :: file_unit = 60 ! |
---|
86 | |
---|
87 | integer :: jfile ! ASCII file scan index |
---|
88 | logical :: endwhile ! Reading loop boolean |
---|
89 | character(len=132) :: scanline ! ASCII file scanning line |
---|
90 | integer :: read_ok ! to know if the line is well read |
---|
91 | |
---|
92 | integer :: nwvl ! Number of wavelengths in the domain (VIS or IR) |
---|
93 | integer :: nsize ! Number of particle sizes stored in the file |
---|
94 | integer :: iwvl,isize ! Wavelength and Particle size loop indices |
---|
95 | real,dimension(:),allocatable :: wvl ! Wavelength axis [m] |
---|
96 | real,dimension(:),allocatable :: radiusdyn ! Particle size axis [m] |
---|
97 | real,dimension(:,:),allocatable :: Qext ! Extinction efficiency coefficient [/] |
---|
98 | |
---|
99 | |
---|
100 | ! Opacity computation |
---|
101 | integer :: iwvl1,iwvl2,isize1,isize2 ! Wavelength and Particle size indices for the interpolation |
---|
102 | real :: coeff ! Interpolation coeficient |
---|
103 | real,dimension(2) :: reffQext ! Qext after reff interpolation |
---|
104 | real :: Qext_val ! Qext after both interpolations |
---|
105 | real,dimension(:),allocatable :: rho_aer ! Density of the aerosols [kg.m-3] |
---|
106 | |
---|
107 | !=============================================================================== |
---|
108 | ! 0. USER INPUTS |
---|
109 | !=============================================================================== |
---|
110 | write(*,*) "Welcome in the aerosol opacities' computation program !" |
---|
111 | write(*,*) |
---|
112 | |
---|
113 | ! Ask the GCM file name |
---|
114 | WRITE(*,*) "Enter an input file name (GCM diagfi/stats/concat) :" |
---|
115 | READ(*,*) gcmfile |
---|
116 | |
---|
117 | ! Ask the reffwice to the user |
---|
118 | write(*,*)"" |
---|
119 | write(*,*) "The water ice effective radius computed by the GCM is known to be a bit inaccurate." |
---|
120 | write(*,*) "Which value do you want to use for it (in meters) ?" |
---|
121 | write(*,*) "(put 0 if you still want the program to read the value in "//trim(gcmfile)//")" |
---|
122 | READ(*,*) reffwice_val |
---|
123 | |
---|
124 | ! Ask the wavelength to the user |
---|
125 | write(*,*)"" |
---|
126 | WRITE(*,*) "Value of the reference wavelength for the opacities' computation (in meters) ?" |
---|
127 | READ(*,*) wvl_val |
---|
128 | |
---|
129 | ! Ask the directory containing the optical properties files |
---|
130 | write(*,*)"" |
---|
131 | WRITE(*,*) "In which directory should we look for the optical properties' files ?" |
---|
132 | READ(*,'(a)') datadir |
---|
133 | |
---|
134 | !=============================================================================== |
---|
135 | ! 1. GCM FILE & OUTPUT FILE INITIALIZATIONS |
---|
136 | !=============================================================================== |
---|
137 | |
---|
138 | !========================================================================== |
---|
139 | ! 1.1 GCM file opening & dimensions - Output file initializations |
---|
140 | !========================================================================== |
---|
141 | |
---|
142 | !========================================================================== |
---|
143 | ! --> 1.1.1 Open the netcdf GCM file given by the user |
---|
144 | |
---|
145 | ierr=nf90_open(gcmfile,nf90_nowrite,gcmfid) ! nowrite mode=the program can only read the opened file |
---|
146 | error_text="Error: could not open file "//trim(gcmfile) |
---|
147 | call status_check(ierr,error_text) |
---|
148 | |
---|
149 | !========================================================================== |
---|
150 | ! --> 1.1.2 Creation of the outfile |
---|
151 | outfile=gcmfile(1:index(gcmfile, ".nc")-1)//"_OPA.nc" |
---|
152 | |
---|
153 | ierr=NF90_CREATE(outfile,nf90_clobber,outfid) ! NB: clobber mode=overwrite any dataset with the same file name ! |
---|
154 | error_text="Error: could not create file "//trim(outfile) |
---|
155 | call status_check(ierr,error_text) |
---|
156 | write(*,*)"";WRITE(*,*)"Output file is: ",trim(outfile);write(*,*)"" |
---|
157 | |
---|
158 | !========================================================================== |
---|
159 | ! --> 1.1.3 Get the dimensions and create them in the output file |
---|
160 | |
---|
161 | ! Initialize output file's lat,lon,alt,time and controle dimensions |
---|
162 | call inidims(gcmfile,gcmfid,outfile,outfid,& |
---|
163 | lonlen,latlen,altlen,timelen,& |
---|
164 | latdimout,londimout,altdimout,timedimout,& |
---|
165 | GCM_layers,layerdimout,interlayerdimout) |
---|
166 | |
---|
167 | ! Initialize output file's aps,bps,ap,bp and phisinit variables |
---|
168 | call init2(gcmfid,lonlen,latlen,altlen,GCM_layers,& |
---|
169 | outfid,londimout,latdimout,altdimout,& |
---|
170 | layerdimout,interlayerdimout) |
---|
171 | |
---|
172 | !========================================================================== |
---|
173 | ! 1.2 GCM file variables |
---|
174 | !========================================================================== |
---|
175 | |
---|
176 | ! Number of aerosols |
---|
177 | naerkind = 2 ! dust & water ice ! + NEW AER? |
---|
178 | |
---|
179 | ! Length allocation of the variables |
---|
180 | allocate(mmr(naerkind,lonlen,latlen,altlen,timelen)) |
---|
181 | allocate(reff(naerkind,lonlen,latlen,altlen,timelen)) |
---|
182 | |
---|
183 | ! Initialize aerok to .true. for all aerosols |
---|
184 | allocate(aerok(naerkind)) |
---|
185 | aerok(:)=.true. |
---|
186 | |
---|
187 | !========================================================================== |
---|
188 | ! --> 1.2.1 DUST |
---|
189 | |
---|
190 | ! DUST MASS MIXING RATIO |
---|
191 | ! Check that the GCM file contains that variable |
---|
192 | ierr=nf90_inq_varid(gcmfid,"dustq",tmpvarid) |
---|
193 | error_text="Failed to find variable dustq in "//trim(gcmfile)& |
---|
194 | //" We'll skip the dust aerosol." |
---|
195 | if (ierr.ne.nf90_noerr) then |
---|
196 | write(*,*)trim(error_text) |
---|
197 | aerok(1)=.false. |
---|
198 | else |
---|
199 | ! Load datasets |
---|
200 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,mmr(1,:,:,:,:)) |
---|
201 | error_text="Failed to load dust mass mixing ratio" |
---|
202 | call status_check(ierr,error_text) |
---|
203 | write(*,*) "Dust mass mixing ratio loaded from "//trim(gcmfile) |
---|
204 | endif |
---|
205 | |
---|
206 | ! DUST EFFECTIVE RADIUS |
---|
207 | if (aerok(1)) then |
---|
208 | ! Check that the GCM file contains that variable |
---|
209 | ierr=nf90_inq_varid(gcmfid,"reffdust",tmpvarid) |
---|
210 | error_text="Failed to find variable reffdust in "//trim(gcmfile)& |
---|
211 | //" We'll skip the dust aerosol." |
---|
212 | if (ierr.ne.nf90_noerr) then |
---|
213 | write(*,*)trim(error_text) |
---|
214 | aerok(1)=.false. |
---|
215 | else |
---|
216 | ! Load datasets |
---|
217 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,reff(1,:,:,:,:)) |
---|
218 | error_text="Failed to load dust effective radius" |
---|
219 | call status_check(ierr,error_text) |
---|
220 | write(*,*) "Dust effective radius loaded from "//trim(gcmfile) |
---|
221 | endif |
---|
222 | endif |
---|
223 | |
---|
224 | !========================================================================== |
---|
225 | ! --> 1.2.2 WATER ICE |
---|
226 | |
---|
227 | ! WATER ICE MASS MIXING RATIO |
---|
228 | ! Check that the GCM file contains that variable |
---|
229 | ierr=nf90_inq_varid(gcmfid,"h2o_ice",tmpvarid) |
---|
230 | error_text="Failed to find variable h2o_ice in "//trim(gcmfile)& |
---|
231 | //" We'll skip the water ice aerosol." |
---|
232 | if (ierr.ne.nf90_noerr) then |
---|
233 | write(*,*)trim(error_text) |
---|
234 | aerok(2)=.false. |
---|
235 | else |
---|
236 | ! Load datasets |
---|
237 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,mmr(2,:,:,:,:)) |
---|
238 | error_text="Failed to load water ice mass mixing ratio" |
---|
239 | call status_check(ierr,error_text) |
---|
240 | write(*,*) "Water ice mass mixing ratio loaded from "//trim(gcmfile) |
---|
241 | endif |
---|
242 | |
---|
243 | ! WATER ICE EFFECTIVE RADIUS |
---|
244 | if (aerok(2)) then |
---|
245 | IF (reffwice_val.eq.0) THEN |
---|
246 | ! Check that the GCM file contains that variable |
---|
247 | ierr=nf90_inq_varid(gcmfid,"reffice",tmpvarid) |
---|
248 | error_text="Failed to find variable reffice in "//trim(gcmfile)& |
---|
249 | //" We'll skip the water ice aerosol." |
---|
250 | if (ierr.ne.nf90_noerr) then |
---|
251 | write(*,*)trim(error_text) |
---|
252 | aerok(2)=.false. |
---|
253 | else |
---|
254 | ! Load datasets |
---|
255 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,reff(2,:,:,1,:)) ! reffice is actually a GCM |
---|
256 | ! surface (column-averaged) variable |
---|
257 | error_text="Failed to load water ice effective radius" |
---|
258 | call status_check(ierr,error_text) |
---|
259 | do ialt=2,altlen |
---|
260 | reff(2,:,:,ialt,:) = reff(2,:,:,1,:) ! build up along altitude axis |
---|
261 | enddo |
---|
262 | write(*,*) "Water ice effective radius loaded from "//trim(gcmfile) |
---|
263 | endif |
---|
264 | ELSE ! if reffwice_val/=0 |
---|
265 | reff(2,:,:,:,:) = reffwice_val |
---|
266 | write(*,*) "Water ice effective radius loaded from the user input" |
---|
267 | ENDIF |
---|
268 | endif |
---|
269 | |
---|
270 | ! Check if there is still at least one aerosol to compute |
---|
271 | IF (.NOT.ANY(aerok)) THEN |
---|
272 | write(*,*) "Neither dust nor water ice are found in the file. Better stop now..." |
---|
273 | stop |
---|
274 | ENDIF |
---|
275 | |
---|
276 | !========================================================================== |
---|
277 | ! --> 1.2.3 + NEW AER? |
---|
278 | |
---|
279 | !========================================================================== |
---|
280 | ! --> 1.2.3 ATMOSPHERIC DENSITY |
---|
281 | |
---|
282 | ! Check that the GCM file contains that variable |
---|
283 | ierr=nf90_inq_varid(gcmfid,"rho",tmpvarid) |
---|
284 | error_text="Failed to find variable rho in "//trim(gcmfile)& |
---|
285 | //" We need it to compute the opacity [1/km]." |
---|
286 | call status_check(ierr,error_text) |
---|
287 | ! Length allocation for each dimension of the 4D variable |
---|
288 | allocate(rho(lonlen,latlen,altlen,timelen)) |
---|
289 | ! Load datasets |
---|
290 | ierr=nf90_get_var(gcmfid,tmpvarid,rho) |
---|
291 | error_text="Failed to load atmospheric density" |
---|
292 | call status_check(ierr,error_text) |
---|
293 | write(*,*) "Loaded atmospheric density rho from "//trim(gcmfile) |
---|
294 | |
---|
295 | !========================================================================== |
---|
296 | ! 1.3 Some output variable's initializations |
---|
297 | !========================================================================== |
---|
298 | ! Define the ID shape of the output variables |
---|
299 | varshape(1)=londimout |
---|
300 | varshape(2)=latdimout |
---|
301 | varshape(3)=altdimout |
---|
302 | varshape(4)=timedimout |
---|
303 | |
---|
304 | ! Fill the aerosol density array |
---|
305 | allocate(rho_aer(naerkind)) |
---|
306 | rho_aer(1)=2500. ! dust |
---|
307 | rho_aer(2)=920. ! water ice |
---|
308 | ! + NEW AER? |
---|
309 | |
---|
310 | allocate(opa_aer(naerkind,lonlen,latlen,altlen,timelen)) |
---|
311 | |
---|
312 | |
---|
313 | !=============================================================================== |
---|
314 | ! 2. OPTICAL PROPERTIES |
---|
315 | !=============================================================================== |
---|
316 | |
---|
317 | DO iaer = 1, naerkind ! Loop on aerosol kind |
---|
318 | if (aerok(iaer)) then ! check if this aerosol can be computed |
---|
319 | |
---|
320 | !========================================================================== |
---|
321 | ! 2.1 Open the ASCII file |
---|
322 | !========================================================================== |
---|
323 | IF (wvl_val.lt.5.e-6) THEN |
---|
324 | ! "VISIBLE" DOMAIN |
---|
325 | |
---|
326 | SELECT CASE(iaer) |
---|
327 | CASE(1) ! DUST |
---|
328 | ! Dust file |
---|
329 | optpropfile = "optprop_dustvis_TM_n50.dat" |
---|
330 | CASE(2) ! WATER ICE |
---|
331 | ! Water ice file |
---|
332 | optpropfile = "optprop_icevis_n30.dat" |
---|
333 | ! + NEW AER? |
---|
334 | END SELECT ! iaer |
---|
335 | |
---|
336 | ELSE ! wvl_val.ge.5.e-6 |
---|
337 | ! "INFRARED" DOMAIN |
---|
338 | |
---|
339 | SELECT CASE(iaer) |
---|
340 | CASE(1) ! DUST |
---|
341 | ! Dust file |
---|
342 | optpropfile = "optprop_dustir_n50.dat" |
---|
343 | CASE(2) ! WATER ICE |
---|
344 | ! Water ice file |
---|
345 | optpropfile = "optprop_iceir_n30.dat" |
---|
346 | ! + NEW AER? |
---|
347 | END SELECT ! iaer |
---|
348 | |
---|
349 | ENDIF ! wvl_val |
---|
350 | |
---|
351 | INQUIRE(FILE=trim(datadir)//'/'//trim(optpropfile),EXIST=file_ok) |
---|
352 | if(.not.file_ok) then |
---|
353 | write(*,*)'Problem opening ',trim(optpropfile) |
---|
354 | write(*,*)'It should be in: ',trim(datadir) |
---|
355 | stop |
---|
356 | endif |
---|
357 | OPEN(UNIT=file_unit,FILE=trim(datadir)//'/'//trim(optpropfile),FORM='formatted') |
---|
358 | |
---|
359 | !========================================================================== |
---|
360 | ! 2.2 Allocate the optical property table |
---|
361 | !========================================================================== |
---|
362 | jfile = 1 |
---|
363 | endwhile = .false. |
---|
364 | DO WHILE (.NOT.endwhile) |
---|
365 | READ(file_unit,*,iostat=read_ok) scanline |
---|
366 | if (read_ok.ne.0) then |
---|
367 | write(*,*)'Error reading file ',& |
---|
368 | trim(datadir)//'/'//trim(optpropfile) |
---|
369 | stop |
---|
370 | endif |
---|
371 | IF ((scanline(1:1).ne.'#').and.(scanline(1:1).ne.' ')) THEN |
---|
372 | BACKSPACE(file_unit) |
---|
373 | reading1_seq: SELECT CASE (jfile) ! FIRST READING SEQUENCE |
---|
374 | CASE(1) reading1_seq ! nwvl ---------------------------- |
---|
375 | read(file_unit,*,iostat=read_ok) nwvl |
---|
376 | if (read_ok.ne.0) then |
---|
377 | write(*,*)'reading1_seq: Error while reading line: ',& |
---|
378 | trim(scanline) |
---|
379 | write(*,*)' of file ',& |
---|
380 | trim(datadir)//'/'//trim(optpropfile) |
---|
381 | stop |
---|
382 | endif |
---|
383 | jfile = jfile+1 |
---|
384 | CASE(2) reading1_seq ! nsize --------------------------- |
---|
385 | read(file_unit,*,iostat=read_ok) nsize |
---|
386 | if (read_ok.ne.0) then |
---|
387 | write(*,*)'reading1_seq: Error while reading line: ',& |
---|
388 | trim(scanline) |
---|
389 | write(*,*)' of file ',& |
---|
390 | trim(datadir)//'/'//trim(optpropfile) |
---|
391 | stop |
---|
392 | endif |
---|
393 | endwhile = .true. |
---|
394 | CASE DEFAULT reading1_seq ! default -------------------- |
---|
395 | write(*,*) 'reading1_seq: ',& |
---|
396 | 'Error while loading optical properties.' |
---|
397 | stop |
---|
398 | END SELECT reading1_seq ! ============================== |
---|
399 | ENDIF |
---|
400 | ENDDO ! DO WHILE(.NOT.endwhile) |
---|
401 | |
---|
402 | ALLOCATE(wvl(nwvl)) ! Wavelength axis |
---|
403 | ALLOCATE(radiusdyn(nsize)) ! Aerosol radius axis |
---|
404 | ALLOCATE(Qext(nwvl,nsize)) ! Extinction efficiency coeff |
---|
405 | |
---|
406 | !========================================================================== |
---|
407 | ! 2.3 Read the data |
---|
408 | !========================================================================== |
---|
409 | jfile = 1 |
---|
410 | endwhile = .false. |
---|
411 | DO WHILE (.NOT.endwhile) |
---|
412 | READ(file_unit,*) scanline |
---|
413 | IF ((scanline(1:1).ne.'#').and.(scanline(1:1).ne.' ')) THEN |
---|
414 | BACKSPACE(file_unit) |
---|
415 | reading2_seq: SELECT CASE (jfile) ! SECOND READING SEQUENCE |
---|
416 | CASE(1) reading2_seq ! wvl ----------------------------- |
---|
417 | read(file_unit,*,iostat=read_ok) wvl |
---|
418 | if (read_ok.ne.0) then |
---|
419 | write(*,*)'reading2_seq: Error while reading line: ',& |
---|
420 | trim(scanline) |
---|
421 | write(*,*)' of file ',& |
---|
422 | trim(datadir)//'/'//trim(optpropfile) |
---|
423 | stop |
---|
424 | endif |
---|
425 | jfile = jfile+1 |
---|
426 | CASE(2) reading2_seq ! radiusdyn ----------------------- |
---|
427 | read(file_unit,*,iostat=read_ok) radiusdyn |
---|
428 | if (read_ok.ne.0) then |
---|
429 | write(*,*)'reading2_seq: Error while reading line: ',& |
---|
430 | trim(scanline) |
---|
431 | write(*,*)' of file ',& |
---|
432 | trim(datadir)//'/'//trim(optpropfile) |
---|
433 | stop |
---|
434 | endif |
---|
435 | jfile = jfile+1 |
---|
436 | CASE(3) reading2_seq ! Qext ---------------------------- |
---|
437 | isize = 1 |
---|
438 | DO WHILE (isize .le. nsize) |
---|
439 | READ(file_unit,*,iostat=read_ok) scanline |
---|
440 | if (read_ok.ne.0) then |
---|
441 | write(*,*)'reading2_seq: Error while reading line: ',& |
---|
442 | trim(scanline) |
---|
443 | write(*,*)' of file ',& |
---|
444 | trim(datadir)//'/'//trim(optpropfile) |
---|
445 | stop |
---|
446 | endif |
---|
447 | IF ((scanline(1:1).ne.'#').and.(scanline(1:1).ne.' ')) THEN |
---|
448 | BACKSPACE(file_unit) |
---|
449 | read(file_unit,*) Qext(:,isize) |
---|
450 | isize = isize + 1 |
---|
451 | ENDIF |
---|
452 | ENDDO |
---|
453 | endwhile = .true. |
---|
454 | CASE DEFAULT reading2_seq ! default -------------------- |
---|
455 | write(*,*) 'reading2_seq: ',& |
---|
456 | 'Error while loading optical properties.' |
---|
457 | stop |
---|
458 | END SELECT reading2_seq ! ============================== |
---|
459 | ENDIF |
---|
460 | ENDDO |
---|
461 | |
---|
462 | ! Close the file |
---|
463 | CLOSE(file_unit) |
---|
464 | |
---|
465 | write(*,*)"" |
---|
466 | write(*,*) "Wavelength array loaded from file ",trim(datadir)//'/'//trim(optpropfile) |
---|
467 | write(*,*) "ranging from ",wvl(1)," to ",wvl(nwvl)," meters" |
---|
468 | write(*,*) "Effective radius array loaded from file ",trim(datadir)//'/'//trim(optpropfile) |
---|
469 | write(*,*) "ranging from ",radiusdyn(1)," to ",radiusdyn(nsize)," meters" |
---|
470 | |
---|
471 | ! + NEW AER? : one may adapt this part to handle the properties' file of the new aerosol |
---|
472 | |
---|
473 | !========================================================================== |
---|
474 | ! 2.4 Get the optpropfile wavelength values encompassing the ref wavelength |
---|
475 | !========================================================================== |
---|
476 | iwvl=1 |
---|
477 | DO WHILE ((iwvl.le.nwvl).and.(wvl(iwvl).lt.wvl_val)) |
---|
478 | iwvl=iwvl+1 |
---|
479 | ENDDO |
---|
480 | if ((iwvl.gt.nwvl).or.(wvl(1).gt.wvl_val)) then |
---|
481 | write(*,*) "ERROR: the reference wavelength is out of the bounds" |
---|
482 | write(*,*) "of the file ",trim(datadir)//'/'//trim(optpropfile) |
---|
483 | write(*,*) "(wvl_inf=",wvl(1),"m ; wvl_sup=",wvl(nwvl),"m)" |
---|
484 | write(*,*) "Ensure you wrote it well (in meters)," |
---|
485 | write(*,*) "or supply a new optical properties' file (change in aeroptical.F90 directly)" |
---|
486 | stop |
---|
487 | endif |
---|
488 | if (wvl(iwvl).eq.wvl_val) then |
---|
489 | ! if the ref wavelength is in the optpropfile |
---|
490 | iwvl1 = iwvl |
---|
491 | iwvl2 = iwvl |
---|
492 | else ! wvl(iwvl)>wvl_val and wvl(iwvl-1)<wvl_val |
---|
493 | iwvl1 = iwvl-1 |
---|
494 | iwvl2 = iwvl |
---|
495 | endif |
---|
496 | |
---|
497 | |
---|
498 | !=============================================================================== |
---|
499 | ! 3. OUTPUT FILE VARIABLES |
---|
500 | !=============================================================================== |
---|
501 | !========================================================================== |
---|
502 | ! 3.1 Creation of the output individual aerosol opacities |
---|
503 | !========================================================================== |
---|
504 | ! Switch to netcdf define mode |
---|
505 | ierr=nf90_redef(outfid) |
---|
506 | write(*,*)"" |
---|
507 | SELECT CASE (iaer) |
---|
508 | CASE(1) ! DUST |
---|
509 | ! Insert the definition of the variable |
---|
510 | tmpvarname="opadust" |
---|
511 | ierr=NF90_DEF_VAR(outfid,trim(tmpvarname),nf90_float,varshape,tmpvaridout) |
---|
512 | error_text="ERROR: Couldn't create "//trim(tmpvarname)//" in the outfile" |
---|
513 | call status_check(ierr,error_text) |
---|
514 | write(*,*) trim(tmpvarname)//" has been created in the outfile" |
---|
515 | |
---|
516 | ! Write the attributes |
---|
517 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name","Dust extinction opacity at reference wavelength") |
---|
518 | |
---|
519 | CASE(2) ! WATER ICE |
---|
520 | ! Insert the definition of the variable |
---|
521 | tmpvarname="opawice" |
---|
522 | ierr=NF90_DEF_VAR(outfid,trim(tmpvarname),nf90_float,varshape,tmpvaridout) |
---|
523 | error_text="ERROR: Couldn't create "//trim(tmpvarname)//" in the outfile" |
---|
524 | call status_check(ierr,error_text) |
---|
525 | write(*,*) trim(tmpvarname)//" has been created in the outfile" |
---|
526 | |
---|
527 | ! Write the attributes |
---|
528 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name","Water ice extinction opacity at reference wavelength") |
---|
529 | |
---|
530 | ! + NEW AER? |
---|
531 | END SELECT ! iaer |
---|
532 | |
---|
533 | ierr=nf90_put_att(outfid,tmpvaridout,"units","opacity/km") |
---|
534 | ierr=nf90_put_att(outfid,tmpvaridout,"refwavelength",wvl_val) |
---|
535 | ierr=nf90_put_att(outfid,tmpvaridout,"missing_value",missval) |
---|
536 | write(*,*)"with missing value = ",missval |
---|
537 | |
---|
538 | ! End netcdf define mode |
---|
539 | ierr=nf90_enddef(outfid) |
---|
540 | |
---|
541 | !========================================================================== |
---|
542 | ! 3.2 Computation of the opacities |
---|
543 | !========================================================================== |
---|
544 | do ilon=1,lonlen |
---|
545 | do ilat=1,latlen |
---|
546 | do ialt=1,altlen |
---|
547 | do it=1,timelen |
---|
548 | ! Get the optpropfile reff values encompassing the GCM reff |
---|
549 | isize=1 |
---|
550 | do while((isize.le.nsize).and.(radiusdyn(isize).lt.reff(iaer,ilon,ilat,ialt,it))) |
---|
551 | isize=isize+1 |
---|
552 | enddo |
---|
553 | if ((isize.gt.nsize).or.(radiusdyn(1).gt.reff(iaer,ilon,ilat,ialt,it))) then |
---|
554 | ! write(*,*) "WARNING: the GCM reff (",reff(iaer,ilon,ilat,ialt,it),"m) is out of the bounds" |
---|
555 | ! write(*,*) "of the file ",trim(datadir)//'/'//trim(optpropfile) |
---|
556 | ! write(*,*) "(reff_inf=",radiusdyn(1),"m ; reff_sup=",radiusdyn(nsize),"m)" |
---|
557 | ! write(*,*) "A missing value will be written for opacity" |
---|
558 | |
---|
559 | ! NB: this should especially handle cases when reff=0 |
---|
560 | opa_aer(iaer,ilon,ilat,ialt,it)=missval |
---|
561 | else |
---|
562 | if (radiusdyn(isize).eq.reff(iaer,ilon,ilat,ialt,it)) then |
---|
563 | ! if the GCM reff is in the optpropfile |
---|
564 | isize1 = isize |
---|
565 | isize2 = isize |
---|
566 | else ! radius(isize)>reff and radiusdyn(isize-1)<reff |
---|
567 | isize1 = isize-1 |
---|
568 | isize2 = isize |
---|
569 | endif |
---|
570 | |
---|
571 | ! Linear interpolation in effective radius |
---|
572 | if (isize2.ne.isize1) then |
---|
573 | coeff = (reff(iaer,ilon,ilat,ialt,it)-radiusdyn(isize1))/(radiusdyn(isize2)-radiusdyn(isize1)) |
---|
574 | else |
---|
575 | coeff = 0. |
---|
576 | endif |
---|
577 | reffQext(1) = Qext(iwvl1,isize1)+coeff*(Qext(iwvl1,isize2)-Qext(iwvl1,isize1)) |
---|
578 | reffQext(2) = Qext(iwvl2,isize1)+coeff*(Qext(iwvl2,isize2)-Qext(iwvl2,isize1)) |
---|
579 | |
---|
580 | ! Linear interpolation in wavelength |
---|
581 | if (iwvl2.ne.iwvl1) then |
---|
582 | coeff = (wvl_val-wvl(iwvl1))/(wvl(iwvl2)-wvl(iwvl1)) |
---|
583 | else |
---|
584 | coeff = 0. |
---|
585 | endif |
---|
586 | Qext_val = reffQext(1)+coeff*(reffQext(2)-reffQext(1)) |
---|
587 | |
---|
588 | ! Computation of the opacity [1/km] |
---|
589 | opa_aer(iaer,ilon,ilat,ialt,it)= 750.*Qext_val*mmr(iaer,ilon,ilat,ialt,it)*rho(ilon,ilat,ialt,it)& |
---|
590 | / ( rho_aer(iaer) * reff(iaer,ilon,ilat,ialt,it) ) |
---|
591 | |
---|
592 | endif |
---|
593 | enddo ! it |
---|
594 | enddo ! ialt |
---|
595 | enddo ! ilat |
---|
596 | enddo ! ilon |
---|
597 | |
---|
598 | !========================================================================== |
---|
599 | ! 3.3 Writing in the output file |
---|
600 | !========================================================================== |
---|
601 | ierr = NF90_PUT_VAR(outfid, tmpvaridout, opa_aer(iaer,:,:,:,:)) |
---|
602 | error_text="Error: could not write "//trim(tmpvarname)//" data in the outfile" |
---|
603 | call status_check(ierr,error_text) |
---|
604 | |
---|
605 | !========================================================================== |
---|
606 | ! 3.4 Prepare next loop |
---|
607 | !========================================================================== |
---|
608 | DEALLOCATE(wvl) |
---|
609 | DEALLOCATE(radiusdyn) |
---|
610 | DEALLOCATE(Qext) |
---|
611 | |
---|
612 | endif ! if aerok(iaer) |
---|
613 | |
---|
614 | ENDDO ! iaer |
---|
615 | |
---|
616 | !=============================================================================== |
---|
617 | ! 4. Close the files and end the program |
---|
618 | !=============================================================================== |
---|
619 | ierr = nf90_close(gcmfid) |
---|
620 | error_text="Error: could not close file "//trim(gcmfile) |
---|
621 | call status_check(ierr,error_text) |
---|
622 | |
---|
623 | ierr = nf90_close(outfid) |
---|
624 | error_text="Error: could not close file "//trim(outfile) |
---|
625 | call status_check(ierr,error_text) |
---|
626 | |
---|
627 | write(*,*)"";write(*,*)"Program aeroptical completed!" |
---|
628 | |
---|
629 | end program aeroptical |
---|
630 | |
---|
631 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
632 | !! SUBROUTINES |
---|
633 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
634 | |
---|
635 | subroutine status_check(ierr,error_text) |
---|
636 | |
---|
637 | use netcdf |
---|
638 | implicit none |
---|
639 | !================================================================ |
---|
640 | ! Arguments: |
---|
641 | !================================================================ |
---|
642 | integer,intent(in) :: ierr ! NetCDF status number |
---|
643 | character(len=256),intent(in) :: error_text |
---|
644 | |
---|
645 | if (ierr.ne.nf90_noerr) then |
---|
646 | write(*,*)trim(error_text) |
---|
647 | write(*,*)trim(nf90_strerror(ierr)) |
---|
648 | stop |
---|
649 | endif |
---|
650 | |
---|
651 | end subroutine status_check |
---|
652 | |
---|
653 | |
---|
654 | !******************************************************************************* |
---|
655 | |
---|
656 | subroutine inidims(gcmfile,gcmfid,outfile,outfid,lonlen,latlen,altlen,timelen,& |
---|
657 | latdimout,londimout,altdimout,timedimout,& |
---|
658 | GCM_layers,layerdimout,interlayerdimout) |
---|
659 | |
---|
660 | !============================================================================== |
---|
661 | ! Purpose: |
---|
662 | ! Read the dimensions of the input file |
---|
663 | ! and write them in the output file |
---|
664 | !============================================================================== |
---|
665 | ! Remarks: |
---|
666 | ! The NetCDF files must be open |
---|
667 | !============================================================================== |
---|
668 | use netcdf |
---|
669 | |
---|
670 | implicit none |
---|
671 | |
---|
672 | !============================================================================== |
---|
673 | ! Arguments: |
---|
674 | !============================================================================== |
---|
675 | character(len=128),intent(in) :: gcmfile ! name of the netcdf GCM input file |
---|
676 | integer,intent(in) :: gcmfid ! [netcdf] GCM input file ID number |
---|
677 | character(len=128),intent(in) :: outfile ! name of the netcdf output file |
---|
678 | integer,intent(in) :: outfid ! [netcdf] file ID number |
---|
679 | |
---|
680 | integer,intent(out) :: lonlen,latlen,altlen,timelen |
---|
681 | ! nb of grid points along longitude,latitude,altitude,time |
---|
682 | integer,intent(out) :: latdimout,londimout,altdimout,timedimout |
---|
683 | ! latdimout: stores latitude dimension ID number |
---|
684 | ! londimout: stores longitude dimension ID number |
---|
685 | ! altdimout: stores altitude dimension ID number |
---|
686 | ! timedimout: stores time dimension ID number |
---|
687 | integer,intent(out) :: GCM_layers ! number of GCM layers |
---|
688 | integer,intent(out) :: layerdimout ! "GCM_layers" ID |
---|
689 | integer,intent(out) :: interlayerdimout ! "GCM_layers+1" ID |
---|
690 | |
---|
691 | !============================================================================== |
---|
692 | ! Local variables: |
---|
693 | !============================================================================== |
---|
694 | integer :: ierr ! [netcdf] subroutine returned error code |
---|
695 | character(len=256) :: error_text ! text to display in case of error |
---|
696 | |
---|
697 | integer :: tmpdimid,tmpvarid ! temporary store a dimension/variable ID number |
---|
698 | character (len=64) :: long_name,units,positive |
---|
699 | ! long_name(): [netcdf] long_name attribute |
---|
700 | ! units(): [netcdf] units attribute |
---|
701 | ! positive(): [netcdf] positive direction attribute (for altitude) |
---|
702 | real, dimension(:), allocatable:: lat,lon,alt,time,ctl |
---|
703 | ! lat(): array, stores latitude coordinates |
---|
704 | ! lon(): array, stores longitude coordinates |
---|
705 | ! alt(): array, stores altitude coordinates |
---|
706 | ! time(): array, stores time coordinates |
---|
707 | ! ctl(): array, stores controle variable |
---|
708 | integer :: ctllen ! nb of elements in the controle array |
---|
709 | integer :: tmpdimidout,tmpvaridout ! temporary stores a dimension/variable ID number |
---|
710 | |
---|
711 | !============================================================================== |
---|
712 | ! LONGITUDE |
---|
713 | !============================================================================== |
---|
714 | ! Get the dimension in GCM file |
---|
715 | ierr=nf90_inq_dimid(gcmfid,"longitude",tmpdimid) |
---|
716 | error_text="Error: Dimension <longitude> is missing in file "//trim(gcmfile) |
---|
717 | call status_check(ierr,error_text) |
---|
718 | |
---|
719 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=lonlen) |
---|
720 | allocate(lon(lonlen)) |
---|
721 | |
---|
722 | ! Create the dimension in output file |
---|
723 | ierr=NF90_DEF_DIM(outfid,"longitude",lonlen,londimout) |
---|
724 | error_text="Error: could not define the longitude dimension in the outfile" |
---|
725 | call status_check(ierr,error_text) |
---|
726 | |
---|
727 | ! Get the field in GCM file |
---|
728 | ierr=nf90_inq_varid(gcmfid,"longitude",tmpvarid) |
---|
729 | error_text="Error: Field <longitude> is missing in file "//trim(gcmfile) |
---|
730 | call status_check(ierr,error_text) |
---|
731 | |
---|
732 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,lon) |
---|
733 | error_text="Failed to load longitude" |
---|
734 | call status_check(ierr,error_text) |
---|
735 | |
---|
736 | ! Create the field in the output file |
---|
737 | ierr=NF90_DEF_VAR(outfid,"longitude",nf90_float,(/londimout/),tmpvaridout) |
---|
738 | error_text="Error: could not define the longitude variable in the outfile" |
---|
739 | call status_check(ierr,error_text) |
---|
740 | |
---|
741 | ! Get the field attributes in the GCM file |
---|
742 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
743 | if (ierr.ne.nf90_noerr) then |
---|
744 | ! if no attribute "long_name", try "title" |
---|
745 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
746 | endif |
---|
747 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
748 | |
---|
749 | ! Put the field attributes in the output file |
---|
750 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
751 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
752 | |
---|
753 | ! Write the field values in the output file |
---|
754 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
755 | error_text="Error: could not end the define mode of the outfile" |
---|
756 | call status_check(ierr,error_text) |
---|
757 | |
---|
758 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,lon) |
---|
759 | error_text="Error: could not write the longitude field in the outfile" |
---|
760 | call status_check(ierr,error_text) |
---|
761 | |
---|
762 | !============================================================================== |
---|
763 | ! LATITUDE |
---|
764 | !============================================================================== |
---|
765 | ! Switch to netcdf define mode |
---|
766 | ierr=nf90_redef(outfid) |
---|
767 | error_text="Error: could not switch to define mode in the outfile" |
---|
768 | call status_check(ierr,error_text) |
---|
769 | |
---|
770 | ! Get the dimension in GCM file |
---|
771 | ierr=nf90_inq_dimid(gcmfid,"latitude",tmpdimid) |
---|
772 | error_text="Error: Dimension <latitude> is missing in file "//trim(gcmfile) |
---|
773 | call status_check(ierr,error_text) |
---|
774 | |
---|
775 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=latlen) |
---|
776 | allocate(lat(latlen)) |
---|
777 | |
---|
778 | ! Create the dimension in output file |
---|
779 | ierr=NF90_DEF_DIM(outfid,"latitude",latlen,latdimout) |
---|
780 | error_text="Error: could not define the latitude dimension in the outfile" |
---|
781 | call status_check(ierr,error_text) |
---|
782 | |
---|
783 | ! Get the field in GCM file |
---|
784 | ierr=nf90_inq_varid(gcmfid,"latitude",tmpvarid) |
---|
785 | error_text="Error: Field <latitude> is missing in file "//trim(gcmfile) |
---|
786 | call status_check(ierr,error_text) |
---|
787 | |
---|
788 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,lat) |
---|
789 | error_text="Failed to load latitude" |
---|
790 | call status_check(ierr,error_text) |
---|
791 | |
---|
792 | ! Create the field in the output file |
---|
793 | ierr=NF90_DEF_VAR(outfid,"latitude",nf90_float,(/latdimout/),tmpvaridout) |
---|
794 | error_text="Error: could not define the latitude variable in the outfile" |
---|
795 | call status_check(ierr,error_text) |
---|
796 | |
---|
797 | ! Get the field attributes in the GCM file |
---|
798 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
799 | if (ierr.ne.nf90_noerr) then |
---|
800 | ! if no attribute "long_name", try "title" |
---|
801 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
802 | endif |
---|
803 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
804 | |
---|
805 | ! Put the field attributes in the output file |
---|
806 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
807 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
808 | |
---|
809 | ! Write the field values in the output file |
---|
810 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
811 | error_text="Error: could not end the define mode of the outfile" |
---|
812 | call status_check(ierr,error_text) |
---|
813 | |
---|
814 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,lat) |
---|
815 | error_text="Error: could not write the latitude field in the outfile" |
---|
816 | call status_check(ierr,error_text) |
---|
817 | |
---|
818 | !============================================================================== |
---|
819 | ! ALTITUDE |
---|
820 | !============================================================================== |
---|
821 | ! Switch to netcdf define mode |
---|
822 | ierr=nf90_redef(outfid) |
---|
823 | error_text="Error: could not switch to define mode in the outfile" |
---|
824 | call status_check(ierr,error_text) |
---|
825 | |
---|
826 | ! Get the dimension in GCM file |
---|
827 | ierr=nf90_inq_dimid(gcmfid,"altitude",tmpdimid) |
---|
828 | error_text="Error: Dimension <altitude> is missing in file "//trim(gcmfile) |
---|
829 | call status_check(ierr,error_text) |
---|
830 | |
---|
831 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=altlen) |
---|
832 | allocate(alt(altlen)) |
---|
833 | |
---|
834 | ! Create the dimension in output file |
---|
835 | ierr=NF90_DEF_DIM(outfid,"altitude",altlen,altdimout) |
---|
836 | error_text="Error: could not define the altitude dimension in the outfile" |
---|
837 | call status_check(ierr,error_text) |
---|
838 | |
---|
839 | ! Get the field in GCM file |
---|
840 | ierr=nf90_inq_varid(gcmfid,"altitude",tmpvarid) |
---|
841 | error_text="Error: Field <altitude> is missing in file "//trim(gcmfile) |
---|
842 | call status_check(ierr,error_text) |
---|
843 | |
---|
844 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,alt) |
---|
845 | error_text="Failed to load altitude" |
---|
846 | call status_check(ierr,error_text) |
---|
847 | |
---|
848 | ! Create the field in the output file |
---|
849 | ierr=NF90_DEF_VAR(outfid,"altitude",nf90_float,(/altdimout/),tmpvaridout) |
---|
850 | error_text="Error: could not define the altitude variable in the outfile" |
---|
851 | call status_check(ierr,error_text) |
---|
852 | |
---|
853 | ! Get the field attributes in the GCM file |
---|
854 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
855 | if (ierr.ne.nf90_noerr) then |
---|
856 | ! if no attribute "long_name", try "title" |
---|
857 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
858 | endif |
---|
859 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
860 | ierr=nf90_get_att(gcmfid,tmpvarid,"positive",positive) |
---|
861 | |
---|
862 | ! Put the field attributes in the output file |
---|
863 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
864 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
865 | ierr=nf90_put_att(outfid,tmpvaridout,"positive",positive) |
---|
866 | |
---|
867 | ! Write the field values in the output file |
---|
868 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
869 | error_text="Error: could not end the define mode of the outfile" |
---|
870 | call status_check(ierr,error_text) |
---|
871 | |
---|
872 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,alt) |
---|
873 | error_text="Error: could not write the altitude field in the outfile" |
---|
874 | call status_check(ierr,error_text) |
---|
875 | |
---|
876 | !============================================================================== |
---|
877 | ! TIME |
---|
878 | !============================================================================== |
---|
879 | ! Switch to netcdf define mode |
---|
880 | ierr=nf90_redef(outfid) |
---|
881 | error_text="Error: could not switch to define mode in the outfile" |
---|
882 | call status_check(ierr,error_text) |
---|
883 | |
---|
884 | ! Get the dimension in GCM file |
---|
885 | ierr=nf90_inq_dimid(gcmfid,"Time",tmpdimid) |
---|
886 | error_text="Error: Dimension <Time> is missing in file "//trim(gcmfile) |
---|
887 | call status_check(ierr,error_text) |
---|
888 | |
---|
889 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=timelen) |
---|
890 | allocate(time(timelen)) |
---|
891 | |
---|
892 | ! Create the dimension in output file |
---|
893 | ierr=NF90_DEF_DIM(outfid,"Time",timelen,timedimout) |
---|
894 | error_text="Error: could not define the time dimension in the outfile" |
---|
895 | call status_check(ierr,error_text) |
---|
896 | |
---|
897 | ! Get the field in GCM file |
---|
898 | ierr=nf90_inq_varid(gcmfid,"Time",tmpvarid) |
---|
899 | error_text="Error: Field <Time> is missing in file "//trim(gcmfile) |
---|
900 | call status_check(ierr,error_text) |
---|
901 | |
---|
902 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,time) |
---|
903 | error_text="Failed to load Time" |
---|
904 | call status_check(ierr,error_text) |
---|
905 | |
---|
906 | ! Create the field in the output file |
---|
907 | ierr=NF90_DEF_VAR(outfid,"Time",nf90_float,(/timedimout/),tmpvaridout) |
---|
908 | error_text="Error: could not define the Time variable in the outfile" |
---|
909 | call status_check(ierr,error_text) |
---|
910 | |
---|
911 | ! Get the field attributes in the GCM file |
---|
912 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
913 | if (ierr.ne.nf90_noerr) then |
---|
914 | ! if no attribute "long_name", try "title" |
---|
915 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
916 | endif |
---|
917 | ierr=nf90_get_att(gcmfid,tmpvarid,"units",units) |
---|
918 | |
---|
919 | ! Put the field attributes in the output file |
---|
920 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
921 | ierr=nf90_put_att(outfid,tmpvaridout,"units",units) |
---|
922 | |
---|
923 | ! Write the field values in the output file |
---|
924 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
925 | error_text="Error: could not end the define mode of the outfile" |
---|
926 | call status_check(ierr,error_text) |
---|
927 | |
---|
928 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,time) |
---|
929 | error_text="Error: could not write the Time field in the outfile" |
---|
930 | call status_check(ierr,error_text) |
---|
931 | |
---|
932 | !============================================================================== |
---|
933 | ! CONTROLE |
---|
934 | !============================================================================== |
---|
935 | ! Switch to netcdf define mode |
---|
936 | ierr=nf90_redef(outfid) |
---|
937 | error_text="Error: could not switch to define mode in the outfile" |
---|
938 | call status_check(ierr,error_text) |
---|
939 | |
---|
940 | ! Get the dimension in GCM file |
---|
941 | ierr=nf90_inq_dimid(gcmfid,"index",tmpdimid) |
---|
942 | error_text="Dimension <index> is missing in file "//trim(gcmfile)& |
---|
943 | //". We'll skip that one." |
---|
944 | if (ierr.ne.nf90_noerr) then |
---|
945 | write(*,*)trim(error_text) |
---|
946 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
947 | error_text="Error: could not end the define mode of the outfile" |
---|
948 | call status_check(ierr,error_text) |
---|
949 | else |
---|
950 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=ctllen) |
---|
951 | allocate(ctl(ctllen)) |
---|
952 | |
---|
953 | ! Create the dimension in output file |
---|
954 | ierr=NF90_DEF_DIM(outfid,"index",ctllen,tmpdimidout) |
---|
955 | error_text="Error: could not define the index dimension in the outfile" |
---|
956 | call status_check(ierr,error_text) |
---|
957 | |
---|
958 | ! Get the field in GCM file |
---|
959 | ierr=nf90_inq_varid(gcmfid,"controle",tmpvarid) |
---|
960 | error_text="Error: Field <controle> is missing in file "//trim(gcmfile) |
---|
961 | call status_check(ierr,error_text) |
---|
962 | |
---|
963 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,ctl) |
---|
964 | error_text="Failed to load ctl" |
---|
965 | call status_check(ierr,error_text) |
---|
966 | |
---|
967 | ! Create the field in the output file |
---|
968 | ierr=NF90_DEF_VAR(outfid,"controle",nf90_float,(/tmpdimidout/),tmpvaridout) |
---|
969 | error_text="Error: could not define the controle variable in the outfile" |
---|
970 | call status_check(ierr,error_text) |
---|
971 | |
---|
972 | ! Get the field attributes in the GCM file |
---|
973 | ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name) |
---|
974 | if (ierr.ne.nf90_noerr) then |
---|
975 | ! if no attribute "long_name", try "title" |
---|
976 | ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name) |
---|
977 | endif |
---|
978 | |
---|
979 | ! Put the field attributes in the output file |
---|
980 | ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name) |
---|
981 | |
---|
982 | ! Write the field values in the output file |
---|
983 | ierr=nf90_enddef(outfid) ! end netcdf define mode |
---|
984 | error_text="Error: could not end the define mode of the outfile" |
---|
985 | call status_check(ierr,error_text) |
---|
986 | |
---|
987 | ierr=NF90_PUT_VAR(outfid,tmpvaridout,ctl) |
---|
988 | error_text="Error: could not write the controle field in the outfile" |
---|
989 | call status_check(ierr,error_text) |
---|
990 | endif |
---|
991 | |
---|
992 | |
---|
993 | !============================================================================== |
---|
994 | ! Load size of aps() or sigma() (in case it is not altlen) |
---|
995 | !============================================================================== |
---|
996 | ! Switch to netcdf define mode |
---|
997 | ierr=nf90_redef(outfid) |
---|
998 | error_text="Error: could not switch to define mode in the outfile" |
---|
999 | call status_check(ierr,error_text) |
---|
1000 | |
---|
1001 | ! Default is that GCM_layers=altlen |
---|
1002 | ! but for outputs of zrecast, it may be a different value |
---|
1003 | ierr=nf90_inq_dimid(gcmfid,"GCM_layers",tmpdimid) |
---|
1004 | if (ierr.ne.nf90_noerr) then |
---|
1005 | ! didn't find a GCM_layers dimension; therefore we have: |
---|
1006 | GCM_layers=altlen |
---|
1007 | else |
---|
1008 | ! load value of GCM_layers |
---|
1009 | ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=GCM_layers) |
---|
1010 | endif |
---|
1011 | |
---|
1012 | ! Create the dimensions in output file |
---|
1013 | ierr = NF90_DEF_DIM(outfid,"GCM_layers",GCM_layers,layerdimout) |
---|
1014 | error_text="Error: could not define the GCM_layers dimension in the outfile" |
---|
1015 | call status_check(ierr,error_text) |
---|
1016 | ierr = NF90_DEF_DIM(outfid,"GCM_interlayers",GCM_layers+1,interlayerdimout) |
---|
1017 | error_text="Error: could not define the GCM_interlayers dimension in the outfile" |
---|
1018 | call status_check(ierr,error_text) |
---|
1019 | |
---|
1020 | ! End netcdf define mode |
---|
1021 | ierr=nf90_enddef(outfid) |
---|
1022 | error_text="Error: could not end the define mode of the outfile" |
---|
1023 | call status_check(ierr,error_text) |
---|
1024 | |
---|
1025 | end subroutine inidims |
---|
1026 | |
---|
1027 | |
---|
1028 | !******************************************************************************* |
---|
1029 | |
---|
1030 | subroutine init2(gcmfid,lonlen,latlen,altlen,GCM_layers, & |
---|
1031 | outfid,londimout,latdimout,altdimout, & |
---|
1032 | layerdimout,interlayerdimout) |
---|
1033 | !============================================================================== |
---|
1034 | ! Purpose: |
---|
1035 | ! Copy ap() , bp(), aps(), bps(), aire() and phisinit() |
---|
1036 | ! from input file to outpout file |
---|
1037 | !============================================================================== |
---|
1038 | ! Remarks: |
---|
1039 | ! The NetCDF files must be open |
---|
1040 | !============================================================================== |
---|
1041 | |
---|
1042 | use netcdf |
---|
1043 | |
---|
1044 | implicit none |
---|
1045 | |
---|
1046 | !============================================================================== |
---|
1047 | ! Arguments: |
---|
1048 | !============================================================================== |
---|
1049 | integer, intent(in) :: gcmfid ! NetCDF output file ID |
---|
1050 | integer, intent(in) :: lonlen ! # of grid points along longitude |
---|
1051 | integer, intent(in) :: latlen ! # of grid points along latitude |
---|
1052 | integer, intent(in) :: altlen ! # of grid points along latitude |
---|
1053 | integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers |
---|
1054 | integer, intent(in) :: outfid ! NetCDF output file ID |
---|
1055 | integer, intent(in) :: londimout ! longitude dimension ID |
---|
1056 | integer, intent(in) :: latdimout ! latitude dimension ID |
---|
1057 | integer, intent(in) :: altdimout ! altitude dimension ID |
---|
1058 | integer, intent(in) :: layerdimout ! GCM_layers dimension ID |
---|
1059 | integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID |
---|
1060 | !============================================================================== |
---|
1061 | ! Local variables: |
---|
1062 | !============================================================================== |
---|
1063 | real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates |
---|
1064 | real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates |
---|
1065 | real,dimension(:),allocatable :: sigma ! sigma levels |
---|
1066 | real,dimension(:,:),allocatable :: aire ! mesh areas |
---|
1067 | real,dimension(:,:),allocatable :: phisinit ! Ground geopotential |
---|
1068 | integer :: ierr |
---|
1069 | integer :: tmpvarid ! temporary variable ID |
---|
1070 | logical :: area ! is "aire" available ? |
---|
1071 | logical :: phis ! is "phisinit" available ? |
---|
1072 | logical :: hybrid ! are "aps" and "bps" available ? |
---|
1073 | logical :: apbp ! are "ap" and "bp" available ? |
---|
1074 | |
---|
1075 | !============================================================================== |
---|
1076 | ! 1. Read data from input file |
---|
1077 | !============================================================================== |
---|
1078 | |
---|
1079 | ! hybrid coordinate aps |
---|
1080 | !write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
1081 | allocate(aps(GCM_layers),stat=ierr) |
---|
1082 | if (ierr.ne.0) then |
---|
1083 | write(*,*) "init2: failed to allocate aps!" |
---|
1084 | stop |
---|
1085 | endif |
---|
1086 | ierr=nf90_inq_varid(gcmfid,"aps",tmpvarid) |
---|
1087 | if (ierr.ne.nf90_noerr) then |
---|
1088 | write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord." |
---|
1089 | hybrid=.false. |
---|
1090 | else |
---|
1091 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,aps) |
---|
1092 | hybrid=.true. |
---|
1093 | if (ierr.ne.nf90_noerr) then |
---|
1094 | stop "init2 Error: Failed reading aps" |
---|
1095 | endif |
---|
1096 | |
---|
1097 | ! hybrid coordinate bps |
---|
1098 | ! write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
1099 | allocate(bps(GCM_layers),stat=ierr) |
---|
1100 | if (ierr.ne.0) then |
---|
1101 | write(*,*) "init2: failed to allocate bps!" |
---|
1102 | stop |
---|
1103 | endif |
---|
1104 | ierr=nf90_inq_varid(gcmfid,"bps",tmpvarid) |
---|
1105 | if (ierr.ne.nf90_noerr) then |
---|
1106 | stop "init2 Error: Failed to get bps ID." |
---|
1107 | endif |
---|
1108 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,bps) |
---|
1109 | if (ierr.ne.nf90_noerr) then |
---|
1110 | stop "init2 Error: Failed reading bps" |
---|
1111 | endif |
---|
1112 | endif |
---|
1113 | |
---|
1114 | ! hybrid coordinate ap |
---|
1115 | allocate(ap(GCM_layers+1),stat=ierr) |
---|
1116 | if (ierr.ne.0) then |
---|
1117 | write(*,*) "init2: failed to allocate ap!" |
---|
1118 | stop |
---|
1119 | else |
---|
1120 | ierr=nf90_inq_varid(gcmfid,"ap",tmpvarid) |
---|
1121 | if (ierr.ne.nf90_noerr) then |
---|
1122 | write(*,*) "Ooops. Failed to get ap ID. OK." |
---|
1123 | apbp=.false. |
---|
1124 | else |
---|
1125 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,ap) |
---|
1126 | apbp=.true. |
---|
1127 | if (ierr.ne.nf90_noerr) then |
---|
1128 | stop "Error: Failed reading ap" |
---|
1129 | endif |
---|
1130 | endif |
---|
1131 | endif |
---|
1132 | |
---|
1133 | ! hybrid coordinate bp |
---|
1134 | allocate(bp(GCM_layers+1),stat=ierr) |
---|
1135 | if (ierr.ne.0) then |
---|
1136 | write(*,*) "init2: failed to allocate bp!" |
---|
1137 | stop |
---|
1138 | else |
---|
1139 | ierr=nf90_inq_varid(gcmfid,"bp",tmpvarid) |
---|
1140 | if (ierr.ne.nf90_noerr) then |
---|
1141 | write(*,*) "Ooops. Failed to get bp ID. OK." |
---|
1142 | apbp=.false. |
---|
1143 | else |
---|
1144 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,bp) |
---|
1145 | apbp=.true. |
---|
1146 | if (ierr.ne.nf90_noerr) then |
---|
1147 | stop "Error: Failed reading bp" |
---|
1148 | endif |
---|
1149 | endif |
---|
1150 | endif |
---|
1151 | |
---|
1152 | ! sigma levels (if any) |
---|
1153 | if (.not.hybrid) then |
---|
1154 | allocate(sigma(GCM_layers),stat=ierr) |
---|
1155 | if (ierr.ne.0) then |
---|
1156 | write(*,*) "init2: failed to allocate sigma" |
---|
1157 | stop |
---|
1158 | endif |
---|
1159 | ierr=nf90_inq_varid(gcmfid,"sigma",tmpvarid) |
---|
1160 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,sigma) |
---|
1161 | if (ierr.ne.nf90_noerr) then |
---|
1162 | stop "init2 Error: Failed reading sigma" |
---|
1163 | endif |
---|
1164 | endif ! of if (.not.hybrid) |
---|
1165 | |
---|
1166 | ! mesh area |
---|
1167 | allocate(aire(lonlen,latlen),stat=ierr) |
---|
1168 | if (ierr.ne.0) then |
---|
1169 | write(*,*) "init2: failed to allocate aire!" |
---|
1170 | stop |
---|
1171 | endif |
---|
1172 | ierr=nf90_inq_varid(gcmfid,"aire",tmpvarid) |
---|
1173 | if (ierr.ne.nf90_noerr) then |
---|
1174 | write(*,*)"init2 warning: Failed to get aire ID." |
---|
1175 | area = .false. |
---|
1176 | else |
---|
1177 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,aire) |
---|
1178 | if (ierr.ne.nf90_noerr) then |
---|
1179 | stop "init2 Error: Failed reading aire" |
---|
1180 | endif |
---|
1181 | area = .true. |
---|
1182 | endif |
---|
1183 | |
---|
1184 | ! ground geopotential phisinit |
---|
1185 | allocate(phisinit(lonlen,latlen),stat=ierr) |
---|
1186 | if (ierr.ne.0) then |
---|
1187 | write(*,*) "init2: failed to allocate phisinit!" |
---|
1188 | stop |
---|
1189 | endif |
---|
1190 | ierr=nf90_inq_varid(gcmfid,"phisinit",tmpvarid) |
---|
1191 | if (ierr.ne.nf90_noerr) then |
---|
1192 | write(*,*)"init2 warning: Failed to get phisinit ID." |
---|
1193 | phis = .false. |
---|
1194 | else |
---|
1195 | ierr=NF90_GET_VAR(gcmfid,tmpvarid,phisinit) |
---|
1196 | if (ierr.ne.nf90_noerr) then |
---|
1197 | stop "init2 Error: Failed reading phisinit" |
---|
1198 | endif |
---|
1199 | phis = .true. |
---|
1200 | endif |
---|
1201 | |
---|
1202 | !============================================================================== |
---|
1203 | ! 2. Write |
---|
1204 | !============================================================================== |
---|
1205 | |
---|
1206 | !============================================================================== |
---|
1207 | ! 2.2. Hybrid coordinates ap() , bp(), aps() and bps() |
---|
1208 | !============================================================================== |
---|
1209 | if(hybrid) then |
---|
1210 | ! define aps |
---|
1211 | ! Switch to netcdf define mode |
---|
1212 | ierr=nf90_redef(outfid) |
---|
1213 | ! Insert the definition of the variable |
---|
1214 | ierr=NF90_DEF_VAR(outfid,"aps",nf90_float,(/layerdimout/),tmpvarid) |
---|
1215 | if (ierr.ne.nf90_noerr) then |
---|
1216 | stop "init2 Error: Failed to define the variable aps" |
---|
1217 | endif |
---|
1218 | ! Write the attributes |
---|
1219 | ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid pressure at midlayers") |
---|
1220 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1221 | ! End netcdf define mode |
---|
1222 | ierr=nf90_enddef(outfid) |
---|
1223 | |
---|
1224 | ! write aps |
---|
1225 | ierr=NF90_PUT_VAR(outfid,tmpvarid,aps) |
---|
1226 | if (ierr.ne.nf90_noerr) then |
---|
1227 | stop "init2 Error: Failed to write aps" |
---|
1228 | endif |
---|
1229 | |
---|
1230 | ! define bps |
---|
1231 | ! Switch to netcdf define mode |
---|
1232 | ierr=nf90_redef(outfid) |
---|
1233 | ! Insert the definition of the variable |
---|
1234 | ierr=NF90_DEF_VAR(outfid,"bps",nf90_float,(/layerdimout/),tmpvarid) |
---|
1235 | if (ierr.ne.nf90_noerr) then |
---|
1236 | stop "init2 Error: Failed to define the variable bps" |
---|
1237 | endif |
---|
1238 | ! Write the attributes |
---|
1239 | ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at midlayers") |
---|
1240 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1241 | ! End netcdf define mode |
---|
1242 | ierr=nf90_enddef(outfid) |
---|
1243 | |
---|
1244 | ! write bps |
---|
1245 | ierr=NF90_PUT_VAR(outfid,tmpvarid,bps) |
---|
1246 | if (ierr.ne.nf90_noerr) then |
---|
1247 | stop "init2 Error: Failed to write bps" |
---|
1248 | endif |
---|
1249 | |
---|
1250 | if (apbp) then |
---|
1251 | ! define ap |
---|
1252 | |
---|
1253 | ! Switch to netcdf define mode |
---|
1254 | ierr=nf90_redef(outfid) |
---|
1255 | ! Insert the definition of the variable |
---|
1256 | ierr=NF90_DEF_VAR(outfid,"ap",nf90_float,(/interlayerdimout/),tmpvarid) |
---|
1257 | if (ierr.ne.nf90_noerr) then |
---|
1258 | stop "init2 Error: Failed to define the variable ap" |
---|
1259 | endif |
---|
1260 | ! Write the attributes |
---|
1261 | ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at interlayers") |
---|
1262 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1263 | ! End netcdf define mode |
---|
1264 | ierr=nf90_enddef(outfid) |
---|
1265 | |
---|
1266 | ! write ap |
---|
1267 | ierr=NF90_PUT_VAR(outfid,tmpvarid,ap) |
---|
1268 | if (ierr.ne.nf90_noerr) then |
---|
1269 | stop "Error: Failed to write ap" |
---|
1270 | endif |
---|
1271 | |
---|
1272 | ! define bp |
---|
1273 | |
---|
1274 | ! Switch to netcdf define mode |
---|
1275 | ierr=nf90_redef(outfid) |
---|
1276 | ! Insert the definition of the variable |
---|
1277 | ierr=NF90_DEF_VAR(outfid,"bp",nf90_float,(/interlayerdimout/),tmpvarid) |
---|
1278 | if (ierr.ne.nf90_noerr) then |
---|
1279 | stop "init2 Error: Failed to define the variable bp" |
---|
1280 | endif |
---|
1281 | ! Write the attributes |
---|
1282 | ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at interlayers") |
---|
1283 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1284 | ! End netcdf define mode |
---|
1285 | ierr=nf90_enddef(outfid) |
---|
1286 | |
---|
1287 | ! write bp |
---|
1288 | ierr=NF90_PUT_VAR(outfid,tmpvarid,bp) |
---|
1289 | if (ierr.ne.nf90_noerr) then |
---|
1290 | stop "Error: Failed to write bp" |
---|
1291 | endif |
---|
1292 | endif ! of if (apbp) |
---|
1293 | |
---|
1294 | else |
---|
1295 | |
---|
1296 | ! Switch to netcdf define mode |
---|
1297 | ierr=nf90_redef(outfid) |
---|
1298 | ! Insert the definition of the variable |
---|
1299 | ierr=NF90_DEF_VAR(outfid,"sigma",nf90_float,(/layerdimout/),tmpvarid) |
---|
1300 | if (ierr.ne.nf90_noerr) then |
---|
1301 | stop "init2 Error: Failed to define the variable sigma" |
---|
1302 | endif |
---|
1303 | ! Write the attributes |
---|
1304 | ierr=nf90_put_att(outfid,tmpvarid,"title","sigma at midlayers") |
---|
1305 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1306 | ! End netcdf define mode |
---|
1307 | ierr=nf90_enddef(outfid) |
---|
1308 | |
---|
1309 | ! write sigma |
---|
1310 | ierr=NF90_PUT_VAR(outfid,tmpvarid,sigma) |
---|
1311 | if (ierr.ne.nf90_noerr) then |
---|
1312 | stop "init2 Error: Failed to write sigma" |
---|
1313 | endif |
---|
1314 | endif ! of if (hybrid) |
---|
1315 | |
---|
1316 | !============================================================================== |
---|
1317 | ! 2.2. aire() and phisinit() |
---|
1318 | !============================================================================== |
---|
1319 | |
---|
1320 | if (area) then |
---|
1321 | |
---|
1322 | ! Switch to netcdf define mode |
---|
1323 | ierr=nf90_redef(outfid) |
---|
1324 | ! Insert the definition of the variable |
---|
1325 | ierr=NF90_DEF_VAR(outfid,"aire",nf90_float,(/londimout,latdimout/),tmpvarid) |
---|
1326 | if (ierr.ne.nf90_noerr) then |
---|
1327 | stop "init2 Error: Failed to define the variable aire" |
---|
1328 | endif |
---|
1329 | ! Write the attributes |
---|
1330 | ierr=nf90_put_att(outfid,tmpvarid,"title","Mesh area") |
---|
1331 | ierr=nf90_put_att(outfid,tmpvarid,"units","m2") |
---|
1332 | ! End netcdf define mode |
---|
1333 | ierr=nf90_enddef(outfid) |
---|
1334 | |
---|
1335 | ! write aire |
---|
1336 | ierr=NF90_PUT_VAR(outfid,tmpvarid,aire) |
---|
1337 | if (ierr.ne.nf90_noerr) then |
---|
1338 | stop "init2 Error: Failed to write aire" |
---|
1339 | endif |
---|
1340 | endif ! of if (area) |
---|
1341 | |
---|
1342 | IF (phis) THEN |
---|
1343 | |
---|
1344 | ! Switch to netcdf define mode |
---|
1345 | ierr=nf90_redef(outfid) |
---|
1346 | ! Insert the definition of the variable |
---|
1347 | ierr=NF90_DEF_VAR(outfid,"phisinit",nf90_float,(/londimout,latdimout/),tmpvarid) |
---|
1348 | if (ierr.ne.nf90_noerr) then |
---|
1349 | stop "init2 Error: Failed to define the variable phisinit" |
---|
1350 | endif |
---|
1351 | ! Write the attributes |
---|
1352 | ierr=nf90_put_att(outfid,tmpvarid,"title","Ground level geopotential") |
---|
1353 | ierr=nf90_put_att(outfid,tmpvarid,"units"," ") |
---|
1354 | ! End netcdf define mode |
---|
1355 | ierr=nf90_enddef(outfid) |
---|
1356 | |
---|
1357 | ! write phisinit |
---|
1358 | ierr=NF90_PUT_VAR(outfid,tmpvarid,phisinit) |
---|
1359 | if (ierr.ne.nf90_noerr) then |
---|
1360 | stop "init2 Error: Failed to write phisinit" |
---|
1361 | endif |
---|
1362 | |
---|
1363 | ENDIF ! of IF (phis) |
---|
1364 | |
---|
1365 | |
---|
1366 | ! Cleanup |
---|
1367 | if (allocated(aps)) deallocate(aps) |
---|
1368 | if (allocated(bps)) deallocate(bps) |
---|
1369 | if (allocated(ap)) deallocate(ap) |
---|
1370 | if (allocated(bp)) deallocate(bp) |
---|
1371 | if (allocated(sigma)) deallocate(sigma) |
---|
1372 | if (allocated(phisinit)) deallocate(phisinit) |
---|
1373 | if (allocated(aire)) deallocate(aire) |
---|
1374 | |
---|
1375 | end subroutine init2 |
---|