source: trunk/LMDZ.TITAN/libf/dynphy_lonlat/phytitan/datareadnc.F @ 1815

Last change on this file since 1815 was 1470, checked in by emillour, 9 years ago

Generic GCM:

  • reorganizing the "datadir" structure: aerosol properties should now be in subdirectory 'aerosol_properties' of datadir, and surface.nc files should be in subdirectory 'surface_data'. These subdirectory names are stored in module datafile_mod.
  • Made things retro-compatible so that using an 'old' datadir structure (ie. aerosol properties files and surface files in datadir) still works.

EM

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=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!#include "fxyprim.h"
97
98      pi=2.*ASIN(1.)
99
100c=======================================================================
101c    rlonud, rlatvd
102c=======================================================================
103
104c-----------------------------------------------------------------------
105c    Lecture NetCDF des donnees latitude et longitude
106c-----------------------------------------------------------------------
107      ierr = NF_OPEN (trim(datadir)//'/'//trim(surfdir)//'/'//
108     &                trim(adjustl(filename)),
109     &                NF_NOWRITE,unit)
110      IF (ierr.NE.NF_NOERR) THEN
111        ! In ye old days this file was stored in datadir;
112        ! let's be retro-compatible
113        ierr = NF_OPEN (trim(datadir)//'/'//
114     &                trim(adjustl(filename)),
115     &                NF_NOWRITE,unit)
116       
117      ENDIF
118      IF (ierr.NE.NF_NOERR) THEN
119        write(*,*)'Error : cannot open file '//trim(filename)
120        write(*,*)'(in phystd/datareadnc.F)'
121        write(*,*)'It should be in :',trim(datadir),'/',trim(surfdir)
122        write(*,*)'Check that your path to datagcm:',trim(datadir)
123        write(*,*)' is correct. You can change it in callphys.def with:'
124        write(*,*)' datadir = /absolute/path/to/datagcm'
125        write(*,*)'If necessary surface.nc (and other datafiles)'
126        write(*,*)' can be obtained online on:'
127        write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/'//
128     &             'LMDZ.GENERIC/datagcm/'
129        STOP
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(1) = 'albedo'
181      string(2) = 'thermal'
182      if (relief.ne.'pla') then
183        write(*,*) ' La topographie est celle de MOLA'
184        relief = 'MOL'
185          string(3) = 'z'//relief
186      else
187          string(3) = 'zMOL'  ! pour qu''il lise qqchose sur le fichier
188                            ! remise a 0 derriere
189      endif
190      string(4) = 'zMOL'    ! lecture pour calcul topog. sous-maille
191 
192
193      DO k=1,4
194          write(*,*) 'string',k,string(k)
195         
196c-----------------------------------------------------------------------
197c    initialisation
198c-----------------------------------------------------------------------
199      pfield(1:iimp1*jjp1)=0
200      zdata(1:imd*jmdp1)=0
201      zdataS(1:iimp1*jjp1)=0
202
203c-----------------------------------------------------------------------
204c    Lecture NetCDF 
205c-----------------------------------------------------------------------
206
207      ierr = NF_INQ_VARID (unit, string(k), nvarid)
208#ifdef NC_DOUBLE
209      ierr = NF_GET_VAR_DOUBLE(unit, nvarid, zdata)
210#else
211      ierr = NF_GET_VAR_REAL(unit, nvarid, zdata)
212#endif
213
214c-----------------------------------------------------------------------
215c        Cas particulier "Francois Lott" ( k=4 ) (relief sous-maille)
216c-----------------------------------------------------------------------
217      if (k.eq.4) then
218
219          zdata(:)=1000.*zdata(:)
220          longitude(:)=(pi/180.)*longitude(:)
221          latitude(:)=(pi/180.)*latitude(:)
222
223          call grid_noro1(360, 180, longitude, latitude, zdata,
224     .         iim, jjp1, rlonv, rlatu, zmea,zstd,zsig,zgam,zthe)
225
226          !CALL dump2d(iip1,jjp1,zmea,'zmea')
227          !CALL dump2d(iip1,jjp1,zstd,'zstd')
228          !CALL dump2d(iip1,jjp1,zsig,'zsig')
229          !CALL dump2d(iip1,jjp1,zgam,'zgam')
230          !CALL dump2d(iip1,jjp1,zthe,'zthe')
231
232      endif
233
234c-----------------------------------------------------------------------
235c   Passage de zdata en grille (imdp1 jmdp1)
236c-----------------------------------------------------------------------
237      do j=1,jmdp1
238          do i=1,imd
239              zdataS(i+imdp1*(j-1)) = zdata(i+ngridmixgdat*(j-1))
240          enddo
241          zdataS(imdp1+imdp1*(j-1)) = zdata(1+ngridmixgdat*(j-1))
242      enddo
243
244c-----------------------------------------------------------------------
245c    Interpolation
246c-----------------------------------------------------------------------
247      call interp_horiz(zdataS,pfield,imd,jmd,
248     .    iim, jjm,1,rlonud,rlatvd,rlonu,rlatv)
249
250c-----------------------------------------------------------------------
251c    Periodicite   
252c-----------------------------------------------------------------------
253
254      do j=1,jjp1
255         pfield(iimp1*j) =  pfield(1+iimp1*(j-1))
256      enddo
257 
258c-----------------------------------------------------------------------
259c    Sauvegarde des champs   
260c-----------------------------------------------------------------------
261
262      if (k.eq.1) then                    ! albedo
263         do i=1,iimp1*jjp1
264              alb(i) = pfield(i)
265          enddo
266      elseif (k.eq.2) then                ! thermal
267         do i=1,iimp1*jjp1
268              ith(i) = pfield(i)
269          enddo
270      elseif (k.eq.3) then                ! relief
271        if (relief.eq.'pla') then
272             phisinit(1:iimp1*jjp1)=0
273        else
274             phisinit(1:iimp1*jjp1)=pfield(1:iimp1*jjp1)
275        endif
276      endif
277
278      ENDDO
279
280c-----------------------------------------------------------------------
281c    Traitement Phisinit
282c-----------------------------------------------------------------------
283
284      phisinit(1:iimp1*jjp1)=1000.*phisinit(1:iimp1*jjp1)
285      !CALL dump2d(iimp1,jjp1,phisinit,'Altitude en m')
286      phisinit(:)=g*phisinit(:)
287
288c-----------------------------------------------------------------------
289c    FIN
290c-----------------------------------------------------------------------
291
292      END
Note: See TracBrowser for help on using the repository browser.