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

Last change on this file since 4029 was 3969, checked in by jmauxion, 6 months ago

Mars PCM:
Bug fix: threadsave statement missing on variable save "isample". Could fail
when writting diagfi.
JM

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