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

Last change on this file since 3493 was 3119, checked in by mwolff, 14 months ago

Small updates to aeroptical and simu_MCS (in LMDZ.MARS/util) to make
Mac OSX gfortran happy.

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