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

Last change on this file since 4027 was 4027, checked in by debatzbr, 42 hours ago

Pluto PCM: Adding a slow_diagfi flag for 1D model from Generic PCM (revision 3928).
BBT

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