source: trunk/LMDZ.GENERIC/libf/phystd/datareadnc.F @ 486

Last change on this file since 486 was 374, checked in by emillour, 13 years ago

Generic GCM: Upgrade: The location of the 'datagcm' directory can now be given in the callphys.def file ( datadir = /absolute/path/to/datagcm ). Changed "datafile.h" into a F90 module "datafile_mod.F90" and spread this change to all routines that used to use "datafile.h".
EM

File size: 9.6 KB
Line 
1c=======================================================================
2      SUBROUTINE datareadnc(relief,phisinit,alb,ith,
3     .                    zmea,zstd,zsig,zgam,zthe)
4c=======================================================================
5c
6c
7c   Author: F. Hourdin      01/1997
8c   -------
9c
10c   Object: To read data from Martian surface to use in a GCM
11c   ------                from NetCDF file "surface.nc"
12c
13c
14c   Arguments:
15c   ----------
16c
17c     Inputs:
18c     ------
19c
20c     Outputs:
21c     --------
22c
23c=======================================================================
24c   donnees ALBEDO, INERTIE THERMIQUE, RELIEF:
25c
26c       Ces donnees sont au format NetCDF dans le fichier "surface.nc"
27c
28c   360 valeurs en longitude (de -179.5 a 179.5)
29c   180 valeurs en latitudes (de 89.5 a -89.5)
30c
31c   Pour les passer au format de la grille, on utilise "interp_horiz.F"
32c
33c   Il faut donc que ces donnees soient au format grille scalaire
34c               (imold+1 jmold+1)
35c       avec passage des coordonnees de la "boite" (rlonu, rlatv)
36c
37c   On prend imd (d pour donnees!)
38c           imd = 360 avec copie de la 1ere valeur sur la imd+1
39c                   (rlonud de -179 a -181)
40c           jmd = 179
41c                   (rlatvd de 89 a -89)
42c=======================================================================
43
44      use datafile_mod, only: datadir
45      implicit none
46
47#include "dimensions.h"
48#include "paramet.h"
49#include "comgeom.h"
50#include "comconst.h"
51#include "netcdf.inc"
52
53c=======================================================================
54c   Declarations:
55C=======================================================================
56
57      INTEGER    imd,jmd,imdp1,jmdp1
58      parameter    (imd=360,jmd=179,imdp1=361,jmdp1=180)
59
60      INTEGER    iimp1
61      parameter    (iimp1=iim+1-1/iim)
62
63      CHARACTER relief*3
64
65      REAL        zdata(imd*jmdp1)
66      REAL        zdataS(imdp1*jmdp1)
67      REAL        pfield(iimp1*jjp1)
68
69      REAL        alb(iimp1*jjp1)
70      REAL        ith(iimp1*jjp1)
71      REAL        phisinit(iimp1*jjp1)
72
73      REAL        zmea(imdp1*jmdp1)
74      REAL        zstd(imdp1*jmdp1)
75      REAL        zsig(imdp1*jmdp1)
76      REAL        zgam(imdp1*jmdp1)
77      REAL        zthe(imdp1*jmdp1)
78
79      INTEGER   lnblnk, ierr
80      EXTERNAL    lnblnk
81
82      INTEGER   unit,nvarid
83
84      INTEGER    i,j,k
85
86      INTEGER klatdat,ngridmxgdat
87      PARAMETER (klatdat=180,ngridmxgdat=360)
88
89c    on passe une grille en rlonu rlatv et im+1 jm a interp_horiz)
90
91      REAL longitude(imd),latitude(jmdp1) ! Pour lecture des donnees
92      REAL rlonud(imdp1),rlatvd(jmd)
93
94      CHARACTER*20 string
95      DIMENSION string(4)
96
97      CHARACTER*20 filename
98      CHARACTER*50 filestring
99
100#include "lmdstd.h"
101#include "fxyprim.h"
102
103      pi=2.*ASIN(1.)
104
105c=======================================================================
106c    rlonud, rlatvd
107c=======================================================================
108
109c-----------------------------------------------------------------------
110c    Lecture NetCDF des donnees latitude et longitude
111c-----------------------------------------------------------------------
112
113      write(*,*) 'Choice of surface data is:'
114      filestring='ls '//trim(datadir)//'/*.nc'
115      call system(filestring)
116
117      write(*,*) ''
118      write(*,*) 'Please choose the relevant file:'
119      write(*,*) '(e.g. type "surface_mars.nc")'
120      read(*,fmt='(a20)') filename
121
122      ierr = NF_OPEN (trim(datadir)//'/'//trim(adjustl(filename)),
123     &  NF_NOWRITE,unit)
124      IF (ierr.NE.NF_NOERR) THEN
125        write(*,*)'Error : cannot open file '//trim(filename)
126        write(*,*)'(in phystd/datareadnc.F)'
127        write(*,*)'It should be in :',trim(datadir),'/'
128        write(*,*)'Check that your path to datagcm:',trim(datadir)
129        write(*,*)' is correct. You can change it in callphys.def with:'
130        write(*,*)' datadir = /absolute/path/to/datagcm'
131        write(*,*)'If necessary surface.nc (and other datafiles)'
132        write(*,*)' can be obtained online on:'
133        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
134        CALL ABORT
135      ENDIF
136
137c
138c Lecture des latitudes (coordonnees):
139c
140      ierr = NF_INQ_VARID (unit, "latitude", nvarid)
141#ifdef NC_DOUBLE
142      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, latitude)
143#else
144      ierr = NF_GET_VAR_REAL(unit, nvarid, latitude)
145#endif
146c
147c Lecture des longitudes (coordonnees):
148c
149      ierr = NF_INQ_VARID (unit, "longitude", nvarid)
150#ifdef NC_DOUBLE
151      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, longitude)
152#else
153      ierr = NF_GET_VAR_REAL(unit, nvarid, longitude)
154#endif
155
156c-----------------------------------------------------------------------
157c    Passage au format boites scalaires
158c-----------------------------------------------------------------------
159
160c-----------------------------------------------------------------------
161c       longitude(imd)        -->      rlonud(imdp1)
162c-----------------------------------------------------------------------
163
164c Passage en coordonnees boites scalaires et en radian
165      do i=1,imd
166          rlonud(i)=(longitude(i)+.5)*pi/180.
167      enddo
168
169c Repetition de la valeur im+1
170      rlonud(imdp1)=rlonud(1) + 2*pi
171
172c-----------------------------------------------------------------------
173c        latitude(jmdp1)         -->        rlonvd(jmd)
174c-----------------------------------------------------------------------
175
176c Passage en coordonnees boites scalaires et en radian
177      do j=1,jmd
178          rlatvd(j)=(latitude(j)-.5)*pi/180.
179      enddo
180
181c=======================================================================
182c   lecture NetCDF de albedo, thermal, relief, zdtm (pour francois Lott)
183c=======================================================================
184
185      string(1) = 'albedo'
186      string(2) = 'thermal'
187      if (relief.ne.'pla') then
188        write(*,*) ' La topographie est celle de MOLA'
189        relief = 'MOL'
190          string(3) = 'z'//relief
191      else
192          string(3) = 'zMOL'  ! pour qu''il lise qqchose sur le fichier
193                            ! remise a 0 derriere
194      endif
195      string(4) = 'zMOL'    ! lecture pour calcul topog. sous-maille
196 
197
198      DO k=1,4
199          write(*,*) 'string',k,string(k)
200         
201c-----------------------------------------------------------------------
202c    initialisation
203c-----------------------------------------------------------------------
204      call initial0(iimp1*jjp1,pfield)
205      call initial0(imd*jmdp1,zdata)
206      call initial0(imdp1*jmdp1,zdataS)
207
208c-----------------------------------------------------------------------
209c    Lecture NetCDF 
210c-----------------------------------------------------------------------
211
212      ierr = NF_INQ_VARID (unit, string(k), nvarid)
213#ifdef NC_DOUBLE
214      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, zdata)
215#else
216      ierr = NF_GET_VAR_REAL(unit, nvarid, zdata)
217#endif
218
219c-----------------------------------------------------------------------
220c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
221c-----------------------------------------------------------------------
222      if (k.eq.4) then
223
224          call multscal(imd*jmdp1,zdata,1000.,zdata)
225          call multscal(imd,longitude,pi/180.,longitude)
226          call multscal(jmdp1,latitude,pi/180.,latitude)
227
228          call grid_noro1(360, 180, longitude, latitude, zdata,
229     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
230
231          CALL dump2d(iip1,jjp1,zmea,'zmea')
232          CALL dump2d(iip1,jjp1,zstd,'zstd')
233          CALL dump2d(iip1,jjp1,zsig,'zsig')
234          CALL dump2d(iip1,jjp1,zgam,'zgam')
235          CALL dump2d(iip1,jjp1,zthe,'zthe')
236
237      endif
238
239c-----------------------------------------------------------------------
240c   Passage de zdata en grille (imdp1 jmdp1)
241c-----------------------------------------------------------------------
242      do j=1,jmdp1
243          do i=1,imd
244              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmxgdat*(j-1))
245          enddo
246          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmxgdat*(j-1))
247      enddo
248
249c-----------------------------------------------------------------------
250c    Interpolation
251c-----------------------------------------------------------------------
252      call interp_horiz(zdataS,pfield,imd,jmd,
253     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
254
255c-----------------------------------------------------------------------
256c    Periodicite   
257c-----------------------------------------------------------------------
258
259      do j=1,jjp1
260         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
261      enddo
262 
263c-----------------------------------------------------------------------
264c    Sauvegarde des champs   
265c-----------------------------------------------------------------------
266
267      if (k.eq.1) then                    ! albedo
268         do i=1,iimp1*jjp1
269              alb(i) = pfield(i)
270          enddo
271      elseif (k.eq.2) then                ! thermal
272         do i=1,iimp1*jjp1
273              ith(i) = pfield(i)
274          enddo
275      elseif (k.eq.3) then                ! relief
276        if (relief.eq.'pla') then
277              call initial0(iimp1*jjp1,phisinit)
278        else
279             do i=1,iimp1*jjp1
280                  phisinit(i) = pfield(i)
281              enddo
282        endif
283      endif
284
285      ENDDO
286
287c-----------------------------------------------------------------------
288c    Traitement Phisinit
289c-----------------------------------------------------------------------
290
291      DO i=1,iimp1*jjp1
292            phisinit(i)=1000.*phisinit(i)
293      ENDDO
294      CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
295      CALL multscal(iimp1*jjp1,phisinit,g,phisinit)
296
297c-----------------------------------------------------------------------
298c    FIN
299c-----------------------------------------------------------------------
300
301      END
Note: See TracBrowser for help on using the repository browser.