source: trunk/LMDZ.MARS/libf/phymars/read_dust_scenario_mod.F90 @ 3807

Last change on this file since 3807 was 3726, checked in by emillour, 2 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
EM

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