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

Last change on this file since 1242 was 1130, checked in by emillour, 11 years ago

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

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