source: trunk/WRF.COMMON/WRFV2/mars_lmd/libf/phymars/readtesassim.F90.modif @ 3428

Last change on this file since 3428 was 341, checked in by aslmd, 13 years ago

MESOSCALE: tests pour faire marcher le modele en parallele sur la ferme. toujours infructueux... toutes les notes incluses et options explorees en commentaire. les options par defaut restent les memes en attendant. ajout de scripts pour compiler NETCDF et MPI. correction d un probleme de Registry et de makemeso pour les runs LES ancienne physique. ajout d un cas test LES phoenix.

File size: 7.5 KB
Line 
1subroutine readtesassim(ngrid,nlayer,zday,pplev,tauref)
2
3! Reading of the dust assimilation file
4
5implicit none
6
7include "dimensions.h"
8include "dimphys.h"
9include "comgeomfi.h"
10include "netcdf.inc"
11include "datafile.h"
12
13integer, intent(in) :: ngrid,nlayer
14real, intent(in) :: zday ! date in martian days
15real, dimension(ngrid,nlayer+1), intent(in) :: pplev
16real, dimension(ngrid), intent(out) :: tauref
17
18! Local variables
19
20real :: realday
21integer nid,nvarid,ierr
22integer ltloop,lsloop,iloop,jloop,varloop,ig
23real, dimension(2) :: taubuf
24real tau1(4),tau
25real alt(4)
26integer latp(4),lonp(4)
27integer yinf,ysup,xinf,xsup,tinf,tsup
28real latinf,latsup,loninf,lonsup
29real latintmp,lonintmp,latdeg,londeg
30real colat,dlat,dlon,colattmp
31logical, save :: firstcall=.true.
32logical :: timeflag
33real radeg,pi
34integer :: timedim,londim,latdim
35!real, dimension(:), allocatable, save :: lat,lon,time
36!real, dimension(:,:,:), allocatable, save :: tautes
37!integer, save :: timelen,lonlen,latlen
38real, dimension(:), allocatable :: lat,lon,time
39real, dimension(:,:,:), allocatable :: tautes
40integer :: timelen,lonlen,latlen
41INTEGER, external:: lnblnk
42
43pi=acos(-1.)
44radeg=180/pi
45realday=mod(zday,669.)
46
47if (firstcall) then
48   !firstcall=.false.
49
50   ! Note: datafile() is defined in "datafile.h"
51   ierr=NF_OPEN (datafile(1:lnblnk(datafile))//"/dust_tes.nc",NF_NOWRITE,nid)
52   IF (ierr.NE.NF_NOERR) THEN
53     write(*,*)'Problem opening dust.nc (in phymars/readtesassim.F90)'
54     write(*,*)'It should be in :',datafile(1:lnblnk(datafile)),'/'
55     write(*,*)'1) You can change this directory address in '
56     write(*,*)'   file phymars/datafile.h'
57     write(*,*)'2) If necessary, dust.nc (and other datafiles)'
58     write(*,*)'   can be obtained online on:'
59     write(*,*)'   http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
60     CALL ABORT
61   ENDIF
62
63   ierr=NF_INQ_DIMID(nid,"Time",timedim)
64   ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
65   ierr=NF_INQ_DIMID(nid,"latitude",latdim)
66   ierr=NF_INQ_DIMLEN(nid,latdim,latlen)
67   ierr=NF_INQ_DIMID(nid,"longitude",londim)
68   ierr=NF_INQ_DIMLEN(nid,londim,lonlen)
69
70
71   allocate(tautes(lonlen,latlen,timelen))
72   allocate(lat(latlen), lon(lonlen), time(timelen))
73
74   ierr = NF_INQ_VARID (nid, "dustop",nvarid)
75#ifdef NC_DOUBLE
76   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tautes)
77#else
78   ierr = NF_GET_VAR_REAL(nid, nvarid, tautes)
79#endif
80   IF (ierr .NE. NF_NOERR) THEN
81      PRINT*, "Error: Readtesassim <dustop> not found"
82      stop
83   ENDIF
84   ierr = NF_INQ_VARID (nid, "Time",nvarid)
85#ifdef NC_DOUBLE
86   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, time)
87#else
88   ierr = NF_GET_VAR_REAL(nid, nvarid, time)
89#endif
90   IF (ierr .NE. NF_NOERR) THEN
91      PRINT*, "Error: Readtesassim <Time> not found"
92      stop
93   ENDIF
94
95   ierr = NF_INQ_VARID (nid, "latitude",nvarid)
96#ifdef NC_DOUBLE
97   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, lat)
98#else
99   ierr = NF_GET_VAR_REAL(nid, nvarid, lat)
100#endif
101   IF (ierr .NE. NF_NOERR) THEN
102      PRINT*, "Error: Readtesassim <latitude> not found"
103      stop
104   ENDIF
105
106   ierr = NF_INQ_VARID (nid, "longitude",nvarid)
107#ifdef NC_DOUBLE
108   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, lon)
109#else
110   ierr = NF_GET_VAR_REAL(nid, nvarid, lon)
111#endif
112   IF (ierr .NE. NF_NOERR) THEN
113      PRINT*, "Error: Readtesassim <longitude> not found"
114      stop
115   ENDIF
116
117   ierr=NF_CLOSE(nid)
118
119endif ! of if (firstcall)
120
121do ig=1,ngrid
122
123! Find the four nearest points, arranged as follows:
124!                               1 2
125!                               3 4
126
127   latdeg=lati(ig)*radeg  ! latitude, in degrees
128   londeg=long(ig)*radeg  ! longitude, in degrees east
129   colat=90-latdeg        ! colatitude, in degrees
130
131! Ehouarn: rounding effects and/or specific compiler issues
132!          sometime cause londeg to be sligthly below -180 ...
133  if (londeg.lt.-180) then
134    ! check if it is by a large amount
135    if ((londeg+180).lt.-1.e-3) then
136      write(*,*) 'reattesassim: error!!'
137      write(*,*) ' ig=',ig,' londeg=',londeg
138      stop
139    else
140      londeg=-180
141    endif
142  endif
143
144 ! Find encompassing latitudes
145   if (colat<(90-lat(1))) then ! between north pole and lat(1)
146      ysup=1
147      yinf=1
148   else if (colat>=90-(lat(latlen))) then ! between south pole and lat(laten)
149      ysup=latlen
150      yinf=latlen
151   else ! general case
152      do iloop=2,latlen
153         if(colat<(90-lat(iloop))) then
154           ysup=iloop-1
155           yinf=iloop
156           exit
157         endif
158      enddo
159   endif
160
161   latinf=lat(yinf)
162   latsup=lat(ysup)
163
164
165 ! Find encompassing longitudes
166   ! Note: in input file, lon(1)=-180.
167   if (londeg>lon(lonlen)) then
168      xsup=1
169      xinf=lonlen
170      loninf=lon(xsup)
171      lonsup=180.0 ! since lon(1)=-180.0
172   else
173      do iloop=2,lonlen
174         if(londeg<lon(iloop)) then
175              xsup=iloop
176              xinf=iloop-1
177              exit
178         endif
179      enddo
180      loninf=lon(xinf)
181      lonsup=lon(xsup)
182   endif
183
184   if ((xsup.gt.lonlen).OR.(yinf.gt.latlen).OR.(xinf.lt.1)&
185     .OR.(ysup.lt.1)) then
186      write (*,*) "Readtesassim: SYSTEM ERROR on x or y in ts_gcm"
187      write (*,*) "xinf: ",xinf
188      write (*,*) "xsup: ",xsup
189      write (*,*) "yinf: ",yinf
190      write (*,*) "ysup: ",ysup
191      stop
192   endif
193
194!   loninf=lon(xinf)
195!   lonsup=lon(xsup)
196!   latinf=lat(yinf)
197!   latsup=lat(ysup)
198
199! The four neighbouring points are arranged as follows:
200!                               1 2
201!                               3 4
202
203   latp(1)=ysup
204   latp(2)=ysup
205   latp(3)=yinf
206   latp(4)=yinf
207
208   lonp(1)=xinf
209   lonp(2)=xsup
210   lonp(3)=xinf
211   lonp(4)=xsup
212
213! Linear interpolation on time, for all four neighbouring points
214
215   if ((realday<time(1)).or.(realday>time(timelen))) then
216      tinf=timelen
217      tsup=1
218      timeflag=.true.
219   else
220      timeflag=.false.
221      do iloop=2,timelen
222         if (realday<time(iloop)) then
223            tinf=iloop-1
224            tsup=iloop
225            exit
226         endif
227      enddo
228   endif
229
230! Bilinear interpolation on the four nearest points
231
232   if ((colat<(90-lat(1))).OR.(colat>(90-lat(latlen))).OR.(latsup==latinf)) then
233      dlat=0
234   else
235      dlat=((90-latinf)-colat)/((90-latinf)-(90-latsup))
236   endif
237
238   if (lonsup==loninf) then
239      dlon=0
240   else
241      dlon=(londeg-loninf)/(lonsup-loninf)
242   endif
243
244   do iloop=1,4
245      taubuf(1)=tautes(lonp(iloop),latp(iloop),tinf)
246      taubuf(2)=tautes(lonp(iloop),latp(iloop),tsup)
247      if (timeflag) then
248         if (realday>time(timelen)) then
249            tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)+(669-time(tinf)))
250         else
251            tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*realday/(time(tsup)+(669-time(tinf)))
252         endif
253      else
254         tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)-time(tinf))
255      endif
256      if (tau1(iloop)<0) then
257          write (*,*) "Readtesassim: SYSTEM ERROR on tau"
258          write (*,*) "utime ",realday
259          write (*,*) "time(tinf) ",time(tinf)
260          write (*,*) "time(tsup) ",time(tsup)
261          write (*,*) "tau1 ",taubuf(1)
262          write (*,*) "tau2 ",taubuf(2)
263          write (*,*) "tau ",tau1(iloop)
264          stop
265      endif
266   enddo
267
268   if ((dlat>1).OR.(dlon>1) .OR. (dlat<0) .OR. (dlon<0)) then
269      write (*,*) "Readtesassim: SYSTEM ERROR on dlat or dlon in ts_gcm"
270      write (*,*) "dlat: ",dlat
271      write (*,*) "lat: ",latdeg
272      write (*,*) "dlon: ",dlon
273      write (*,*) "lon: ",londeg
274      stop
275   endif
276
277   tau= dlat*(dlon*(tau1(2)+tau1(3)-tau1(1)-tau1(4))+tau1(1)-tau1(3)) +dlon*(tau1(4)-tau1(3))+tau1(3)
278
279!   tauref(ig)=tau*700/pplev(ig,1)*0.8
280    tauref(ig)=tau*0.825
281enddo
282
283   deallocate(tautes)
284   deallocate(lat)
285   deallocate(lon)
286   deallocate(time)
287
288end
Note: See TracBrowser for help on using the repository browser.