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

Last change on this file since 3026 was 2900, checked in by romain.vande, 23 months ago

Mars PCM:
Following r2896, further implementation of subslope parametrisation.
Carefull ! This is a devolpment revision and it still need improvements and tests.
However, this commit should not change anything for nslope=1.
Only nslope=1 is possible for now!
Utilitaries are not adapted yet.
Dimension of some variables (30 listed below) is changed to add the nslope dimension.
Outputs (diagfi and restartfi) are not changed yet.
nslope is now read and written in the startfi.nc

Changed variables:
_ In surfdat_h:

  • qsurf
  • tsurf
  • watercap
  • emis
  • capcal
  • fluxgrd

_ In comsoil_h

  • tsoil
  • mthermdiff
  • thermdiff
  • coefd
  • alph
  • beta

_ In dimradmars_mod

  • albedo
  • fluxrad_sky
  • fluxrad

_ In physiq_mod

  • inertiesoil
  • fluxsurf_lw
  • fluxsurf_dn_sw
  • dqsurf
  • zdtsurf
  • zdtsdif
  • zdtsurfc
  • zdqsdif
  • zdqsc
  • dwatercap
  • dwatercap_dif
  • zflubid
  • fluxsurf_dn_sw_tot
  • inertiedat_tmp

New variables called VAR_meshavg correspond to the mesh average over all the subslope distribution of the variable

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