source: trunk/LMDZ.PLUTO.old/libf/phypluto/datareadnc.F @ 3436

Last change on this file since 3436 was 3175, checked in by emillour, 11 months ago

Pluto PCM:
Add the old Pluto LMDZ for reference (required prior step to making
an LMDZ.PLUTO using the same framework as the other physics packages).
TB+EM

File size: 9.8 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 Plutonian surface to use in a GCM
11c   ------                from NetCDF file "surface.nc"
12c
13c   Arguments:
14c   ----------
15c
16c     Inputs:
17c     ------
18c
19c     Outputs:
20c     --------
21c
22c=======================================================================
23c   donnees ALBEDO, INERTIE THERMIQUE, RELIEF:
24c
25c       Ces donnees sont au format NetCDF dans le fichier "surface.nc"
26c
27c   360 valeurs en longitude (de -179.5 a 179.5)
28c   180 valeurs en latitudes (de 89.5 a -89.5)
29c
30c   Pour les passer au format de la grille, on utilise "interp_horiz.F"
31c
32c   Il faut donc que ces donnees soient au format grille scalaire
33c               (imold+1 jmold+1)
34c       avec passage des coordonnees de la "boite" (rlonu, rlatv)
35c
36c   On prend imd (d pour donnees!)
37c           imd = 360 avec copie de la 1ere valeur sur la imd+1
38c                   (rlonud de -179 a -181)
39c           jmd = 179
40c                   (rlatvd de 89 a -89)
41c=======================================================================
42
43      use datafile_mod, only: datadir2
44
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=1440,jmd=719,imdp1=1441,jmdp1=720)
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=720,ngridmxgdat=1440)
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*40 filename
98      CHARACTER*250 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     
115      filestring='ls '//TRIM(datadir2)//
116     &                                       '/surface_data/*.nc'
117      call system(filestring)
118
119      write(*,*) ''
120      write(*,*) 'Please choose the relevant file:'
121      write(*,*) '(e.g. type "surface_pluto.nc")'
122      read(*,fmt='(a40)') filename
123
124      ierr = NF_OPEN (TRIM(datadir2)//
125     &                                     '/surface_data/'//filename,
126     &  NF_NOWRITE,unit)
127      IF (ierr.NE.NF_NOERR) THEN
128        write(*,*)'Error : cannot open file '//filename
129        write(*,*)'(in phypluto/datareadnc.F)'
130        write(*,*)'It should be in:',TRIM(datadir2),
131     &                                                '/surface_data/'
132        write(*,*)'1) You can change this directory address in '
133        write(*,*)'   file phypluto/datafile_mod or callphys.def'
134        write(*,*)'2) If necessary surface.nc (and other datafiles)'
135        write(*,*)'   can be obtained in datagcm/surface_data'
136        CALL ABORT
137      ENDIF
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      write(*,*)' TBlongitude',longitude
157
158c-----------------------------------------------------------------------
159c    Passage au format boites scalaires
160c-----------------------------------------------------------------------
161
162c-----------------------------------------------------------------------
163c       longitude(imd)        -->      rlonud(imdp1)
164c-----------------------------------------------------------------------
165
166c Passage en coordonnees boites scalaires et en radian
167      do i=1,imd
168          rlonud(i)=(longitude(i)+.5)*pi/180.
169      enddo
170
171c Repetition de la valeur im+1
172      rlonud(imdp1)=rlonud(1) + 2*pi
173
174c-----------------------------------------------------------------------
175c        latitude(jmdp1)         -->        rlonvd(jmd)
176c-----------------------------------------------------------------------
177
178c Passage en coordonnees boites scalaires et en radian
179      do j=1,jmd
180          rlatvd(j)=(latitude(j)-.5)*pi/180.
181      enddo
182
183c=======================================================================
184c   lecture NetCDF de albedo, thermal, relief, zdtm (pour francois Lott)
185c=======================================================================
186
187      string(1) = 'albedo'
188      string(2) = 'thermal'
189      !if (relief.ne.'pla') then
190      !  write(*,*) ' La topographie est celle de MOLA'
191      !  relief = 'MOL'
192      !    string(3) = 'z'//relief
193      !else
194      string(3) = 'zMOL'  ! pour qu''il lise qqchose sur le fichier
195                            ! remise a 0 derriere
196      !endif
197      string(4) = 'zMOL'    ! lecture pour calcul topog. sous-maille
198 
199
200      DO k=1,3
201          write(*,*) 'TB string',k,string(k)
202
203c-----------------------------------------------------------------------
204c    initialisation
205c-----------------------------------------------------------------------
206      call initial0(iimp1*jjp1,pfield)
207      call initial0(imd*jmdp1,zdata)
208      call initial0(imdp1*jmdp1,zdataS)
209          write(*,*) 'TB initi ok'
210
211c-----------------------------------------------------------------------
212c    Lecture NetCDF 
213c-----------------------------------------------------------------------
214
215      ierr = NF_INQ_VARID (unit, string(k), nvarid)
216#ifdef NC_DOUBLE
217      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, zdata)
218#else
219      ierr = NF_GET_VAR_REAL(unit, nvarid, zdata)
220#endif
221      write(*,*) 'TB lect ok'
222
223c-----------------------------------------------------------------------
224c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
225c-----------------------------------------------------------------------
226      if (k.eq.4) then
227
228          call multscal(imd*jmdp1,zdata,1000.,zdata)
229          call multscal(imd,longitude,pi/180.,longitude)
230          call multscal(jmdp1,latitude,pi/180.,latitude)
231
232          call grid_noro1(1440, 720, longitude, latitude, zdata,
233     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
234
235          CALL dump2d(iip1,jjp1,zmea,'zmea')
236          CALL dump2d(iip1,jjp1,zstd,'zstd')
237          CALL dump2d(iip1,jjp1,zsig,'zsig')
238          CALL dump2d(iip1,jjp1,zgam,'zgam')
239          CALL dump2d(iip1,jjp1,zthe,'zthe')
240
241      endif
242
243c-----------------------------------------------------------------------
244c   Passage de zdata en grille (imdp1 jmdp1)
245c-----------------------------------------------------------------------
246      do j=1,jmdp1
247          do i=1,imd
248              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmxgdat*(j-1))
249          enddo
250          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmxgdat*(j-1))
251      enddo
252      write(*,*) 'TB zdata ok'
253
254c-----------------------------------------------------------------------
255c    Interpolation
256c-----------------------------------------------------------------------
257      call interp_horiz(zdataS,pfield,imd,jmd,
258     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
259
260c-----------------------------------------------------------------------
261c    Periodicite   
262c-----------------------------------------------------------------------
263
264      do j=1,jjp1
265         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
266      enddo
267 
268c-----------------------------------------------------------------------
269c    Sauvegarde des champs   
270c-----------------------------------------------------------------------
271
272      if (k.eq.1) then                    ! albedo
273         do i=1,iimp1*jjp1
274              alb(i) = pfield(i)
275          enddo
276      elseif (k.eq.2) then                ! thermal
277         do i=1,iimp1*jjp1
278              ith(i) = pfield(i)
279          enddo
280      elseif (k.eq.3) then                ! relief
281        if (relief.eq.'pla') then
282              call initial0(iimp1*jjp1,phisinit)
283        else
284             do i=1,iimp1*jjp1
285                  phisinit(i) = pfield(i)
286              enddo
287        endif
288      endif
289
290      ENDDO
291
292c-----------------------------------------------------------------------
293c    Traitement Phisinit
294c-----------------------------------------------------------------------
295
296      DO i=1,iimp1*jjp1
297            phisinit(i)=1000.*phisinit(i)
298      ENDDO
299      CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
300      CALL multscal(iimp1*jjp1,phisinit,g,phisinit)
301
302c-----------------------------------------------------------------------
303c    FIN
304c-----------------------------------------------------------------------
305
306      END
Note: See TracBrowser for help on using the repository browser.