source: trunk/LMDZ.MARS/libf/phymars/writediagfi.F

Last change on this file was 3369, checked in by emillour, 6 months ago

Mars PCM:
Change the way the rate of outputs for diagfi.nc files is specified:
IMPORTANT: Specifying "ecritphy" no longer possible and will trigger an error.
Use "outputs_per_sol" to specify output rate instead.
This should makes things (hopefully) clearer for users and also better
enforces a cleaner and clearer separation between dynamics and physics.
EM

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