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

Last change on this file since 3026 was 2616, checked in by romain.vande, 3 years ago

LMDZ_MARS RV : Open_MP;
Put all the "save" variables as "!$OMP THREADPRIVATE" in phymars.
The code can now be tested, see README for more info

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