source: trunk/LMDZ.MARS/libf/phymars/read_dust_scenario.F90 @ 636

Last change on this file since 636 was 607, checked in by emillour, 13 years ago

Mars GCM:

Updated model so that 1) reference pressure for opacity is now 610 Pa (and
not 700 Pa anymore) and 2) the new MY24-MY30 dust scenarios (input files
dust_MY24.nc, dust_MY25.nc, ..., dust_MY30.nc) can be used
(when setting iaervar=24,25,...,30 in callphys.def): changed
"readtesassim.F90" into "read_dust_scenario.F90" and apadpted aeropacity.F.
one can still use the "old" MY24-MY25-MY26 (input files dust_tes_MY24.nc,
dust_tes_MY25.nc and dust_tes_MY26.nc) scenarios by setting
iaervar=124, 125 or 126 in callphys.def.

EM

File size: 9.2 KB
Line 
1subroutine read_dust_scenario(ngrid,nlayer,zday,pplev,tauref)
2
3! Reading of the dust scenario file
4
5use netcdf
6implicit none
7
8#include "dimensions.h"
9#include "dimphys.h"
10#include "comgeomfi.h"
11#include "datafile.h"
12#include "callkeys.h"
13
14integer, intent(in) :: ngrid,nlayer
15real, intent(in) :: zday ! date in martian days
16real, dimension(ngrid,nlayer+1), intent(in) :: pplev
17real, dimension(ngrid), intent(out) :: tauref
18
19! Local variables
20
21real :: realday
22integer nid,nvarid,ierr
23integer ltloop,lsloop,iloop,jloop,varloop,ig
24real, dimension(2) :: taubuf
25real tau1(4),tau
26real alt(4)
27integer latp(4),lonp(4)
28integer yinf,ysup,xinf,xsup,tinf,tsup
29real latinf,latsup,loninf,lonsup
30real latintmp,lonintmp,latdeg,londeg
31real colat,dlat,dlon,colattmp
32logical, save :: firstcall=.true.
33logical :: timeflag
34real,save :: radeg,pi
35integer :: timedim,londim,latdim
36real, dimension(:), allocatable, save :: lat,lon,time
37real, dimension(:,:,:), allocatable, save :: tautes
38integer, save :: timelen,lonlen,latlen
39character(len=33),save :: filename
40
41realday=mod(zday,669.)
42
43if (firstcall) then
44   firstcall=.false.
45
46   pi=acos(-1.)
47   radeg=180/pi
48   
49   ! assimilated dust file: (NB: iaervar is a common in "callkeys.h")
50   ! iaervar=4 means read dust_tes.nc file
51   ! iaervar=24 means read dust_MY24.nc file
52   ! iaervar=25 means read dust_MY25.nc file
53   ! iaervar=26 means read dust_MY26.nc file, etc.
54   if (iaervar.eq.4) then
55   ! NB: 4: old TES assimilated MY24 dust scenarios (at 700Pa ref pressure!)
56     filename="dust_tes.nc"
57   else if (iaervar.eq.24) then
58     filename="dust_MY24.nc"
59   else if (iaervar.eq.25) then
60     filename="dust_MY25.nc"
61   else if (iaervar.eq.26) then
62     filename="dust_MY26.nc"
63   else if (iaervar.eq.27) then
64     filename="dust_MY27.nc"
65   else if (iaervar.eq.28) then
66     filename="dust_MY28.nc"
67   else if (iaervar.eq.29) then
68     filename="dust_MY29.nc"   
69   else if (iaervar.eq.30) then
70     filename="dust_MY30.nc"   
71   ! 124,125,126: old TES assimilated dust scenarios (at 700Pa ref pressure!)
72   else if  (iaervar.eq.124) then
73     filename="dust_tes_MY24.nc"
74   else if (iaervar.eq.125) then
75     filename="dust_tes_MY25.nc"
76   else if (iaervar.eq.126) then
77     filename="dust_tes_MY26.nc"
78   endif
79   
80   ! Note: datafile() is defined in "datafile.h"
81   !ierr=NF_OPEN(trim(datafile)//"/"//trim(filename),NF_NOWRITE,nid)
82   ierr=nf90_open(trim(datafile)//"/"//trim(filename),NF90_NOWRITE,nid)
83   IF (ierr.NE.nf90_noerr) THEN
84     write(*,*)'Problem opening ',trim(filename),' (in phymars/read_dust_scenario.F90)'
85     write(*,*)'It should be in :',trim(datafile),'/'
86     write(*,*)'1) You can change this directory address in callfis.def with'
87     write(*,*)'   datadir=/path/to/datafiles'
88     write(*,*)'2) If necessary, dust*.nc (and other datafiles)'
89     write(*,*)'   can be obtained online on:'
90     write(*,*)'   http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
91     CALL ABORT
92   ENDIF
93
94   ierr=nf90_inq_dimid(nid,"Time",timedim)
95   if (ierr.ne.nf90_noerr) then
96    ! 'Time' dimension not found, look for 'time'
97    ierr=nf90_inq_dimid(nid,"time",timedim)
98    if (ierr.ne.nf90_noerr) then
99      write(*,*)"Error: read_dust_scenario <time> not found"
100    endif
101   endif
102   ierr=nf90_inquire_dimension(nid,timedim,len=timelen)
103   ierr=nf90_inq_dimid(nid,"latitude",latdim)
104   ierr=nf90_inquire_dimension(nid,latdim,len=latlen)
105   ierr=nf90_inq_dimid(nid,"longitude",londim)
106   ierr=nf90_inquire_dimension(nid,londim,len=lonlen)
107
108
109   allocate(tautes(lonlen,latlen,timelen))
110   allocate(lat(latlen), lon(lonlen), time(timelen))
111
112   ierr=nf90_inq_varid(nid,"dustop",nvarid)
113   ierr=nf90_get_var(nid,nvarid,tautes)
114   IF (ierr .NE. nf90_noerr) THEN
115      PRINT*, "Error: read_dust_scenario <dustop> not found"
116      write(*,*)trim(nf90_strerror(ierr))
117      stop
118   ENDIF
119
120   ierr=nf90_inq_varid(nid,"Time",nvarid)
121   ierr=nf90_get_var(nid,nvarid,time)
122   IF (ierr .NE. nf90_noerr) THEN
123      PRINT*, "Error: read_dust_scenario <Time> not found"
124      write(*,*)trim(nf90_strerror(ierr))
125      stop
126   ENDIF
127
128   ierr=nf90_inq_varid(nid,"latitude",nvarid)
129   ierr=nf90_get_var(nid,nvarid,lat)
130   IF (ierr .NE. nf90_noerr) THEN
131      PRINT*, "Error: read_dust_scenario <latitude> not found"
132      write(*,*)trim(nf90_strerror(ierr))
133      stop
134   ENDIF
135
136   ierr=nf90_inq_varid(nid,"longitude",nvarid)
137   ierr=nf90_get_var(nid,nvarid,lon)
138   IF (ierr .NE. nf90_noerr) THEN
139      PRINT*, "Error: read_dust_scenario <longitude> not found"
140      write(*,*)trim(nf90_strerror(ierr))
141      stop
142   ENDIF
143
144   ierr=nf90_close(nid)
145
146endif ! of if (firstcall)
147
148do ig=1,ngrid
149
150! Find the four nearest points, arranged as follows:
151!                               1 2
152!                               3 4
153
154   latdeg=lati(ig)*radeg  ! latitude, in degrees
155   londeg=long(ig)*radeg  ! longitude, in degrees east
156   colat=90-latdeg        ! colatitude, in degrees
157
158! Ehouarn: rounding effects and/or specific compiler issues
159!          sometime cause londeg to be sligthly below -180 ...
160  if (londeg.lt.-180) then
161    ! check if it is by a large amount
162    if ((londeg+180).lt.-1.e-3) then
163      write(*,*) 'reattesassim: error!!'
164      write(*,*) ' ig=',ig,' londeg=',londeg
165      stop
166    else
167      londeg=-180
168    endif
169  endif
170
171 ! Find encompassing latitudes
172   if (colat<(90-lat(1))) then ! between north pole and lat(1)
173      ysup=1
174      yinf=1
175   else if (colat>=90-(lat(latlen))) then ! between south pole and lat(laten)
176      ysup=latlen
177      yinf=latlen
178   else ! general case
179      do iloop=2,latlen
180         if(colat<(90-lat(iloop))) then
181           ysup=iloop-1
182           yinf=iloop
183           exit
184         endif
185      enddo
186   endif
187
188   latinf=lat(yinf)
189   latsup=lat(ysup)
190
191
192 ! Find encompassing longitudes
193   ! Note: in input file, lon(1)=-180.
194   if (londeg>lon(lonlen)) then
195      xsup=1
196      xinf=lonlen
197      loninf=lon(xsup)
198      lonsup=180.0 ! since lon(1)=-180.0
199   else
200      do iloop=2,lonlen
201         if(londeg<lon(iloop)) then
202              xsup=iloop
203              xinf=iloop-1
204              exit
205         endif
206      enddo
207      loninf=lon(xinf)
208      lonsup=lon(xsup)
209   endif
210
211   if ((xsup.gt.lonlen).OR.(yinf.gt.latlen).OR.(xinf.lt.1)&
212     .OR.(ysup.lt.1)) then
213      write (*,*) "read_dust_scenario: SYSTEM ERROR on x or y in ts_gcm"
214      write (*,*) "xinf: ",xinf
215      write (*,*) "xsup: ",xsup
216      write (*,*) "yinf: ",yinf
217      write (*,*) "ysup: ",ysup
218      stop
219   endif
220
221!   loninf=lon(xinf)
222!   lonsup=lon(xsup)
223!   latinf=lat(yinf)
224!   latsup=lat(ysup)
225
226! The four neighbouring points are arranged as follows:
227!                               1 2
228!                               3 4
229
230   latp(1)=ysup
231   latp(2)=ysup
232   latp(3)=yinf
233   latp(4)=yinf
234
235   lonp(1)=xinf
236   lonp(2)=xsup
237   lonp(3)=xinf
238   lonp(4)=xsup
239
240! Linear interpolation on time, for all four neighbouring points
241
242   if ((realday<time(1)).or.(realday>time(timelen))) then
243      tinf=timelen
244      tsup=1
245      timeflag=.true.
246   else
247      timeflag=.false.
248      do iloop=2,timelen
249         if (realday<time(iloop)) then
250            tinf=iloop-1
251            tsup=iloop
252            exit
253         endif
254      enddo
255   endif
256
257! Bilinear interpolation on the four nearest points
258
259   if ((colat<(90-lat(1))).OR.(colat>(90-lat(latlen))).OR.(latsup==latinf)) then
260      dlat=0
261   else
262      dlat=((90-latinf)-colat)/((90-latinf)-(90-latsup))
263   endif
264
265   if (lonsup==loninf) then
266      dlon=0
267   else
268      dlon=(londeg-loninf)/(lonsup-loninf)
269   endif
270
271   do iloop=1,4
272      taubuf(1)=tautes(lonp(iloop),latp(iloop),tinf)
273      taubuf(2)=tautes(lonp(iloop),latp(iloop),tsup)
274      if (timeflag) then
275         if (realday>time(timelen)) then
276            tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)+(669-time(tinf)))
277         else
278            tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*realday/(time(tsup)+(669-time(tinf)))
279         endif
280      else
281         tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)-time(tinf))
282      endif
283      if (tau1(iloop)<0) then
284          write (*,*) "read_dust_scenario: SYSTEM ERROR on tau"
285          write (*,*) "utime ",realday
286          write (*,*) "time(tinf) ",time(tinf)
287          write (*,*) "time(tsup) ",time(tsup)
288          write (*,*) "tau1 ",taubuf(1)
289          write (*,*) "tau2 ",taubuf(2)
290          write (*,*) "tau ",tau1(iloop)
291          stop
292      endif
293   enddo
294
295   if ((dlat>1).OR.(dlon>1) .OR. (dlat<0) .OR. (dlon<0)) then
296      write (*,*) "read_dust_scenario: SYSTEM ERROR on dlat or dlon in ts_gcm"
297      write (*,*) "dlat: ",dlat
298      write (*,*) "lat: ",latdeg
299      write (*,*) "dlon: ",dlon
300      write (*,*) "lon: ",londeg
301      stop
302   endif
303
304   tau= dlat*(dlon*(tau1(2)+tau1(3)-tau1(1)-tau1(4))+tau1(1)-tau1(3)) +dlon*(tau1(4)-tau1(3))+tau1(3)
305
306   tauref(ig)=tau
307!
308enddo ! of ig=1,ngrid
309
310if (filename(1:8)=="dust_tes") then
311  ! when using old TES files, correction for:
312  ! - the reference pressure was 700Pa (unlike 610Pa now)
313  ! - the 1.3 conversion factor from IR absorbtion opacity to
314  !   IR extinction opacity
315  tauref(:)=tauref(:)*1.3*(610./700.)
316endif
317
318if (swrtype.eq.1) then ! Fouquart (NB: swrtype is set in callkeys.h)
319 ! when using old radiative transfer (like in MCD 4.x)
320 ! needed to decrease opacity (*0.825) to compensate overestimation of
321 ! heating rates
322  tauref(:)=tauref(:)*0.825/1.3
323endif
324
325end
Note: See TracBrowser for help on using the repository browser.