source: trunk/LMDZ.GENERIC/libf/phystd/phyetat0.F @ 223

Last change on this file since 223 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

File size: 17.5 KB
Line 
1      SUBROUTINE phyetat0 (fichnom,tab0,Lmodif,nsoil,nq,
2     .           day_ini,time,
3     .           tsurf,tsoil,emis,q2,qsurf)
4      implicit none
5c======================================================================
6c Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
7c  Adaptation à Mars : Yann Wanherdrick
8c Objet: Lecture de l etat initial pour la physique
9c======================================================================
10#include "netcdf.inc"
11#include "dimensions.h"
12#include "dimphys.h"
13#include "comgeomfi.h"
14#include "surfdat.h"
15#include "planete.h"
16#include "comcstfi.h"
17#include"advtrac.h"
18
19c======================================================================
20      INTEGER nbsrf !Mars nbsrf a 1 au lieu de 4
21      PARAMETER (nbsrf=1) ! nombre de sous-fractions pour une maille
22!======================================================================
23!  Arguments:
24!  ---------
25!  inputs:
26      character*(*) fichnom ! "startfi.nc" file
27      integer tab0
28      integer Lmodif
29      integer nsoil ! # of soil layers
30      integer nq
31      integer day_ini
32      real time
33
34!  outputs:
35      real tsurf(ngridmx,nbsrf) ! surface temperature
36      real tsoil(ngridmx,nsoil,nbsrf) ! soil temperature
37      real emis(ngridmx) ! surface emissivity
38      real q2(ngridmx, llm+1) !
39      real qsurf(ngridmx,nq) ! tracers on surface
40!      real co2ice(ngridmx) ! co2 ice cover
41
42!======================================================================
43!  Local variables:
44
45!      INTEGER radpas
46!      REAL co2_ppm
47!      REAL solaire
48
49      real xmin,xmax ! to display min and max of a field
50c
51      INTEGER ig,iq,lmax
52      INTEGER nid, nvarid
53      INTEGER ierr, i, nsrf
54!      integer isoil
55!      INTEGER length
56!      PARAMETER (length=100)
57      CHARACTER*7 str7
58      CHARACTER*2 str2
59      CHARACTER*1 yes
60c
61      REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec
62      INTEGER nqold
63
64! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...)
65      logical :: oldtracernames=.false.
66      integer :: count
67      character(len=30) :: txt ! to store some text
68
69 
70c
71c Ouvrir le fichier contenant l etat initial:
72c
73
74      ierr = NF_OPEN (fichnom, NF_NOWRITE,nid)
75      IF (ierr.NE.NF_NOERR) THEN
76        write(6,*)' Pb d''ouverture du fichier '//fichnom
77        CALL ABORT
78      ENDIF
79
80! Preliminary stuff: check if tracers follow old naming convention (qsurf01,
81!                    qsurf02, ...)
82      count=0
83      do iq=1,nqmx
84        txt= " "
85        write(txt,'(a5,i2.2)')'qsurf',iq
86        ierr=NF_INQ_VARID(nid,txt,nvarid)
87        if (ierr.ne.NF_NOERR) then
88          ! did not find old tracer name
89          exit ! might as well stop here
90        else
91          ! found old tracer name
92          count=count+1
93        endif
94      enddo
95      if (count.eq.nqmx) then
96        write(*,*) "phyetat0:tracers seem to follow old naming ",
97     &             "convention (qsurf01,qsurf02,...)"
98        write(*,*) "   => will work for now ... "
99        write(*,*) "      but you should run newstart to rename them"
100        oldtracernames=.true.
101      endif
102
103c modifications possibles des variables de tab_cntrl
104      PRINT*
105      write(*,*) 'TABFI de phyeta0',Lmodif,tab0
106      call tabfi (nid,Lmodif,tab0,day_ini,lmax,p_rad,
107     .              p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
108c
109c Lecture des latitudes (coordonnees):
110c
111      ierr = NF_INQ_VARID (nid, "latitude", nvarid)
112      IF (ierr.NE.NF_NOERR) THEN
113         PRINT*, 'phyetat0: Le champ <latitude> est absent'
114         CALL abort
115      ENDIF
116#ifdef NC_DOUBLE
117      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, lati)
118#else
119      ierr = NF_GET_VAR_REAL(nid, nvarid, lati)
120#endif
121      IF (ierr.NE.NF_NOERR) THEN
122         PRINT*, 'phyetat0: Lecture echouee pour <latitude>'
123         CALL abort
124      ENDIF
125c
126c Lecture des longitudes (coordonnees):
127c
128      ierr = NF_INQ_VARID (nid, "longitude", nvarid)
129      IF (ierr.NE.NF_NOERR) THEN
130         PRINT*, 'phyetat0: Le champ <longitude> est absent'
131         CALL abort
132      ENDIF
133#ifdef NC_DOUBLE
134      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, long)
135#else
136      ierr = NF_GET_VAR_REAL(nid, nvarid, long)
137#endif
138      IF (ierr.NE.NF_NOERR) THEN
139         PRINT*, 'phyetat0: Lecture echouee pour <longitude>'
140         CALL abort
141      ENDIF
142c
143c Lecture des aires des mailles:
144c
145      ierr = NF_INQ_VARID (nid, "area", nvarid)
146      IF (ierr.NE.NF_NOERR) THEN
147         PRINT*, 'phyetat0: Le champ <area> est absent'
148         CALL abort
149      ENDIF
150#ifdef NC_DOUBLE
151      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, area)
152#else
153      ierr = NF_GET_VAR_REAL(nid, nvarid, area)
154#endif
155      IF (ierr.NE.NF_NOERR) THEN
156         PRINT*, 'phyetat0: Lecture echouee pour <area>'
157         CALL abort
158      ENDIF
159      xmin = 1.0E+20
160      xmax = -1.0E+20
161      xmin = MINVAL(area)
162      xmax = MAXVAL(area)
163      PRINT*,'Aires des mailles <area>:', xmin, xmax
164c
165c Lecture du geopotentiel au sol:
166c
167      ierr = NF_INQ_VARID (nid, "phisfi", nvarid)
168      IF (ierr.NE.NF_NOERR) THEN
169         PRINT*, 'phyetat0: Le champ <phisfi> est absent'
170         CALL abort
171      ENDIF
172#ifdef NC_DOUBLE
173      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, phisfi)
174#else
175      ierr = NF_GET_VAR_REAL(nid, nvarid, phisfi)
176#endif
177      IF (ierr.NE.NF_NOERR) THEN
178         PRINT*, 'phyetat0: Lecture echouee pour <phisfi>'
179         CALL abort
180      ENDIF
181      xmin = 1.0E+20
182      xmax = -1.0E+20
183      xmin = MINVAL(phisfi)
184      xmax = MAXVAL(phisfi)
185      PRINT*,'Geopotentiel au sol <phisfi>:', xmin, xmax
186c
187c Lecture de l''albedo du sol nu:
188c
189      ierr = NF_INQ_VARID (nid, "albedodat", nvarid)
190      IF (ierr.NE.NF_NOERR) THEN
191         PRINT*, 'phyetat0: Le champ <albedodat> est absent'
192         CALL abort
193      ENDIF
194#ifdef NC_DOUBLE
195      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, albedodat)
196#else
197      ierr = NF_GET_VAR_REAL(nid, nvarid, albedodat)
198#endif
199      IF (ierr.NE.NF_NOERR) THEN
200         PRINT*, 'phyetat0: Lecture echouee pour <albedodat>'
201         CALL abort
202      ENDIF
203      xmin = 1.0E+20
204      xmax = -1.0E+20
205      xmin = MINVAL(albedodat)
206      xmax = MAXVAL(albedodat)
207      PRINT*,'Albedo du sol nu <albedodat>:', xmin, xmax
208c
209c Lecture de l''inertie thermique du sol:
210c
211!      ierr = NF_INQ_VARID (nid, "inertiedat", nvarid)
212!      IF (ierr.NE.NF_NOERR) THEN
213!         PRINT*, 'phyetat0: Le champ <inertiedat> est absent'
214!         CALL abort
215!      ENDIF
216!#ifdef NC_DOUBLE
217!      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, inertiedat)
218!#else
219!      ierr = NF_GET_VAR_REAL(nid, nvarid, inertiedat)
220!#endif
221!      IF (ierr.NE.NF_NOERR) THEN
222!         PRINT*, 'phyetat0: Lecture echouee pour <inertiedat>'
223!         CALL abort
224!      ENDIF
225!      xmin = 1.0E+20
226!      xmax = -1.0E+20
227!      xmin = MINVAL(inertiedat)
228!      xmax = MAXVAL(inertiedat)
229!      PRINT*,'Inertie thermique du sol <inertiedat>:', xmin, xmax
230c
231c ZMEA
232c
233      ierr = NF_INQ_VARID (nid, "ZMEA", nvarid)
234      IF (ierr.NE.NF_NOERR) THEN
235         PRINT*, 'phyetat0: Le champ <ZMEA> est absent'
236         CALL abort
237      ENDIF
238#ifdef NC_DOUBLE
239      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, zmea)
240#else
241      ierr = NF_GET_VAR_REAL(nid, nvarid, zmea)
242#endif
243      IF (ierr.NE.NF_NOERR) THEN
244         PRINT*, 'phyetat0: Lecture echouee pour <ZMEA>'
245         CALL abort
246      ENDIF
247      xmin = 1.0E+20
248      xmax = -1.0E+20
249      DO i = 1, ngridmx
250         xmin = MIN(zmea(i),xmin)
251         xmax = MAX(zmea(i),xmax)
252      ENDDO
253      PRINT*,'<zmea>:', xmin, xmax
254c
255c ZSTD
256c
257      ierr = NF_INQ_VARID (nid, "ZSTD", nvarid)
258      IF (ierr.NE.NF_NOERR) THEN
259         PRINT*, 'phyetat0: Le champ <ZSTD> est absent'
260         CALL abort
261      ENDIF
262#ifdef NC_DOUBLE
263      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, zstd)
264#else
265      ierr = NF_GET_VAR_REAL(nid, nvarid, zstd)
266#endif
267      IF (ierr.NE.NF_NOERR) THEN
268         PRINT*, 'phyetat0: Lecture echouee pour <ZSTD>'
269         CALL abort
270      ENDIF
271      xmin = 1.0E+20
272      xmax = -1.0E+20
273      DO i = 1, ngridmx
274         xmin = MIN(zstd(i),xmin)
275         xmax = MAX(zstd(i),xmax)
276      ENDDO
277      PRINT*,'<zstd>:', xmin, xmax
278c
279c ZSIG
280c
281      ierr = NF_INQ_VARID (nid, "ZSIG", nvarid)
282      IF (ierr.NE.NF_NOERR) THEN
283         PRINT*, 'phyetat0: Le champ <ZSIG> est absent'
284         CALL abort
285      ENDIF
286#ifdef NC_DOUBLE
287      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, zsig)
288#else
289      ierr = NF_GET_VAR_REAL(nid, nvarid, zsig)
290#endif
291      IF (ierr.NE.NF_NOERR) THEN
292         PRINT*, 'phyetat0: Lecture echouee pour <ZSIG>'
293         CALL abort
294      ENDIF
295      xmin = 1.0E+20
296      xmax = -1.0E+20
297      DO i = 1, ngridmx
298         xmin = MIN(zsig(i),xmin)
299         xmax = MAX(zsig(i),xmax)
300      ENDDO
301      PRINT*,'<zsig>:', xmin, xmax
302c
303c ZGAM
304c
305      ierr = NF_INQ_VARID (nid, "ZGAM", nvarid)
306      IF (ierr.NE.NF_NOERR) THEN
307         PRINT*, 'phyetat0: Le champ <ZGAM> est absent'
308         CALL abort
309      ENDIF
310#ifdef NC_DOUBLE
311      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, zgam)
312#else
313      ierr = NF_GET_VAR_REAL(nid, nvarid, zgam)
314#endif
315      IF (ierr.NE.NF_NOERR) THEN
316         PRINT*, 'phyetat0: Lecture echouee pour <ZGAM>'
317         CALL abort
318      ENDIF
319      xmin = 1.0E+20
320      xmax = -1.0E+20
321      DO i = 1, ngridmx
322         xmin = MIN(zgam(i),xmin)
323         xmax = MAX(zgam(i),xmax)
324      ENDDO
325      PRINT*,'<zgam>:', xmin, xmax
326c
327c ZTHE
328c
329      ierr = NF_INQ_VARID (nid, "ZTHE", nvarid)
330      IF (ierr.NE.NF_NOERR) THEN
331         PRINT*, 'phyetat0: Le champ <ZTHE> est absent'
332         CALL abort
333      ENDIF
334#ifdef NC_DOUBLE
335      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, zthe)
336#else
337      ierr = NF_GET_VAR_REAL(nid, nvarid, zthe)
338#endif
339      IF (ierr.NE.NF_NOERR) THEN
340         PRINT*, 'phyetat0: Lecture echouee pour <ZTHE>'
341         CALL abort
342      ENDIF
343      xmin = 1.0E+20
344      xmax = -1.0E+20
345      DO i = 1, ngridmx
346         xmin = MIN(zthe(i),xmin)
347         xmax = MAX(zthe(i),xmax)
348      ENDDO
349      PRINT*,'<zthe>:', xmin, xmax
350c
351c CO2 ice cover
352c
353! Ehouarn: from now on, there is no "co2ice" standalone field; it is supposed
354! to be with stored in qsurf(i_co2_ice)
355!      ierr = NF_INQ_VARID (nid, "co2ice", nvarid)
356!      IF (ierr.NE.NF_NOERR) THEN
357!         PRINT*, 'phyetat0: Le champ <co2ice> est absent'
358!         CALL abort
359!      ENDIF
360!#ifdef NC_DOUBLE
361!      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, co2ice)
362!#else
363!      ierr = NF_GET_VAR_REAL(nid, nvarid, co2ice)
364!#endif
365!      IF (ierr.NE.NF_NOERR) THEN
366!         PRINT*, 'phyetat0: Lecture echouee pour <co2ice>'
367!         CALL abort
368!      ENDIF
369!      xmin = 1.0E+20
370!      xmax = -1.0E+20
371!      xmin = MINVAL(co2ice)
372!      xmax = MAXVAL(co2ice)
373!      PRINT*,'CO2 ice cover <co2ice>:', xmin, xmax
374c
375c Lecture des temperatures du sol:
376c
377      ierr = NF_INQ_VARID (nid, "tsurf", nvarid)
378      IF (ierr.NE.NF_NOERR) THEN
379         PRINT*, 'phyetat0: Le champ <tsurf> est absent'
380         PRINT*, '          Mais je vais essayer de lire TS**'
381         IF (nbsrf.GT.99) THEN
382            PRINT*, "Trop de sous-mailles"
383            CALL abort
384         ENDIF
385         DO nsrf = 1, nbsrf
386           WRITE(str2,'(i2.2)') nsrf
387           ierr = NF_INQ_VARID (nid, "TS"//str2, nvarid)
388           IF (ierr.NE.NF_NOERR) THEN
389           PRINT*, "phyetat0: Le champ <TS"//str2//"> est absent"
390              CALL abort
391           ENDIF
392#ifdef NC_DOUBLE
393           ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tsurf(1,nsrf))
394#else
395           ierr = NF_GET_VAR_REAL(nid, nvarid, tsurf(1,nsrf))
396#endif
397           IF (ierr.NE.NF_NOERR) THEN
398             PRINT*, "phyetat0: Lecture echouee pour <TS"//str2//">"
399             CALL abort
400           ENDIF
401           xmin = 1.0E+20
402           xmax = -1.0E+20
403           xmin = MINVAL(tsurf)
404           xmax = MAXVAL(tsurf)
405           PRINT*,'Temperature du sol TS**:', nsrf, xmin, xmax
406         ENDDO
407      ELSE
408         PRINT*, 'phyetat0: Le champ <tsurf> est present'
409         PRINT*, '          J ignore donc les autres temperatures TS**'
410#ifdef NC_DOUBLE
411         ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tsurf(1,1))
412#else
413         ierr = NF_GET_VAR_REAL(nid, nvarid, tsurf(1,1))
414#endif
415         IF (ierr.NE.NF_NOERR) THEN
416            PRINT*, "phyetat0: Lecture echouee pour <TSURF>"
417            CALL abort
418         ENDIF
419         xmin = 1.0E+20
420         xmax = -1.0E+20
421         xmin = MINVAL(tsurf)
422         xmax = MAXVAL(tsurf)
423         PRINT*,'Temperature du sol <tsurf>', xmin, xmax
424         IF (nbsrf >= 2) THEN
425            DO nsrf = 2, nbsrf
426               DO i = 1, ngridmx
427                  tsurf(i,nsrf) = tsurf(i,1)
428               ENDDO
429            ENDDO
430         ENDIF
431      ENDIF
432c
433c Lecture des temperatures du sol profond:
434c
435!      IF (nsoil.GT.99 .OR. nbsrf.GT.99) THEN
436!         PRINT*, "Trop de couches ou sous-mailles"
437!         CALL abort
438!      ENDIF
439!      DO nsrf = 1, nbsrf
440!         DO isoil=1, nsoil
441!            WRITE(str7,'(i2.2,"srf",i2.2)') isoil, nsrf
442!            ierr = NF_INQ_VARID (nid, 'tsoil', nvarid)
443!            IF (ierr.NE.NF_NOERR) THEN
444!               PRINT*, "phyetat0: Le champ <tsoil> est absent"
445!               PRINT*, "          Il prend donc la valeur de surface"
446!               DO i=1, ngridmx
447!                  tsoil(i,isoil,nsrf)=tsurf(i,nsrf)
448!               ENDDO
449!            ELSE
450!#ifdef NC_DOUBLE
451!              ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tsoil(1,1,nsrf))
452!#else
453!              ierr = NF_GET_VAR_REAL(nid, nvarid, tsoil(1,1,nsrf))
454!#endif
455!               IF (ierr.NE.NF_NOERR) THEN
456!                  PRINT*, "Lecture echouee pour <tsoil>"
457!                  CALL abort
458!               ENDIF
459!            ENDIF
460!         ENDDO
461!      ENDDO
462!      xmin = 1.0E+20
463!      xmax = -1.0E+20
464!      xmin = MINVAL(tsoil)
465!      xmax = MAXVAL(tsoil)
466!      PRINT*,'Temperatures du sol profond <tsoil>', xmin, xmax
467c
468c Surface emissivity
469c
470      ierr = NF_INQ_VARID (nid, "emis", nvarid)
471      IF (ierr.NE.NF_NOERR) THEN
472         PRINT*, 'phyetat0: Le champ <emis> est absent'
473         CALL abort
474      ENDIF
475#ifdef NC_DOUBLE
476      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, emis)
477#else
478      ierr = NF_GET_VAR_REAL(nid, nvarid, emis)
479#endif
480      IF (ierr.NE.NF_NOERR) THEN
481         PRINT*, 'phyetat0: Lecture echouee pour <emis>'
482         CALL abort
483      ENDIF
484      xmin = 1.0E+20
485      xmax = -1.0E+20
486      xmin = MINVAL(emis)
487      xmax = MAXVAL(emis)
488      PRINT*,'Surface emissivity <emis>:', xmin, xmax
489c
490c pbl wind variance
491c
492      ierr = NF_INQ_VARID (nid, "q2", nvarid)
493      IF (ierr.NE.NF_NOERR) THEN
494         PRINT*, 'phyetat0: Le champ <q2> est absent'
495         CALL abort
496      ENDIF
497#ifdef NC_DOUBLE
498      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, q2)
499#else
500      ierr = NF_GET_VAR_REAL(nid, nvarid, q2)
501#endif
502      IF (ierr.NE.NF_NOERR) THEN
503         PRINT*, 'phyetat0: Lecture echouee pour <q2>'
504         CALL abort
505      ENDIF
506      xmin = 1.0E+20
507      xmax = -1.0E+20
508      xmin = MINVAL(q2)
509      xmax = MAXVAL(q2)
510      PRINT*,'pbl wind variance <q2>:', xmin, xmax
511c
512c tracer on surface
513c
514
515      IF(nq.GE.1) THEN
516         nqold=nq
517         DO iq=1,nq
518!            str7(1:5)='qsurf'
519!            WRITE(str7(6:7),'(i2.2)') iq
520!            ierr = NF_INQ_VARID (nid,str7,nvarid)
521           IF (oldtracernames) THEN
522             txt=" "
523             write(txt,'(a5,i2.2)')'qsurf',iq
524           ELSE
525             txt=tnom(iq)
526             if (txt.eq."h2o_vap") then
527               ! There is no surface tracer for h2o_vap;
528               ! "h2o_ice" should be loaded instead
529               txt="h2o_ice"
530               write(*,*) 'phyetat0: loading surface tracer',
531     &                     ' h2o_ice instead of h2o_vap'
532             endif
533           ENDIF ! of IF (oldtracernames) THEN
534           ierr=NF_INQ_VARID(nid,txt,nvarid)
535           IF (ierr.NE.NF_NOERR) THEN
536             write(*,*) 'PHYETAT0: WARNING : surface tracer',trim(txt),
537     &                  ' not found in file'
538             write(*,*) trim(txt), ' set to 0'
539             do ig=1,ngridmx
540               qsurf(ig,iq)=0.
541             end do
542             nqold=min(iq-1,nqold)
543           ELSE
544#ifdef NC_DOUBLE
545             ierr = NF_GET_VAR_DOUBLE(nid, nvarid,qsurf(1,iq))
546#else
547             ierr = NF_GET_VAR_REAL(nid, nvarid,qsurf(1,iq))
548#endif
549             IF (ierr.NE.NF_NOERR) THEN
550               PRINT*, 'phyetat0: Lecture echouee pour <',trim(txt),'>'
551               CALL abort
552             ENDIF
553           ENDIF
554           xmin = 1.0E+20
555           xmax = -1.0E+20
556           xmin = MINVAL(qsurf(1:ngridmx,iq))
557           xmax = MAXVAL(qsurf(1:ngridmx,iq))
558           PRINT*,'tracer on surface <',trim(txt),'>:',xmin,xmax
559         ENDDO
560         if ((nqold.lt.nq).and.(nqold.ge.1)) then
561c        case when new tracer are added in addition to old ones
562             write(*,*)'qsurf 1 to ', nqold,'were already present'
563             write(*,*)'qsurf ', nqold+1,' to ', nqmx,'are new'
564             write(*,*)' and initialized to zero'
565             qsurf(:,nqold+1:nqmx)=0.0
566!            yes=' '
567!            do while ((yes.ne.'y').and.(yes.ne.'n'))
568!             write(*,*) 'Would you like to reindex qsurf # 1 ->',nqold
569!             write(*,*) 'to #',nqmx-nqold+1,'->', nqmx,'   (y or n) ?'
570!             read(*,fmt='(a)') yes
571!            end do
572!            if (yes.eq.'y') then
573!              write(*,*) 'OK, let s reindex qsurf'
574!                 do ig=1,ngridmx
575!                    do iq=nqmx,nqmx-nqold+1,-1
576!                       qsurf(ig,iq)=qsurf(ig,iq-nqmx+nqold)
577!                    end do
578!                    do iq=nqmx-nqold,1,-1
579!                       qsurf(ig,iq)= 0.
580!                    end do
581!                 end do
582!            end if
583         end if ! of if ((nqold.lt.nq).and.(nqold.ge.1))
584      ENDIF ! of IF(nq.GE.1)
585
586! Call to soil_settings, in order to read soil temperatures,
587! as well as thermal inertia and volumetric heat capacity
588
589      call soil_settings(nid,ngridmx,nsoil,tsurf,tsoil)
590c
591c Fermer le fichier:
592c
593      ierr = NF_CLOSE(nid)
594c
595      RETURN
596      END
Note: See TracBrowser for help on using the repository browser.