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

Last change on this file since 304 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

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      implicit none
45
46#include "dimensions.h"
47#include "paramet.h"
48#include "comgeom.h"
49#include "comconst.h"
50#include "netcdf.inc"
51#include "datafile.h"
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 '//datafile(1:lnblnk(datafile))//'/*.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 (datafile(1:lnblnk(datafile))//'/surface.nc',
123!     &  NF_NOWRITE,unit)
124      ierr = NF_OPEN (datafile(1:lnblnk(datafile))//'/'//filename,
125     &  NF_NOWRITE,unit)
126      IF (ierr.NE.NF_NOERR) THEN
127        write(*,*)'Error : cannot open file '//filename
128        write(*,*)'(in phystd/datareadnc.F)'
129        write(*,*)'It should be in :',datafile(1:lnblnk(datafile)),'/'
130        write(*,*)'1) You can change this directory address in '
131        write(*,*)'   file phystd/datafile.h'
132        write(*,*)'2) If necessary surface.nc (and other datafiles)'
133        write(*,*)'   can be obtained online on:'
134        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
135        CALL ABORT
136      ENDIF
137
138c
139c Lecture des latitudes (coordonnees):
140c
141      ierr = NF_INQ_VARID (unit, "latitude", nvarid)
142#ifdef NC_DOUBLE
143      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, latitude)
144#else
145      ierr = NF_GET_VAR_REAL(unit, nvarid, latitude)
146#endif
147c
148c Lecture des longitudes (coordonnees):
149c
150      ierr = NF_INQ_VARID (unit, "longitude", nvarid)
151#ifdef NC_DOUBLE
152      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, longitude)
153#else
154      ierr = NF_GET_VAR_REAL(unit, nvarid, longitude)
155#endif
156
157c-----------------------------------------------------------------------
158c    Passage au format boites scalaires
159c-----------------------------------------------------------------------
160
161c-----------------------------------------------------------------------
162c       longitude(imd)        -->      rlonud(imdp1)
163c-----------------------------------------------------------------------
164
165c Passage en coordonnees boites scalaires et en radian
166      do i=1,imd
167          rlonud(i)=(longitude(i)+.5)*pi/180.
168      enddo
169
170c Repetition de la valeur im+1
171      rlonud(imdp1)=rlonud(1) + 2*pi
172
173c-----------------------------------------------------------------------
174c        latitude(jmdp1)         -->        rlonvd(jmd)
175c-----------------------------------------------------------------------
176
177c Passage en coordonnees boites scalaires et en radian
178      do j=1,jmd
179          rlatvd(j)=(latitude(j)-.5)*pi/180.
180      enddo
181
182c=======================================================================
183c   lecture NetCDF de albedo, thermal, relief, zdtm (pour francois Lott)
184c=======================================================================
185
186      string(1) = 'albedo'
187      string(2) = 'thermal'
188      if (relief.ne.'pla') then
189        write(*,*) ' La topographie est celle de MOLA'
190        relief = 'MOL'
191          string(3) = 'z'//relief
192      else
193          string(3) = 'zMOL'  ! pour qu''il lise qqchose sur le fichier
194                            ! remise a 0 derriere
195      endif
196      string(4) = 'zMOL'    ! lecture pour calcul topog. sous-maille
197 
198
199      DO k=1,4
200          write(*,*) 'string',k,string(k)
201         
202c-----------------------------------------------------------------------
203c    initialisation
204c-----------------------------------------------------------------------
205      call initial0(iimp1*jjp1,pfield)
206      call initial0(imd*jmdp1,zdata)
207      call initial0(imdp1*jmdp1,zdataS)
208
209c-----------------------------------------------------------------------
210c    Lecture NetCDF 
211c-----------------------------------------------------------------------
212
213      ierr = NF_INQ_VARID (unit, string(k), nvarid)
214#ifdef NC_DOUBLE
215      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, zdata)
216#else
217      ierr = NF_GET_VAR_REAL(unit, nvarid, zdata)
218#endif
219
220c-----------------------------------------------------------------------
221c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
222c-----------------------------------------------------------------------
223      if (k.eq.4) then
224
225          call multscal(imd*jmdp1,zdata,1000.,zdata)
226          call multscal(imd,longitude,pi/180.,longitude)
227          call multscal(jmdp1,latitude,pi/180.,latitude)
228
229          call grid_noro1(360, 180, longitude, latitude, zdata,
230     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
231
232          CALL dump2d(iip1,jjp1,zmea,'zmea')
233          CALL dump2d(iip1,jjp1,zstd,'zstd')
234          CALL dump2d(iip1,jjp1,zsig,'zsig')
235          CALL dump2d(iip1,jjp1,zgam,'zgam')
236          CALL dump2d(iip1,jjp1,zthe,'zthe')
237
238      endif
239
240c-----------------------------------------------------------------------
241c   Passage de zdata en grille (imdp1 jmdp1)
242c-----------------------------------------------------------------------
243      do j=1,jmdp1
244          do i=1,imd
245              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmxgdat*(j-1))
246          enddo
247          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmxgdat*(j-1))
248      enddo
249
250c-----------------------------------------------------------------------
251c    Interpolation
252c-----------------------------------------------------------------------
253      call interp_horiz(zdataS,pfield,imd,jmd,
254     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
255
256c-----------------------------------------------------------------------
257c    Periodicite   
258c-----------------------------------------------------------------------
259
260      do j=1,jjp1
261         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
262      enddo
263 
264c-----------------------------------------------------------------------
265c    Sauvegarde des champs   
266c-----------------------------------------------------------------------
267
268      if (k.eq.1) then                    ! albedo
269         do i=1,iimp1*jjp1
270              alb(i) = pfield(i)
271          enddo
272      elseif (k.eq.2) then                ! thermal
273         do i=1,iimp1*jjp1
274              ith(i) = pfield(i)
275          enddo
276      elseif (k.eq.3) then                ! relief
277        if (relief.eq.'pla') then
278              call initial0(iimp1*jjp1,phisinit)
279        else
280             do i=1,iimp1*jjp1
281                  phisinit(i) = pfield(i)
282              enddo
283        endif
284      endif
285
286      ENDDO
287
288c-----------------------------------------------------------------------
289c    Traitement Phisinit
290c-----------------------------------------------------------------------
291
292      DO i=1,iimp1*jjp1
293            phisinit(i)=1000.*phisinit(i)
294      ENDDO
295      CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
296      CALL multscal(iimp1*jjp1,phisinit,g,phisinit)
297
298c-----------------------------------------------------------------------
299c    FIN
300c-----------------------------------------------------------------------
301
302      END
Note: See TracBrowser for help on using the repository browser.