source: trunk/LMDZ.PLUTO/libf/phypluto/writediagfi.F @ 3627

Last change on this file since 3627 was 3623, checked in by afalco, 5 months ago

Pluto: better message for change in "temperature" written by diagfi (replaced the deprecated "temp" name).
zrecast also reads "temperature" if "temp" fails.
AF

File size: 21.4 KB
Line 
1      subroutine writediagfi(ngrid,nom,titre,unite,dim,px)
2
3!  Ecriture de variables diagnostiques au choix dans la physique
4!  dans un fichier NetCDF nomme  'diagfi'. Ces variables peuvent etre
5!  3d (ex : temperature), 2d (ex : temperature de surface), ou
6!  0d (pour un scalaire qui ne depend que du temps : ex : la longitude
7!  solaire)
8!  (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne)
9!  La periode d'ecriture est donnee par
10!  "ecritphy " regle dans le fichier de controle de run :  run.def
11!
12!    writediagfi peut etre appele de n'importe quelle subroutine
13!    de la physique, plusieurs fois. L'initialisation et la creation du
14!    fichier se fait au tout premier appel.
15!
16! WARNING : les variables dynamique (u,v,t,q,ps)
17!  sauvees par writediagfi avec une
18! date donnee sont legerement differentes que dans le fichier histoire car
19! on ne leur a pas encore ajoute de la dissipation et de la physique !!!
20! IL est  RECOMMANDE d'ajouter les tendance physique a ces variables
21! avant l'ecriture dans diagfi (cf. physiq.F)
22!
23! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4
24!         Oct 2011 Francois: enable having a 'diagfi.def' file to select
25!                            at runtime, which variables to put in file
26!
27!  parametres (input) :
28!  ----------
29!      ngrid : nombres de point ou est calcule la physique
30!                (ngrid = 2+(jjm-1)*iim - 1/jjm)
31!                 (= nlon ou klon dans la physique terrestre)
32!
33!      unit : unite logique du fichier de sortie (toujours la meme)
34!      nom  : nom de la variable a sortir (chaine de caracteres)
35!      titre: titre de la variable (chaine de caracteres)
36!      unite : unite de la variable (chaine de caracteres)
37!      px : variable a sortir (real 0, 1, 2, ou 3d)
38!      dim : dimension de px : 0, 1, 2, ou 3 dimensions
39!
40!=================================================================
41      use surfdat_h, only: phisfi
42      use geometry_mod, only: cell_area
43      use time_phylmdz_mod, only: ecritphy, day_step, iphysiq, day_ini
44      USE mod_phys_lmdz_para, only : is_parallel, is_mpi_root,
45     &                               is_master, gather
46      USE mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo,
47     &                              nbp_lon, nbp_lat, nbp_lev,
48     &                              grid_type, unstructured
49      implicit none
50
51! Commons
52      include "netcdf.inc"
53
54! Arguments on input:
55      integer,intent(in) :: ngrid
56      character (len=*),intent(in) :: nom,titre,unite
57      integer,intent(in) :: dim
58      real,intent(in) :: px(ngrid,nbp_lev)
59
60! Local variables:
61
62      real*4 dx3(nbp_lon+1,nbp_lat,nbp_lev) ! to store a 3D data set
63      real*4 dx2(nbp_lon+1,nbp_lat)     ! to store a 2D (surface) data set
64      real*4 dx1(nbp_lev)           ! to store a 1D (column) data set
65      real*4 dx0
66      real*4 dx3_1d(1,nbp_lev) ! to store a profile with 1D model
67      real*4 dx2_1d ! to store a surface value with 1D model
68
69      real*4,save :: date
70!$OMP THREADPRIVATE(date)
71
72      REAL phis((nbp_lon+1),nbp_lat)
73      REAL area((nbp_lon+1),nbp_lat)
74
75      integer irythme
76      integer ierr,ierr2
77      integer i,j,l, ig0
78
79      integer,save :: zitau=0
80      character(len=40),save :: firstnom='1234567890'
81!$OMP THREADPRIVATE(zitau,firstnom)
82
83! Ajouts
84      integer, save :: ntime=0
85!$OMP THREADPRIVATE(ntime)
86      integer :: idim,varid
87      integer :: nid
88      character(len=*),parameter :: fichnom="diagfi.nc"
89      integer, dimension(4) :: id
90      integer, dimension(4) :: edges,corner
91
92! Added to use diagfi.def to select output variable
93      logical,save :: diagfi_def
94      logical :: getout
95      integer,save :: n_nom_def
96      integer :: n
97      integer,parameter :: n_nom_def_max=199
98      character(len=120),save :: nom_def(n_nom_def_max)
99      logical,save :: firstcall=.true.
100!$OMP THREADPRIVATE(firstcall)  !diagfi_def,n_nom_def,nom_def read in diagfi.def
101
102#ifdef CPP_PARA
103! Added to work in parallel mode
104      real dx3_glop(klon_glo,nbp_lev)
105      real dx3_glo(nbp_lon,nbp_lat,nbp_lev) ! to store a global 3D data set
106      real dx2_glop(klon_glo)
107      real dx2_glo(nbp_lon,nbp_lat)     ! to store a global 2D (surface) data set
108      real px2(ngrid)
109!      real dx1_glo(nbp_lev)          ! to store a 1D (column) data set
110!      real dx0_glo
111      real phisfi_glo(klon_glo) ! surface geopotential on global physics grid
112      real areafi_glo(klon_glo) ! mesh area on global physics grid
113#else
114      real phisfi_glo(ngrid) ! surface geopotential on global physics grid
115      real areafi_glo(ngrid) ! mesh area on global physics grid
116#endif
117
118      if (grid_type==unstructured) then
119           return
120      endif
121
122!***************************************************************
123!Sortie des variables au rythme voulu
124
125      irythme = int(ecritphy) ! output rate set by ecritphy
126
127!***************************************************************
128
129! At very first call, check if there is a "diagfi.def" to use and read it
130! -----------------------------------------------------------------------
131      IF (firstcall) THEN
132         firstcall=.false.
133
134!$OMP MASTER
135  !      Open diagfi.def definition file if there is one:
136         open(99,file="diagfi.def",status='old',form='formatted',
137     s   iostat=ierr2)
138
139         if(ierr2.eq.0) then
140            diagfi_def=.true.
141            write(*,*) "******************"
142            write(*,*) "Reading diagfi.def"
143            write(*,*) "******************"
144            do n=1,n_nom_def_max
145              read(99,fmt='(a)',end=88) nom_def(n)
146              write(*,*) 'Output in diagfi: ', trim(nom_def(n))
147              if (nom_def(n).eq.'temp') then
148                 write(*,*) 'WARNING: "temp" is outdated.'
149                 write(*,*) 'Do you really have a "temp" variable?'
150                 write(*,*) 'Otherwise use "temperature" instead.'
151              end if
152            end do
153 88         continue
154            if (n.ge.n_nom_def_max) then
155               write(*,*)"n_nom_def_max too small in writediagfi.F:",n
156               call abort_physic("writediagfi",
157     &             "n_nom_def_max too small",1)
158            end if
159            n_nom_def=n-1
160            close(99)
161         else
162            diagfi_def=.false.
163         endif
164!$OMP END MASTER
165!$OMP BARRIER
166      END IF ! of IF (firstcall)
167
168! Get out of write_diagfi if there is diagfi.def AND variable not listed
169!  ---------------------------------------------------------------------
170      if (diagfi_def) then
171          getout=.true.
172          do n=1,n_nom_def
173             if(trim(nom_def(n)).eq.nom) getout=.false.
174          end do
175          if (getout) return
176      end if
177
178! Initialisation of 'firstnom' and create/open the "diagfi.nc" NetCDF file
179! ------------------------------------------------------------------------
180! (at very first call to the subroutine, in accordance with diagfi.def)
181
182      if (firstnom.eq.'1234567890') then ! .true. for the very first valid
183      !   call to this subroutine; now set 'firstnom'
184         firstnom = nom
185         ! just to be sure, check that firstnom is large enough to hold nom
186         if (len_trim(firstnom).lt.len_trim(nom)) then
187           write(*,*) "writediagfi: Error !!!"
188           write(*,*) "   firstnom string not long enough!!"
189           write(*,*) "   increase its size to at least ",len_trim(nom)
190           call abort_physic("writediagfi","firstnom too short",1)
191         endif
192
193         zitau = -1 ! initialize zitau
194
195#ifdef CPP_PARA
196          ! Gather phisfi() geopotential on physics grid
197          call Gather(phisfi,phisfi_glo)
198          ! Gather cell_area() mesh area on physics grid
199          call Gather(cell_area,areafi_glo)
200#else
201         phisfi_glo(:)=phisfi(:)
202         areafi_glo(:)=cell_area(:)
203#endif
204
205         !! parallel: we cannot use the usual writediagfi method
206!!         call iophys_ini
207         if (is_master) then
208         ! only the master is required to do this
209
210         ! Create the NetCDF file
211         ierr = NF_CREATE(fichnom, NF_CLOBBER, nid)
212         ! Define the 'Time' dimension
213         ierr = nf_def_dim(nid,"Time",NF_UNLIMITED,idim)
214         ! Define the 'Time' variable
215!#ifdef NC_DOUBLE
216!         ierr = NF_DEF_VAR (nid, "Time", NF_DOUBLE, 1, idim,varid)
217!#else
218         ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, idim,varid)
219!#endif
220         ! Add a long_name attribute
221         ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name",
222     .          4,"Time")
223         ! Add a units attribute
224         ierr = NF_PUT_ATT_TEXT(nid, varid,'units',29,
225     .          "days since 0000-00-0 00:00:00")
226         ! Switch out of NetCDF Define mode
227         ierr = NF_ENDDEF(nid)
228
229         ! Build phis() and area()
230         IF (klon_glo>1) THEN
231          do i=1,nbp_lon+1 ! poles
232           phis(i,1)=phisfi_glo(1)
233           phis(i,nbp_lat)=phisfi_glo(klon_glo)
234           ! for area, divide at the poles by nbp_lon
235           area(i,1)=areafi_glo(1)/nbp_lon
236           area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon
237          enddo
238          do j=2,nbp_lat-1
239           ig0= 1+(j-2)*nbp_lon
240           do i=1,nbp_lon
241              phis(i,j)=phisfi_glo(ig0+i)
242              area(i,j)=areafi_glo(ig0+i)
243           enddo
244           ! handle redundant point in longitude
245           phis(nbp_lon+1,j)=phis(1,j)
246           area(nbp_lon+1,j)=area(1,j)
247          enddo
248         ENDIF
249
250         ! write "header" of file (longitudes, latitudes, geopotential, ...)
251         IF (klon_glo>1) THEN ! general 3D case
252           call iniwrite(nid,day_ini,phis,area,nbp_lon+1,nbp_lat)
253         ELSE ! 1D model
254           call iniwrite(nid,day_ini,phisfi_glo(1),areafi_glo(1),1,1)
255         ENDIF
256
257         endif ! of if (is_master)
258
259      else
260
261         if (is_master) then
262           ! only the master is required to do this
263
264           ! Open the NetCDF file
265           ierr = NF_OPEN(fichnom,NF_WRITE,nid)
266         endif ! of if (is_master)
267
268      endif ! if (firstnom.eq.'1234567890')
269
270! Increment time index 'zitau' if it is the "fist call" (at given time level)
271! to writediagfi
272!------------------------------------------------------------------------
273      if (nom.eq.firstnom) then
274          zitau = zitau + iphysiq
275      end if
276
277!--------------------------------------------------------
278! Write the variables to output file if it's time to do so
279!--------------------------------------------------------
280
281      if ( MOD(zitau+1,irythme) .eq.0.) then
282
283! Compute/write/extend 'Time' coordinate (date given in days)
284! (done every "first call" (at given time level) to writediagfi)
285! Note: date is incremented as 1 step ahead of physics time
286!--------------------------------------------------------
287
288        if (is_master) then
289           ! only the master is required to do this
290        if (nom.eq.firstnom) then
291        ! We have identified a "first call" (at given date)
292           ntime=ntime+1 ! increment # of stored time steps
293           ! compute corresponding date (in days and fractions thereof)
294           date=(zitau +1.)/day_step
295           ! Get NetCDF ID of 'Time' variable
296           ierr= NF_INQ_VARID(nid,"Time",varid)
297           ! Write (append) the new date to the 'Time' array
298!#ifdef NC_DOUBLE
299!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date)
300!#else
301           ierr= NF_PUT_VARA_REAL(nid,varid,ntime,1,date)
302!#endif
303           if (ierr.ne.NF_NOERR) then
304              write(*,*) "***** PUT_VAR matter in writediagfi_nc"
305              write(*,*) "***** with time"
306              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
307c             call abort
308           endif
309
310           write(6,*)'WRITEDIAGFI: date= ', date
311        end if ! of if (nom.eq.firstnom)
312
313        endif ! of if (is_master)
314
315!Case of a 3D variable
316!---------------------
317        if (dim.eq.3) then
318
319#ifdef CPP_PARA
320          ! Gather field on a "global" (without redundant longitude) array
321          call Gather(px,dx3_glop)
322!$OMP MASTER
323          if (is_mpi_root) then
324            call Grid1Dto2D_glo(dx3_glop,dx3_glo)
325            ! copy dx3_glo() to dx3(:) and add redundant longitude
326            dx3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:)
327            dx3(nbp_lon+1,:,:)=dx3(1,:,:)
328          endif
329!$OMP END MASTER
330!$OMP BARRIER
331#else
332!         Passage variable physique -->  variable dynamique
333!         recast (copy) variable from physics grid to dynamics grid
334          IF (klon_glo>1) THEN ! General case
335           DO l=1,nbp_lev
336             DO i=1,nbp_lon+1
337                dx3(i,1,l)=px(1,l)
338                dx3(i,nbp_lat,l)=px(ngrid,l)
339             ENDDO
340             DO j=2,nbp_lat-1
341                ig0= 1+(j-2)*nbp_lon
342                DO i=1,nbp_lon
343                   dx3(i,j,l)=px(ig0+i,l)
344                ENDDO
345                dx3(nbp_lon+1,j,l)=dx3(1,j,l)
346             ENDDO
347           ENDDO
348          ELSE ! 1D model case
349           dx3_1d(1,1:nbp_lev)=px(1,1:nbp_lev)
350          ENDIF
351#endif
352!         Ecriture du champs
353
354          if (is_master) then
355           ! only the master writes to output
356! name of the variable
357           ierr= NF_INQ_VARID(nid,nom,varid)
358           if (ierr /= NF_NOERR) then
359! corresponding dimensions
360              ierr= NF_INQ_DIMID(nid,"longitude",id(1))
361              ierr= NF_INQ_DIMID(nid,"latitude",id(2))
362              ierr= NF_INQ_DIMID(nid,"altitude",id(3))
363              ierr= NF_INQ_DIMID(nid,"Time",id(4))
364
365! Create the variable if it doesn't exist yet
366
367              write (*,*) "=========================="
368              write (*,*) "DIAGFI: creating variable ",trim(nom)
369              call def_var(nid,nom,titre,unite,4,id,varid,ierr)
370
371           else
372             if (ntime==0) then
373              write(*,*) "DIAGFI Error: failed creating variable ",
374     &                   trim(nom)
375              write(*,*) "it seems it already exists!"
376              call abort_physic("writediagfi",
377     &             trim(nom)//" already exists",1)
378             endif
379           endif
380
381           corner(1)=1
382           corner(2)=1
383           corner(3)=1
384           corner(4)=ntime
385
386           IF (klon_glo==1) THEN
387             edges(1)=1
388           ELSE
389             edges(1)=nbp_lon+1
390           ENDIF
391           edges(2)=nbp_lat
392           edges(3)=nbp_lev
393           edges(4)=1
394!#ifdef NC_DOUBLE
395!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx3)
396!#else
397!           write(*,*)"test:  nid=",nid," varid=",varid
398!           write(*,*)"       corner()=",corner
399!           write(*,*)"       edges()=",edges
400!           write(*,*)"       dx3()=",dx3
401           IF (klon_glo>1) THEN ! General case
402             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3)
403           ELSE
404             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3_1d)
405           ENDIF
406!#endif
407
408           if (ierr.ne.NF_NOERR) then
409              write(*,*) "***** PUT_VAR problem in writediagfi"
410              write(*,*) "***** with dx3: ",trim(nom)
411              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
412              call abort_physic("writediagfi",
413     &             "failed writing "//trim(nom),1)
414           endif
415
416          endif !of if (is_master)
417
418!Case of a 2D variable
419!---------------------
420
421        else if (dim.eq.2) then
422
423#ifdef CPP_PARA
424          ! Gather field on a "global" (without redundant longitude) array
425          px2(:)=px(:,1)
426          call Gather(px2,dx2_glop)
427!$OMP MASTER
428          if (is_mpi_root) then
429            call Grid1Dto2D_glo(dx2_glop,dx2_glo)
430            ! copy dx2_glo() to dx2(:) and add redundant longitude
431            dx2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:)
432            dx2(nbp_lon+1,:)=dx2(1,:)
433          endif
434!$OMP END MASTER
435!$OMP BARRIER
436#else
437
438!         Passage variable physique -->  physique dynamique
439!         recast (copy) variable from physics grid to dynamics grid
440          IF (klon_glo>1) THEN ! General case
441             DO i=1,nbp_lon+1
442                dx2(i,1)=px(1,1)
443                dx2(i,nbp_lat)=px(ngrid,1)
444             ENDDO
445             DO j=2,nbp_lat-1
446                ig0= 1+(j-2)*nbp_lon
447                DO i=1,nbp_lon
448                   dx2(i,j)=px(ig0+i,1)
449                ENDDO
450                dx2(nbp_lon+1,j)=dx2(1,j)
451             ENDDO
452          ELSE ! 1D model case
453            dx2_1d=px(1,1)
454          ENDIF
455#endif
456
457          if (is_master) then
458           ! only the master writes to output
459!         write (*,*) 'In  writediagfi, on sauve:  ' , nom
460!         write (*,*) 'In  writediagfi. Estimated date = ' ,date
461           ierr= NF_INQ_VARID(nid,nom,varid)
462           if (ierr /= NF_NOERR) then
463! corresponding dimensions
464              ierr= NF_INQ_DIMID(nid,"longitude",id(1))
465              ierr= NF_INQ_DIMID(nid,"latitude",id(2))
466              ierr= NF_INQ_DIMID(nid,"Time",id(3))
467
468! Create the variable if it doesn't exist yet
469
470              write (*,*) "=========================="
471              write (*,*) "DIAGFI: creating variable ",trim(nom)
472
473              call def_var(nid,nom,titre,unite,3,id,varid,ierr)
474
475           else
476             if (ntime==0) then
477              write(*,*) "DIAGFI Error: failed creating variable ",
478     &                   trim(nom)
479              write(*,*) "it seems it already exists!"
480              call abort_physic("writediagfi",
481     &             trim(nom)//" already exists",1)
482             endif
483           endif
484
485           corner(1)=1
486           corner(2)=1
487           corner(3)=ntime
488           IF (klon_glo==1) THEN
489             edges(1)=1
490           ELSE
491             edges(1)=nbp_lon+1
492           ENDIF
493           edges(2)=nbp_lat
494           edges(3)=1
495
496
497!#ifdef NC_DOUBLE
498!           ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx2)
499!#else
500           IF (klon_glo>1) THEN ! General case
501             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2)
502           ELSE
503             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2_1d)
504           ENDIF
505!#endif
506
507           if (ierr.ne.NF_NOERR) then
508              write(*,*) "***** PUT_VAR matter in writediagfi"
509              write(*,*) "***** with dx2: ",trim(nom)
510              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
511              call abort_physic("writediagfi",
512     &             "failed writing "//trim(nom),1)
513           endif
514
515          endif !of if (is_master)
516
517!Case of a 1D variable (ie: a column)
518!---------------------------------------------------
519
520       else if (dim.eq.1) then
521         if (is_parallel) then
522           write(*,*) "writediagfi error: dim=1 not implemented ",
523     &                 "in parallel mode. Problem for ",trim(nom)
524              call abort_physic("writediagfi",
525     &             "failed writing "//trim(nom),1)
526         endif
527!         Passage variable physique -->  physique dynamique
528!         recast (copy) variable from physics grid to dynamics grid
529          do l=1,nbp_lev
530            dx1(l)=px(1,l)
531          enddo
532
533          ierr= NF_INQ_VARID(nid,nom,varid)
534           if (ierr /= NF_NOERR) then
535! corresponding dimensions
536              ierr= NF_INQ_DIMID(nid,"altitude",id(1))
537              ierr= NF_INQ_DIMID(nid,"Time",id(2))
538
539! Create the variable if it doesn't exist yet
540
541              write (*,*) "=========================="
542              write (*,*) "DIAGFI: creating variable ",trim(nom)
543
544              call def_var(nid,nom,titre,unite,2,id,varid,ierr)
545
546           else
547             if (ntime==0) then
548              write(*,*) "DIAGFI Error: failed creating variable ",
549     &                   trim(nom)
550              write(*,*) "it seems it already exists!"
551              call abort_physic("writediagfi",
552     &             trim(nom)//" already exists",1)
553             endif
554           endif
555
556           corner(1)=1
557           corner(2)=ntime
558
559           edges(1)=nbp_lev
560           edges(2)=1
561!#ifdef NC_DOUBLE
562!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx1)
563!#else
564           ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx1)
565!#endif
566
567           if (ierr.ne.NF_NOERR) then
568              write(*,*) "***** PUT_VAR problem in writediagfi"
569              write(*,*) "***** with dx1: ",trim(nom)
570              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
571              call abort_physic("writediagfi",
572     &             "failed writing "//trim(nom),1)
573           endif
574
575!Case of a 0D variable (ie: a time-dependent scalar)
576!---------------------------------------------------
577
578        else if (dim.eq.0) then
579
580           dx0 = px (1,1)
581
582          if (is_master) then
583           ! only the master writes to output
584           ierr= NF_INQ_VARID(nid,nom,varid)
585           if (ierr /= NF_NOERR) then
586! corresponding dimensions
587              ierr= NF_INQ_DIMID(nid,"Time",id(1))
588
589! Create the variable if it doesn't exist yet
590
591              write (*,*) "=========================="
592              write (*,*) "DIAGFI: creating variable ",trim(nom)
593
594              call def_var(nid,nom,titre,unite,1,id,varid,ierr)
595
596           else
597             if (ntime==0) then
598              write(*,*) "DIAGFI Error: failed creating variable ",
599     &                   trim(nom)
600              write(*,*) "it seems it already exists!"
601              call abort_physic("writediagfi",
602     &             trim(nom)//" already exists",1)
603             endif
604           endif
605
606           corner(1)=ntime
607           edges(1)=1
608
609!#ifdef NC_DOUBLE
610!           ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx0)
611!#else
612           ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx0)
613!#endif
614           if (ierr.ne.NF_NOERR) then
615              write(*,*) "***** PUT_VAR matter in writediagfi"
616              write(*,*) "***** with dx0: ",trim(nom)
617              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
618              call abort_physic("writediagfi",
619     &             "failed writing "//trim(nom),1)
620           endif
621
622          endif !of if (is_master)
623
624        endif ! of if (dim.eq.3) elseif(dim.eq.2)...
625
626      endif ! of if ( MOD(zitau+1,irythme) .eq.0.)
627
628      if (is_master) then
629        ierr= NF_CLOSE(nid)
630      endif
631
632      end
Note: See TracBrowser for help on using the repository browser.