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

Last change on this file since 3813 was 3749, checked in by afalco, 8 months ago

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