source: trunk/LMDZ.MARS/util/aeroptical.F90 @ 2317

Last change on this file since 2317 was 2317, checked in by abierjon, 5 years ago

Mars GCM:
Creation of the utilitary program aeroptical.F90, whose goal is to compute the opacities [1/km]
of dust and water ice aerosols, from their mass mixing ratios and their effective radius
stored in a GCM file, all at a wavelength given by the user. To that end, it reads optical
properties stored in files of the datadir/ directory.
+ associated aeroptical.def

AB

File size: 51.4 KB
Line 
1program 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
34use netcdf
35
36!===============================================================================
37! Variable declarations
38!===============================================================================
39
40implicit none ! for no implicitly typed variables
41
42
43! GCM file
44character(len=128) :: gcmfile                 ! name of the netcdf GCM input file
45integer :: gcmfid                             ! [netcdf] GCM input file ID number
46integer :: ierr                               ! [netcdf] subroutine returned error code
47character(len=256) :: error_text              ! text to display in case of error
48integer :: tmpvarid                           ! temporary store a variable ID number
49integer :: lonlen,latlen,altlen,timelen       ! nb of grid points along longitude,latitude,altitude,time
50integer :: GCM_layers                         ! number of GCM layers
51integer :: layerdimout,interlayerdimout       ! "GCM_layers" and "GCM_layers+1" IDs
52
53logical,dimension(:),allocatable :: aerok     ! to know if the needed fields are in the GCM file
54real,dimension(:,:,:,:,:),allocatable :: mmr  ! aerosols mass mixing ratios [kg/kg]
55real,dimension(:,:,:,:,:),allocatable :: reff ! aerosols effective radii [m]
56real :: reffwice_val                          ! value given by the user for reffwice (0 if read in the GCM file) [m]
57real,dimension(:,:,:,:),allocatable :: rho    ! atmospheric density [kg.m-3]
58integer :: naerkind                           ! nb of aerosols to consider
59integer :: iaer                               ! aerosol kind index (1=dust, 2=water ice)
60integer :: ilon,ilat,ialt,it                  ! Loop indices for coordinates
61
62
63! Output file
64character(len=128) :: outfile                    ! name of the netcdf output file
65integer :: outfid                                ! [netcdf] file ID number
66integer :: 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
71integer :: tmpvaridout                           ! temporary stores a variable ID number
72
73real :: wvl_val                                  ! reference wavelength of the output opacities (given by the user) [m]
74integer :: varshape(4)                           ! stores a variable shape (of dimensions' IDs)
75character(len=16) :: tmpvarname                  ! temporary stores a variable name
76real,dimension(:,:,:,:,:),allocatable :: opa_aer ! Opacity of the aerosols [1/km]
77real :: missval                                  ! Value to put in outfile when the reff is out of bounds
78PARAMETER(missval=1E+20)
79
80
81! Optical properties (read in external ASCII files)
82character(len=256) :: datadir              ! directory containing the ASCII files
83character(len=128) :: optpropfile          ! name of the ASCII optical properties file
84logical :: file_ok                         ! to know if the file can be opened
85integer :: file_unit = 60                  !
86
87integer :: jfile                           ! ASCII file scan index
88logical :: endwhile                        ! Reading loop boolean
89character(len=132) :: scanline             ! ASCII file scanning line
90integer :: read_ok                         ! to know if the line is well read
91                   
92integer :: nwvl                            ! Number of wavelengths in the domain (VIS or IR)
93integer :: nsize                           ! Number of particle sizes stored in the file
94integer :: iwvl,isize                      ! Wavelength and Particle size loop indices
95real,dimension(:),allocatable :: wvl       ! Wavelength axis [m]
96real,dimension(:),allocatable :: radiusdyn ! Particle size axis [m]
97real,dimension(:,:),allocatable :: Qext    ! Extinction efficiency coefficient [/]
98
99
100! Opacity computation
101integer :: iwvl1,iwvl2,isize1,isize2     ! Wavelength and Particle size indices for the interpolation
102real :: coeff                            ! Interpolation coeficient
103real,dimension(2) :: reffQext            ! Qext after reff interpolation
104real :: Qext_val                         ! Qext after both interpolations
105real,dimension(:),allocatable :: rho_aer ! Density of the aerosols [kg.m-3]
106
107!===============================================================================
108! 0. USER INPUTS
109!===============================================================================
110write(*,*) "Welcome in the aerosol opacities' computation program !"
111write(*,*)
112
113! Ask the GCM file name
114WRITE(*,*) "Enter an input file name (GCM diagfi/stats/concat) :"
115READ(*,*) gcmfile
116
117! Ask the reffwice to the user
118write(*,*)""
119write(*,*) "The water ice effective radius computed by the GCM is known to be a bit inaccurate."
120write(*,*) "Which value do you want to use for it (in meters) ?"
121write(*,*) "(put 0 if you still want the program to read the value in "//trim(gcmfile)//")"
122READ(*,*) reffwice_val
123
124! Ask the wavelength to the user
125write(*,*)""
126WRITE(*,*) "Value of the reference wavelength for the opacities' computation (in meters) ?"
127READ(*,*) wvl_val
128
129! Ask the directory containing the optical properties files
130write(*,*)""
131WRITE(*,*) "In which directory should we look for the optical properties' files ?"
132READ(*,'(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
145ierr=nf90_open(gcmfile,nf90_nowrite,gcmfid) ! nowrite mode=the program can only read the opened file
146error_text="Error: could not open file "//trim(gcmfile)
147call status_check(ierr,error_text)
148
149!==========================================================================
150! --> 1.1.2 Creation of the outfile
151outfile=gcmfile(1:index(gcmfile, ".nc")-1)//"_OPA.nc"
152
153ierr=NF90_CREATE(outfile,nf90_clobber,outfid) ! NB: clobber mode=overwrite any dataset with the same file name !
154error_text="Error: could not create file "//trim(outfile)
155call status_check(ierr,error_text)
156write(*,*)"";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
162call 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
168call 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
177naerkind = 2 ! dust & water ice ! + NEW AER?
178
179! Length allocation of the variables
180allocate(mmr(naerkind,lonlen,latlen,altlen,timelen))
181allocate(reff(naerkind,lonlen,latlen,altlen,timelen))
182
183! Initialize aerok to .true. for all aerosols
184allocate(aerok(naerkind))
185aerok(:)=.true.
186
187!==========================================================================
188! --> 1.2.1 DUST
189
190! DUST MASS MIXING RATIO
191! Check that the GCM file contains that variable
192ierr=nf90_inq_varid(gcmfid,"dustq",tmpvarid)
193error_text="Failed to find variable dustq in "//trim(gcmfile)&
194            //" We'll skip the dust aerosol."
195if (ierr.ne.nf90_noerr) then
196  write(*,*)trim(error_text)
197  aerok(1)=.false.
198else
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)
204endif
205
206! DUST EFFECTIVE RADIUS
207if (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
222endif
223
224!==========================================================================
225! --> 1.2.2 WATER ICE
226
227! WATER ICE MASS MIXING RATIO
228! Check that the GCM file contains that variable
229ierr=nf90_inq_varid(gcmfid,"h2o_ice",tmpvarid)
230error_text="Failed to find variable h2o_ice in "//trim(gcmfile)&
231            //" We'll skip the water ice aerosol."
232if (ierr.ne.nf90_noerr) then
233  write(*,*)trim(error_text)
234  aerok(2)=.false.
235else
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)
241endif
242
243! WATER ICE EFFECTIVE RADIUS
244if (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
268endif
269
270! Check if there is still at least one aerosol to compute
271IF (.NOT.ANY(aerok)) THEN
272  write(*,*) "Neither dust nor water ice are found in the file. Better stop now..."
273  stop
274ENDIF
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
283ierr=nf90_inq_varid(gcmfid,"rho",tmpvarid)
284error_text="Failed to find variable rho in "//trim(gcmfile)&
285            //" We need it to compute the opacity [1/km]."
286call status_check(ierr,error_text)
287! Length allocation for each dimension of the 4D variable
288allocate(rho(lonlen,latlen,altlen,timelen))
289! Load datasets
290ierr=nf90_get_var(gcmfid,tmpvarid,rho)
291error_text="Failed to load atmospheric density"
292call status_check(ierr,error_text)
293write(*,*) "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
299varshape(1)=londimout
300varshape(2)=latdimout
301varshape(3)=altdimout
302varshape(4)=timedimout
303
304! Fill the aerosol density array
305allocate(rho_aer(naerkind))
306rho_aer(1)=2500. ! dust
307rho_aer(2)=920.  ! water ice
308! + NEW AER?
309
310allocate(opa_aer(naerkind,lonlen,latlen,altlen,timelen))
311
312
313!===============================================================================
314! 2. OPTICAL PROPERTIES
315!===============================================================================
316
317DO 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)
373reading1_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)
415reading2_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
614ENDDO ! iaer
615
616!===============================================================================
617! 4. Close the files and end the program
618!===============================================================================
619ierr = nf90_close(gcmfid)
620error_text="Error: could not close file "//trim(gcmfile)
621call status_check(ierr,error_text)
622
623ierr = nf90_close(outfid)
624error_text="Error: could not close file "//trim(outfile)
625call status_check(ierr,error_text)
626
627write(*,*)"";write(*,*)"Program aeroptical completed!"
628
629end program aeroptical
630
631!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
632!! SUBROUTINES
633!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
634
635subroutine status_check(ierr,error_text)
636
637use netcdf
638implicit none
639!================================================================
640!  Arguments:
641!================================================================
642integer,intent(in) :: ierr ! NetCDF status number
643character(len=256),intent(in) :: error_text
644
645if (ierr.ne.nf90_noerr) then
646  write(*,*)trim(error_text)
647  write(*,*)trim(nf90_strerror(ierr))
648  stop
649endif
650
651end subroutine status_check
652
653
654!*******************************************************************************
655
656subroutine 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!==============================================================================
668use netcdf
669
670implicit none
671
672!==============================================================================
673! Arguments:
674!==============================================================================
675character(len=128),intent(in) :: gcmfile ! name of the netcdf GCM input file
676integer,intent(in) :: gcmfid ! [netcdf] GCM input file ID number
677character(len=128),intent(in) :: outfile ! name of the netcdf output file
678integer,intent(in) :: outfid ! [netcdf] file ID number
679
680integer,intent(out) ::  lonlen,latlen,altlen,timelen
681! nb of grid points along longitude,latitude,altitude,time
682integer,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
687integer,intent(out) :: GCM_layers ! number of GCM layers
688integer,intent(out) :: layerdimout ! "GCM_layers" ID
689integer,intent(out) :: interlayerdimout ! "GCM_layers+1" ID
690
691!==============================================================================
692! Local variables:
693!==============================================================================
694integer :: ierr ! [netcdf] subroutine returned error code
695character(len=256) :: error_text ! text to display in case of error
696
697integer :: tmpdimid,tmpvarid ! temporary store a dimension/variable ID number
698character (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)
702real, 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
708integer :: ctllen ! nb of elements in the controle array
709integer :: tmpdimidout,tmpvaridout ! temporary stores a dimension/variable ID number
710
711!==============================================================================
712! LONGITUDE
713!==============================================================================
714! Get the dimension in GCM file
715ierr=nf90_inq_dimid(gcmfid,"longitude",tmpdimid)
716error_text="Error: Dimension <longitude> is missing in file "//trim(gcmfile)
717call status_check(ierr,error_text)
718
719ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=lonlen)
720allocate(lon(lonlen))
721
722! Create the dimension in output file
723ierr=NF90_DEF_DIM(outfid,"longitude",lonlen,londimout)
724error_text="Error: could not define the longitude dimension in the outfile"
725call status_check(ierr,error_text)
726
727! Get the field in GCM file
728ierr=nf90_inq_varid(gcmfid,"longitude",tmpvarid)
729error_text="Error: Field <longitude> is missing in file "//trim(gcmfile)
730call status_check(ierr,error_text)
731
732ierr=NF90_GET_VAR(gcmfid,tmpvarid,lon)
733error_text="Failed to load longitude"
734call status_check(ierr,error_text)
735
736! Create the field in the output file
737ierr=NF90_DEF_VAR(outfid,"longitude",nf90_float,(/londimout/),tmpvaridout)
738error_text="Error: could not define the longitude variable in the outfile"
739call status_check(ierr,error_text)
740
741! Get the field attributes in the GCM file
742ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name)
743if (ierr.ne.nf90_noerr) then
744! if no attribute "long_name", try "title"
745  ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name)
746endif
747ierr=nf90_get_att(gcmfid,tmpvarid,"units",units)
748
749! Put the field attributes in the output file
750ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name)
751ierr=nf90_put_att(outfid,tmpvaridout,"units",units)
752
753! Write the field values in the output file
754ierr=nf90_enddef(outfid) ! end netcdf define mode
755error_text="Error: could not end the define mode of the outfile"
756call status_check(ierr,error_text)
757
758ierr=NF90_PUT_VAR(outfid,tmpvaridout,lon)
759error_text="Error: could not write the longitude field in the outfile"
760call status_check(ierr,error_text)
761
762!==============================================================================
763! LATITUDE
764!==============================================================================
765! Switch to netcdf define mode
766ierr=nf90_redef(outfid)
767error_text="Error: could not switch to define mode in the outfile"
768call status_check(ierr,error_text)
769
770! Get the dimension in GCM file
771ierr=nf90_inq_dimid(gcmfid,"latitude",tmpdimid)
772error_text="Error: Dimension <latitude> is missing in file "//trim(gcmfile)
773call status_check(ierr,error_text)
774
775ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=latlen)
776allocate(lat(latlen))
777
778! Create the dimension in output file
779ierr=NF90_DEF_DIM(outfid,"latitude",latlen,latdimout)
780error_text="Error: could not define the latitude dimension in the outfile"
781call status_check(ierr,error_text)
782
783! Get the field in GCM file
784ierr=nf90_inq_varid(gcmfid,"latitude",tmpvarid)
785error_text="Error: Field <latitude> is missing in file "//trim(gcmfile)
786call status_check(ierr,error_text)
787
788ierr=NF90_GET_VAR(gcmfid,tmpvarid,lat)
789error_text="Failed to load latitude"
790call status_check(ierr,error_text)
791
792! Create the field in the output file
793ierr=NF90_DEF_VAR(outfid,"latitude",nf90_float,(/latdimout/),tmpvaridout)
794error_text="Error: could not define the latitude variable in the outfile"
795call status_check(ierr,error_text)
796
797! Get the field attributes in the GCM file
798ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name)
799if (ierr.ne.nf90_noerr) then
800! if no attribute "long_name", try "title"
801  ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name)
802endif
803ierr=nf90_get_att(gcmfid,tmpvarid,"units",units)
804
805! Put the field attributes in the output file
806ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name)
807ierr=nf90_put_att(outfid,tmpvaridout,"units",units)
808
809! Write the field values in the output file
810ierr=nf90_enddef(outfid) ! end netcdf define mode
811error_text="Error: could not end the define mode of the outfile"
812call status_check(ierr,error_text)
813
814ierr=NF90_PUT_VAR(outfid,tmpvaridout,lat)
815error_text="Error: could not write the latitude field in the outfile"
816call status_check(ierr,error_text)
817
818!==============================================================================
819! ALTITUDE
820!==============================================================================
821! Switch to netcdf define mode
822ierr=nf90_redef(outfid)
823error_text="Error: could not switch to define mode in the outfile"
824call status_check(ierr,error_text)
825
826! Get the dimension in GCM file
827ierr=nf90_inq_dimid(gcmfid,"altitude",tmpdimid)
828error_text="Error: Dimension <altitude> is missing in file "//trim(gcmfile)
829call status_check(ierr,error_text)
830
831ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=altlen)
832allocate(alt(altlen))
833
834! Create the dimension in output file
835ierr=NF90_DEF_DIM(outfid,"altitude",altlen,altdimout)
836error_text="Error: could not define the altitude dimension in the outfile"
837call status_check(ierr,error_text)
838
839! Get the field in GCM file
840ierr=nf90_inq_varid(gcmfid,"altitude",tmpvarid)
841error_text="Error: Field <altitude> is missing in file "//trim(gcmfile)
842call status_check(ierr,error_text)
843
844ierr=NF90_GET_VAR(gcmfid,tmpvarid,alt)
845error_text="Failed to load altitude"
846call status_check(ierr,error_text)
847
848! Create the field in the output file
849ierr=NF90_DEF_VAR(outfid,"altitude",nf90_float,(/altdimout/),tmpvaridout)
850error_text="Error: could not define the altitude variable in the outfile"
851call status_check(ierr,error_text)
852
853! Get the field attributes in the GCM file
854ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name)
855if (ierr.ne.nf90_noerr) then
856! if no attribute "long_name", try "title"
857  ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name)
858endif
859ierr=nf90_get_att(gcmfid,tmpvarid,"units",units)
860ierr=nf90_get_att(gcmfid,tmpvarid,"positive",positive)
861
862! Put the field attributes in the output file
863ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name)
864ierr=nf90_put_att(outfid,tmpvaridout,"units",units)
865ierr=nf90_put_att(outfid,tmpvaridout,"positive",positive)
866
867! Write the field values in the output file
868ierr=nf90_enddef(outfid) ! end netcdf define mode
869error_text="Error: could not end the define mode of the outfile"
870call status_check(ierr,error_text)
871
872ierr=NF90_PUT_VAR(outfid,tmpvaridout,alt)
873error_text="Error: could not write the altitude field in the outfile"
874call status_check(ierr,error_text)
875
876!==============================================================================
877! TIME
878!==============================================================================
879! Switch to netcdf define mode
880ierr=nf90_redef(outfid)
881error_text="Error: could not switch to define mode in the outfile"
882call status_check(ierr,error_text)
883
884! Get the dimension in GCM file
885ierr=nf90_inq_dimid(gcmfid,"Time",tmpdimid)
886error_text="Error: Dimension <Time> is missing in file "//trim(gcmfile)
887call status_check(ierr,error_text)
888
889ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=timelen)
890allocate(time(timelen))
891
892! Create the dimension in output file
893ierr=NF90_DEF_DIM(outfid,"Time",timelen,timedimout)
894error_text="Error: could not define the time dimension in the outfile"
895call status_check(ierr,error_text)
896
897! Get the field in GCM file
898ierr=nf90_inq_varid(gcmfid,"Time",tmpvarid)
899error_text="Error: Field <Time> is missing in file "//trim(gcmfile)
900call status_check(ierr,error_text)
901
902ierr=NF90_GET_VAR(gcmfid,tmpvarid,time)
903error_text="Failed to load Time"
904call status_check(ierr,error_text)
905
906! Create the field in the output file
907ierr=NF90_DEF_VAR(outfid,"Time",nf90_float,(/timedimout/),tmpvaridout)
908error_text="Error: could not define the Time variable in the outfile"
909call status_check(ierr,error_text)
910
911! Get the field attributes in the GCM file
912ierr=nf90_get_att(gcmfid,tmpvarid,"long_name",long_name)
913if (ierr.ne.nf90_noerr) then
914! if no attribute "long_name", try "title"
915  ierr=nf90_get_att(gcmfid,tmpvarid,"title",long_name)
916endif
917ierr=nf90_get_att(gcmfid,tmpvarid,"units",units)
918
919! Put the field attributes in the output file
920ierr=nf90_put_att(outfid,tmpvaridout,"long_name",long_name)
921ierr=nf90_put_att(outfid,tmpvaridout,"units",units)
922
923! Write the field values in the output file
924ierr=nf90_enddef(outfid) ! end netcdf define mode
925error_text="Error: could not end the define mode of the outfile"
926call status_check(ierr,error_text)
927
928ierr=NF90_PUT_VAR(outfid,tmpvaridout,time)
929error_text="Error: could not write the Time field in the outfile"
930call status_check(ierr,error_text)
931
932!==============================================================================
933! CONTROLE
934!==============================================================================
935! Switch to netcdf define mode
936ierr=nf90_redef(outfid)
937error_text="Error: could not switch to define mode in the outfile"
938call status_check(ierr,error_text)
939
940! Get the dimension in GCM file
941ierr=nf90_inq_dimid(gcmfid,"index",tmpdimid)
942error_text="Dimension <index> is missing in file "//trim(gcmfile)&
943           //". We'll skip that one."
944if (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)
949else
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)
990endif
991
992
993!==============================================================================
994! Load size of aps() or sigma() (in case it is not altlen)
995!==============================================================================
996! Switch to netcdf define mode
997ierr=nf90_redef(outfid)
998error_text="Error: could not switch to define mode in the outfile"
999call 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
1003ierr=nf90_inq_dimid(gcmfid,"GCM_layers",tmpdimid)
1004if (ierr.ne.nf90_noerr) then
1005  ! didn't find a GCM_layers dimension; therefore we have:
1006  GCM_layers=altlen 
1007else
1008  ! load value of GCM_layers
1009  ierr=nf90_inquire_dimension(gcmfid,tmpdimid,len=GCM_layers)
1010endif
1011
1012! Create the dimensions in output file
1013ierr = NF90_DEF_DIM(outfid,"GCM_layers",GCM_layers,layerdimout)
1014error_text="Error: could not define the GCM_layers dimension in the outfile"
1015call status_check(ierr,error_text)
1016ierr = NF90_DEF_DIM(outfid,"GCM_interlayers",GCM_layers+1,interlayerdimout)
1017error_text="Error: could not define the GCM_interlayers dimension in the outfile"
1018call status_check(ierr,error_text)
1019
1020! End netcdf define mode
1021ierr=nf90_enddef(outfid)
1022error_text="Error: could not end the define mode of the outfile"
1023call status_check(ierr,error_text)
1024
1025end subroutine inidims
1026
1027
1028!*******************************************************************************
1029
1030subroutine 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
1042use netcdf
1043
1044implicit none
1045
1046!==============================================================================
1047! Arguments:
1048!==============================================================================
1049integer, intent(in) :: gcmfid  ! NetCDF output file ID
1050integer, intent(in) :: lonlen ! # of grid points along longitude
1051integer, intent(in) :: latlen ! # of grid points along latitude
1052integer, intent(in) :: altlen ! # of grid points along latitude
1053integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers
1054integer, intent(in) :: outfid ! NetCDF output file ID
1055integer, intent(in) :: londimout ! longitude dimension ID
1056integer, intent(in) :: latdimout ! latitude dimension ID
1057integer, intent(in) :: altdimout ! altitude dimension ID
1058integer, intent(in) :: layerdimout ! GCM_layers dimension ID
1059integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID
1060!==============================================================================
1061! Local variables:
1062!==============================================================================
1063real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
1064real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates
1065real,dimension(:),allocatable :: sigma ! sigma levels
1066real,dimension(:,:),allocatable :: aire ! mesh areas
1067real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
1068integer :: ierr
1069integer :: tmpvarid ! temporary variable ID
1070logical :: area ! is "aire" available ?
1071logical :: phis ! is "phisinit" available ?
1072logical :: hybrid ! are "aps" and "bps" available ?
1073logical :: 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
1081allocate(aps(GCM_layers),stat=ierr)
1082if (ierr.ne.0) then
1083  write(*,*) "init2: failed to allocate aps!"
1084  stop
1085endif
1086ierr=nf90_inq_varid(gcmfid,"aps",tmpvarid)
1087if (ierr.ne.nf90_noerr) then
1088  write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord."
1089  hybrid=.false.
1090else
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
1112endif
1113
1114! hybrid coordinate ap
1115allocate(ap(GCM_layers+1),stat=ierr)
1116if (ierr.ne.0) then
1117  write(*,*) "init2: failed to allocate ap!"
1118  stop
1119else
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
1131endif
1132
1133! hybrid coordinate bp
1134allocate(bp(GCM_layers+1),stat=ierr)
1135if (ierr.ne.0) then
1136  write(*,*) "init2: failed to allocate bp!"
1137  stop
1138else
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
1150endif
1151
1152! sigma levels (if any)
1153if (.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
1164endif ! of if (.not.hybrid)
1165
1166! mesh area
1167allocate(aire(lonlen,latlen),stat=ierr)
1168if (ierr.ne.0) then
1169  write(*,*) "init2: failed to allocate aire!"
1170  stop
1171endif
1172ierr=nf90_inq_varid(gcmfid,"aire",tmpvarid)
1173if (ierr.ne.nf90_noerr) then
1174  write(*,*)"init2 warning: Failed to get aire ID."
1175  area = .false.
1176else
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.
1182endif
1183
1184! ground geopotential phisinit
1185allocate(phisinit(lonlen,latlen),stat=ierr)
1186if (ierr.ne.0) then
1187  write(*,*) "init2: failed to allocate phisinit!"
1188  stop
1189endif
1190ierr=nf90_inq_varid(gcmfid,"phisinit",tmpvarid)
1191if (ierr.ne.nf90_noerr) then
1192  write(*,*)"init2 warning: Failed to get phisinit ID."
1193  phis = .false.
1194else
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.
1200endif
1201
1202!==============================================================================
1203! 2. Write
1204!==============================================================================
1205
1206!==============================================================================
1207! 2.2. Hybrid coordinates ap() , bp(), aps() and bps()
1208!==============================================================================
1209if(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
1294else
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
1314endif ! of if (hybrid)
1315
1316!==============================================================================
1317! 2.2. aire() and phisinit()
1318!==============================================================================
1319
1320if (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
1340endif ! of if (area)
1341
1342IF (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
1363ENDIF ! of IF (phis)
1364
1365
1366! Cleanup
1367if (allocated(aps)) deallocate(aps)
1368if (allocated(bps)) deallocate(bps)
1369if (allocated(ap)) deallocate(ap)
1370if (allocated(bp)) deallocate(bp)
1371if (allocated(sigma)) deallocate(sigma)
1372if (allocated(phisinit)) deallocate(phisinit)
1373if (allocated(aire)) deallocate(aire)
1374
1375end subroutine init2
Note: See TracBrowser for help on using the repository browser.