source: trunk/LMDZ.MARS/libf/phymars/datareadnc.F @ 171

Last change on this file since 171 was 164, checked in by emillour, 13 years ago

Mars GCM:

Updates and corrections (to enable compiling/running in debug mode with

ifort)

  • removed option "-free-form" from makegcm_ifort and set mod_loc_dir="." so that module files (produced in local directory by ifort) are moved to LIBO
  • updated initracer.F, physdem1.F, physiq.F, inichim_newstart.F to avoid referencing out-of-bound array indexes (even if unused)
  • cosmetic updates on inwrite.F, datareadnc.F
  • updated newstart.F to initialize and use 'datadir' when looking for files
  • corrected bug on interpolation of sub-surface temperatures in lect_start_archive.F

EM

File size: 9.5 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! to use  'getin'
45       use ioipsl_getincom
46
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#include "datafile.h"
55
56c=======================================================================
57c   Declarations:
58C=======================================================================
59
60      INTEGER    imd,jmd,imdp1,jmdp1
61      parameter    (imd=360,jmd=179,imdp1=361,jmdp1=180)
62
63      INTEGER    iimp1
64      parameter    (iimp1=iim+1-1/iim)
65
66      CHARACTER relief*3
67
68      REAL        zdata(imd*jmdp1)
69      REAL        zdataS(imdp1*jmdp1)
70      REAL        pfield(iimp1*jjp1)
71
72      REAL        alb(iimp1*jjp1)
73      REAL        ith(iimp1*jjp1)
74      REAL        phisinit(iimp1*jjp1)
75
76      REAL        zmea(imdp1*jmdp1)
77      REAL        zstd(imdp1*jmdp1)
78      REAL        zsig(imdp1*jmdp1)
79      REAL        zgam(imdp1*jmdp1)
80      REAL        zthe(imdp1*jmdp1)
81
82      INTEGER     ierr
83
84      INTEGER   unit,nvarid
85
86      INTEGER    i,j,k
87
88      INTEGER klatdat,ngridmxgdat
89      PARAMETER (klatdat=180,ngridmxgdat=360)
90
91c    on passe une grille en rlonu rlatv et im+1 jm a interp_horiz)
92
93      REAL longitude(imd),latitude(jmdp1) ! Pour lecture des donnees
94      REAL rlonud(imdp1),rlatvd(jmd)
95
96      CHARACTER*20 string
97      DIMENSION string(4)
98
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      write(*,*) 'datareadnc: opening file surface.nc'
113
114      datafile="/u/forget/WWW/datagcm/datafile" ! default path to surface.nc
115      call getin("datadir",datafile) ! but users may specify another path
116     
117      ierr = NF_OPEN (trim(datafile)//'/surface.nc',
118     &  NF_NOWRITE,unit)
119      IF (ierr.NE.NF_NOERR) THEN
120        write(*,*)'Error : cannot open file surface.nc '
121        write(*,*)'(in phymars/datareadnc.F)'
122        write(*,*)'It should be in :',trim(datafile),'/'
123        write(*,*)'1) You can set this path in the callphys.def file:'
124        write(*,*)'   datadir=/path/to/the/datafiles'
125        write(*,*)'   OR specify the path to use in callphys.def, as:'
126        write(*,*)'   datadir=/path/to/the/directory'
127        write(*,*)'2) If necessary, surface.nc (and other datafiles)'
128        write(*,*)'   can be obtained online on:'
129        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
130        CALL ABORT
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(*,*) ' MOLA topography'
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      call initial0(iimp1*jjp1,pfield)
201      call initial0(imd*jmdp1,zdata)
202      call initial0(imdp1*jmdp1,zdataS)
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          call multscal(imd*jmdp1,zdata,1000.,zdata)
221          call multscal(imd,longitude,pi/180.,longitude)
222          call multscal(jmdp1,latitude,pi/180.,latitude)
223
224          call grid_noro1(360, 180, 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+ngridmxgdat*(j-1))
241          enddo
242          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmxgdat*(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              call initial0(iimp1*jjp1,phisinit)
274        else
275             do i=1,iimp1*jjp1
276                  phisinit(i) = pfield(i)
277              enddo
278        endif
279      endif
280
281      ENDDO
282
283c-----------------------------------------------------------------------
284c    Traitement Phisinit
285c-----------------------------------------------------------------------
286
287      DO i=1,iimp1*jjp1
288            phisinit(i)=1000.*phisinit(i)
289      ENDDO
290      CALL dump2d(iimp1,jjp1,phisinit,'Altitude in m')
291      CALL multscal(iimp1*jjp1,phisinit,g,phisinit)
292
293c-----------------------------------------------------------------------
294c    FIN
295c-----------------------------------------------------------------------
296
297      END
Note: See TracBrowser for help on using the repository browser.