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

Last change on this file since 1769 was 1544, checked in by emillour, 9 years ago

Mars GCM:

  • Bug fix in read_dust_scenario.F90 to handle the limit case of identical realday and time() values at the end of the year.

LS

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