source: trunk/LMDZ.PLUTO/libf/dynphy_lonlat/phypluto/datareadnc.F @ 3945

Last change on this file since 3945 was 3878, checked in by tbertrand, 6 months ago

LMDZ.PLUTO:
For newstart : reading in the intial albedo map and the new topography map
TB

File size: 9.7 KB
Line 
1c=======================================================================
2      SUBROUTINE datareadnc(relief,filename,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, surfdir
45! to use  'getin'
46      USE ioipsl_getincom
47      USE comconst_mod, ONLY: g,pi
48      implicit none
49
50#include "dimensions.h"
51#include "paramet.h"
52#include "comgeom.h"
53#include "netcdf.inc"
54
55c=======================================================================
56c   Declarations:
57C=======================================================================
58
59      INTEGER    imd,jmd,imdp1,jmdp1
60      !parameter    (imd=1440,jmd=719,imdp1=1441,jmdp1=720)
61      parameter    (imd=360,jmd=179,imdp1=361,jmdp1=180)
62
63      INTEGER    iimp1
64      parameter    (iimp1=iim+1-1/iim)
65
66      character(len=3),intent(inout) :: relief*3
67      character(len=*),intent(in) :: filename ! surface.nc file
68      real,intent(out) :: phisinit(iimp1*jjp1) ! surface geopotential
69      real,intent(out) :: alb(iimp1*jjp1) ! albedo
70      real,intent(out) :: ith(iimp1*jjp1) ! thermal inertia
71      real,intent(out) :: zmea(imdp1*jmdp1)
72      real,intent(out) :: zstd(imdp1*jmdp1)
73      real,intent(out) :: zsig(imdp1*jmdp1)
74      real,intent(out) :: zgam(imdp1*jmdp1)
75      real,intent(out) :: zthe(imdp1*jmdp1)
76     
77      REAL        zdata(imd*jmdp1)
78      REAL        zdataS(imdp1*jmdp1)
79      REAL        pfield(iimp1*jjp1)
80
81      INTEGER   ierr
82
83      INTEGER   unit,nvarid
84
85      INTEGER    i,j,k
86
87      INTEGER klatdat,ngridmixgdat
88      PARAMETER (klatdat=jmd+1,ngridmixgdat=imd)
89
90c    on passe une grille en rlonu rlatv et im+1 jm a interp_horiz)
91
92      REAL longitude(imd),latitude(jmdp1) ! Pour lecture des donnees
93      REAL rlonud(imdp1),rlatvd(jmd)
94
95      CHARACTER*20 string
96      DIMENSION string(4)
97!#include "fxyprim.h"
98
99      pi=2.*ASIN(1.)
100
101c=======================================================================
102c    rlonud, rlatvd
103c=======================================================================
104
105c-----------------------------------------------------------------------
106c    Lecture NetCDF des donnees latitude et longitude
107c-----------------------------------------------------------------------
108      ierr = NF_OPEN (trim(datadir)//'/'//trim(surfdir)//'/'//
109     &                trim(adjustl(filename)),
110     &                NF_NOWRITE,unit)
111      IF (ierr.NE.NF_NOERR) THEN
112        ! In ye old days this file was stored in datadir;
113        ! let's be retro-compatible
114        ierr = NF_OPEN (trim(datadir)//'/'//
115     &                trim(adjustl(filename)),
116     &                NF_NOWRITE,unit)
117       
118      ENDIF
119      IF (ierr.NE.NF_NOERR) THEN
120        write(*,*)'Error : cannot open file '//trim(filename)
121        write(*,*)'(in phystd/datareadnc.F)'
122        write(*,*)'It should be in :',trim(datadir),'/',trim(surfdir)
123        write(*,*)'Check that your path to datagcm:',trim(datadir)
124        write(*,*)' is correct. You can change it in callphys.def with:'
125        write(*,*)' datadir = /absolute/path/to/datagcm'
126        write(*,*)'If necessary surface.nc (and other datafiles)'
127        write(*,*)' can be obtained online on:'
128        write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/'//
129     &             'LMDZ.GENERIC/datagcm/'
130        STOP
131      ENDIF
132
133c
134c Lecture des latitudes (coordonnees):
135c
136      ierr = NF_INQ_VARID (unit, "latitude", nvarid)
137#ifdef NC_DOUBLE
138      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, latitude)
139#else
140      ierr = NF_GET_VAR_REAL(unit, nvarid, latitude)
141#endif
142c
143c Lecture des longitudes (coordonnees):
144c
145      ierr = NF_INQ_VARID (unit, "longitude", nvarid)
146#ifdef NC_DOUBLE
147      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, longitude)
148#else
149      ierr = NF_GET_VAR_REAL(unit, nvarid, longitude)
150#endif
151
152c-----------------------------------------------------------------------
153c    Passage au format boites scalaires
154c-----------------------------------------------------------------------
155
156c-----------------------------------------------------------------------
157c       longitude(imd)        -->      rlonud(imdp1)
158c-----------------------------------------------------------------------
159
160c Passage en coordonnees boites scalaires et en radian
161      do i=1,imd
162          rlonud(i)=(longitude(i)+.5)*pi/180.
163      enddo
164
165c Repetition de la valeur im+1
166      rlonud(imdp1)=rlonud(1) + 2*pi
167
168c-----------------------------------------------------------------------
169c        latitude(jmdp1)         -->        rlonvd(jmd)
170c-----------------------------------------------------------------------
171
172c Passage en coordonnees boites scalaires et en radian
173      do j=1,jmd
174          rlatvd(j)=(latitude(j)-.5)*pi/180.
175      enddo
176
177c=======================================================================
178c   lecture NetCDF de albedo, thermal, relief, zdtm (pour francois Lott)
179c=======================================================================
180
181      string(1) = 'albedo'
182      string(2) = 'thermal'
183      if (relief.ne.'pla') then
184        write(*,*) ' La topographie est celle de MOLA'
185        relief = 'MOL'
186          string(3) = 'z'//relief
187      else
188          string(3) = 'zMOL'  ! pour qu''il lise qqchose sur le fichier
189                            ! remise a 0 derriere
190      endif
191      string(4) = 'zMOL'    ! lecture pour calcul topog. sous-maille
192 
193
194      DO k=1,4
195          write(*,*) 'string',k,string(k)
196         
197c-----------------------------------------------------------------------
198c    initialisation
199c-----------------------------------------------------------------------
200      pfield(1:iimp1*jjp1)=0
201      zdata(1:imd*jmdp1)=0
202      zdataS(1:iimp1*jjp1)=0
203
204c-----------------------------------------------------------------------
205c    Lecture NetCDF 
206c-----------------------------------------------------------------------
207
208      ierr = NF_INQ_VARID (unit, string(k), nvarid)
209#ifdef NC_DOUBLE
210      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, zdata)
211#else
212      ierr = NF_GET_VAR_REAL(unit, nvarid, zdata)
213#endif
214
215c-----------------------------------------------------------------------
216c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
217c-----------------------------------------------------------------------
218      if (k.eq.4) then
219
220          zdata(:)=1000.*zdata(:)
221          longitude(:)=(pi/180.)*longitude(:)
222          latitude(:)=(pi/180.)*latitude(:)
223
224          call grid_noro1(imd,jmd+1, longitude, latitude, zdata,
225     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
226
227          !CALL dump2d(iip1,jjp1,zmea,'zmea')
228          !CALL dump2d(iip1,jjp1,zstd,'zstd')
229          !CALL dump2d(iip1,jjp1,zsig,'zsig')
230          !CALL dump2d(iip1,jjp1,zgam,'zgam')
231          !CALL dump2d(iip1,jjp1,zthe,'zthe')
232
233      endif
234
235c-----------------------------------------------------------------------
236c   Passage de zdata en grille (imdp1 jmdp1)
237c-----------------------------------------------------------------------
238      do j=1,jmdp1
239          do i=1,imd
240              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmixgdat*(j-1))
241          enddo
242          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmixgdat*(j-1))
243      enddo
244
245c-----------------------------------------------------------------------
246c    Interpolation
247c-----------------------------------------------------------------------
248      call interp_horiz(zdataS,pfield,imd,jmd,
249     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
250
251c-----------------------------------------------------------------------
252c    Periodicite   
253c-----------------------------------------------------------------------
254
255      do j=1,jjp1
256         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
257      enddo
258 
259c-----------------------------------------------------------------------
260c    Sauvegarde des champs   
261c-----------------------------------------------------------------------
262
263      if (k.eq.1) then                    ! albedo
264         do i=1,iimp1*jjp1
265              alb(i) = pfield(i)
266          enddo
267      elseif (k.eq.2) then                ! thermal
268         do i=1,iimp1*jjp1
269              ith(i) = pfield(i)
270          enddo
271      elseif (k.eq.3) then                ! relief
272        if (relief.eq.'pla') then
273             phisinit(1:iimp1*jjp1)=0
274        else
275             phisinit(1:iimp1*jjp1)=pfield(1:iimp1*jjp1)
276        endif
277      endif
278
279      ENDDO
280
281c-----------------------------------------------------------------------
282c    Traitement Phisinit
283c-----------------------------------------------------------------------
284
285      phisinit(1:iimp1*jjp1)=1000.*phisinit(1:iimp1*jjp1)
286      !CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
287      phisinit(:)=g*phisinit(:)
288
289c-----------------------------------------------------------------------
290c    FIN
291c-----------------------------------------------------------------------
292
293      END
Note: See TracBrowser for help on using the repository browser.