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

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

Mars GCM:
Just some corrections of the preface comments in aeroptical.F90 following r2317

AB

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