source: trunk/LMDZ.MARS/libf/phymars/writediagmicrofi.F @ 3493

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