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

Last change on this file since 1930 was 1918, checked in by emillour, 8 years ago

Mars GCM:
Code cleanup:

  • remove "comorbit.h" since it is no longer used.
  • turn "datafile.h" into module datafile_mod.F90 (and rename variable "datafile" as "datadir" since it stores the path to the datafile directory).

EM

File size: 9.7 KB
Line 
1c=======================================================================
2      SUBROUTINE datareadnc(relief,phisinit,alb,ith,z0,
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 ioipsl_getincom, only: getin
45      USE comconst_mod, ONLY: g,pi
46      use datafile_mod, only: datadir
47
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=360,jmd=179,imdp1=361,jmdp1=180)
61
62      INTEGER    iimp1
63      parameter    (iimp1=iim+1-1/iim)
64
65! Arguments:
66      CHARACTER(len=3),intent(inout) :: relief
67      REAL,intent(out) :: phisinit(iimp1*jjp1)
68      REAL,intent(out) :: alb(iimp1*jjp1)
69      REAL,intent(out) :: ith(iimp1*jjp1)
70      REAL,intent(out) :: z0(iimp1*jjp1)
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! Local variables:
78      REAL        zdata(imd*jmdp1)
79      REAL        zdataS(imdp1*jmdp1)
80      REAL        pfield(iimp1*jjp1)
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(0: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      datadir="/u/lmdz/WWW/planets/mars/datadir" ! default path to surface.nc
115      call getin("datadir",datadir) ! but users may specify another path
116     
117      ierr = NF_OPEN (trim(datadir)//'/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(datadir),'/'
123        write(*,*)'1) You can set this path in the
124     & callphys.def file:'
125        write(*,*)'   datadir=/path/to/the/datafiles'
126        write(*,*)'2) If necessary, surface.nc (and other datafiles)'
127        write(*,*)'   can be obtained online on:'
128        write(*,*)'http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
129        CALL ABORT
130      ENDIF
131
132c
133c Lecture des latitudes (coordonnees):
134c
135      ierr = NF_INQ_VARID (unit, "latitude", nvarid)
136#ifdef NC_DOUBLE
137      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, latitude)
138#else
139      ierr = NF_GET_VAR_REAL(unit, nvarid, latitude)
140#endif
141c
142c Lecture des longitudes (coordonnees):
143c
144      ierr = NF_INQ_VARID (unit, "longitude", nvarid)
145#ifdef NC_DOUBLE
146      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, longitude)
147#else
148      ierr = NF_GET_VAR_REAL(unit, nvarid, longitude)
149#endif
150
151c-----------------------------------------------------------------------
152c    Passage au format boites scalaires
153c-----------------------------------------------------------------------
154
155c-----------------------------------------------------------------------
156c       longitude(imd)        -->      rlonud(imdp1)
157c-----------------------------------------------------------------------
158
159c Passage en coordonnees boites scalaires et en radian
160      do i=1,imd
161          rlonud(i)=(longitude(i)+.5)*pi/180.
162      enddo
163
164c Repetition de la valeur im+1
165      rlonud(imdp1)=rlonud(1) + 2*pi
166
167c-----------------------------------------------------------------------
168c        latitude(jmdp1)         -->        rlonvd(jmd)
169c-----------------------------------------------------------------------
170
171c Passage en coordonnees boites scalaires et en radian
172      do j=1,jmd
173          rlatvd(j)=(latitude(j)-.5)*pi/180.
174      enddo
175
176c=======================================================================
177c   lecture NetCDF de albedo, thermal, relief, zdtm (pour francois Lott)
178c=======================================================================
179
180      string(0) = 'z0'
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=0,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      if (ierr.ne.nf_noerr) then
210        write(*,*) 'datareadnc error, cannot find ',trim(string(k))
211        write(*,*) ' in file ',trim(datadir),'/surface.nc'
212        stop
213      endif
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      if (ierr.ne.nf_noerr) then
220        write(*,*) 'datareadnc: error failed loading ',trim(string(k))
221        stop
222      endif
223
224c-----------------------------------------------------------------------
225c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
226c-----------------------------------------------------------------------
227      if (k.eq.4) then
228
229          zdata(:)=1000.*zdata(:)
230          longitude(:)=(pi/180.)*longitude(:)
231          latitude(:)=(pi/180.)*latitude(:)
232
233          call grid_noro1(360, 180, longitude, latitude, zdata,
234     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
235
236      endif
237
238c-----------------------------------------------------------------------
239c   Passage de zdata en grille (imdp1 jmdp1)
240c-----------------------------------------------------------------------
241      do j=1,jmdp1
242          do i=1,imd
243              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmxgdat*(j-1))
244          enddo
245          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmxgdat*(j-1))
246      enddo
247
248c-----------------------------------------------------------------------
249c    Interpolation
250c-----------------------------------------------------------------------
251      call interp_horiz(zdataS,pfield,imd,jmd,
252     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
253
254c-----------------------------------------------------------------------
255c    Periodicite   
256c-----------------------------------------------------------------------
257
258      do j=1,jjp1
259         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
260      enddo
261 
262c-----------------------------------------------------------------------
263c    Sauvegarde des champs   
264c-----------------------------------------------------------------------
265
266      if (k.eq.0) then                    ! z0
267         z0(1:iimp1*jjp1)=pfield(1:iimp1*jjp1)*.01
268         ! multiplied by 0.01 to have z0 in m
269      elseif (k.eq.1) then                    ! albedo
270         do i=1,iimp1*jjp1
271              alb(i) = pfield(i)
272          enddo
273      elseif (k.eq.2) then                ! thermal
274         do i=1,iimp1*jjp1
275              ith(i) = pfield(i)
276          enddo
277      elseif (k.eq.3) then                ! relief
278        if (relief.eq.'pla') then
279              call initial0(iimp1*jjp1,phisinit)
280        else
281             do i=1,iimp1*jjp1
282                  phisinit(i) = pfield(i)
283              enddo
284        endif
285      endif
286
287      ENDDO
288
289c-----------------------------------------------------------------------
290c    Traitement Phisinit
291c-----------------------------------------------------------------------
292
293      phisinit(1:iimp1*jjp1)=1000.*phisinit(1:iimp1*jjp1)
294      phisinit(:)=g*phisinit(:)
295
296c-----------------------------------------------------------------------
297c    FIN
298c-----------------------------------------------------------------------
299
300      END
Note: See TracBrowser for help on using the repository browser.