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

Last change on this file since 3026 was 2904, checked in by abierjon, 23 months ago

Mars GCM:
In aeroptical.F90, handle cases between :

  • "all opacities in the column are missing values" (in which case we set column-integrated tau=missval)

and

  • "at least one opacity is not a missing value but is 0" (in which case we set tau=0)

AB

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