source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/datareadnc.F @ 222

Last change on this file since 222 was 222, checked in by ymipsl, 10 years ago

Creating temporary dynamico/lmdz/saturn branche

YM

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