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

Last change on this file since 374 was 305, checked in by rwordsworth, 13 years ago

Several new files added as part of the climate evolution model
(main program kcm.F90). Some general cleanup in physiq.F90 and
callcorrk.F90. Bugs in dust radiative transfer and H2 Rayleigh
scattering corrected.

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