source: trunk/LMDZ.UNIVERSAL/libf/phygeneric/iotd_ecrit.F90 @ 965

Last change on this file since 965 was 862, checked in by aslmd, 12 years ago

LMDZ.UNIVERSAL
LMDZ.GENERIC

Added calls to Frederic Hourdin's subroutines to output physical fields in parallel. This is under precompiling flags CPP_PARA.
This allows for LMDZ.UNIVERSAL users to use writediagfi with parallel computations.
These lines are not compiled by casual users of LMDZ.GENERIC (or users of LMDZ.UNIVERSAL in sequential mode).

The precompiling flag NOWRITEDIAGFI is obsolete and has been deleted.

NOTES:

  • A better cleaner method might be proposed later
  • Added subroutines

A 0 LMDZ.UNIVERSAL/libf/phygeneric/iotd_ecrit.F90
A 0 LMDZ.UNIVERSAL/libf/phygeneric/iophys.F90
A 0 LMDZ.UNIVERSAL/libf/phygeneric/iotd.h
A 0 LMDZ.UNIVERSAL/libf/phygeneric/iotd_ini.F90
A 0 LMDZ.UNIVERSAL/libf/phygeneric/iotd_fin.F90
The iotd_* subroutines are actually supposed to be put in bibio in LMDZ.COMMON
This will be done later once agreed in the team

File size: 5.6 KB
Line 
1      subroutine iotd_ecrit(nom,llm,titre,unite,px)
2
3!!
4!! Provided by Frederic Hourdin 01/2013
5!!
6
7
8
9
10!  Ecriture de variables diagnostiques au choix dans la physique
11!  dans un fichier NetCDF nomme  'diagfi'. Ces variables peuvent etre
12!  3d (ex : temperature), 2d (ex : temperature de surface), ou
13!  0d (pour un scalaire qui ne depend que du temps : ex : la longitude
14!  solaire)
15!  (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne)
16!  La periode d'ecriture est donnee par
17!  "ecritphy " regle dans le fichier de controle de run :  run.def
18!
19!    writediagfi peut etre appele de n'importe quelle subroutine
20!    de la physique, plusieurs fois. L'initialisation et la creation du
21!    fichier se fait au tout premier appel.
22!
23! WARNING : les variables dynamique (u,v,t,q,ps)
24!  sauvees par writediagfi avec une
25! date donnee sont legerement differentes que dans le fichier histoire car
26! on ne leur a pas encore ajoute de la dissipation et de la physique !!!
27! IL est  RECOMMANDE d'ajouter les tendance physique a ces variables
28! avant l'ecriture dans diagfi (cf. physiq.F)
29
30! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4
31!
32!  parametres (input) :
33!  ----------
34!      unit : unite logique du fichier de sortie (toujours la meme)
35!      nom  : nom de la variable a sortir (chaine de caracteres)
36!      titre: titre de la variable (chaine de caracteres)
37!      unite : unite de la variable (chaine de caracteres)
38!      px : variable a sortir (real 0, 1, 2, ou 3d)
39!
40!=================================================================
41 
42      implicit none
43
44! Commons
45
46#include "netcdf.inc"
47#include "iotd.h"
48
49
50! Arguments on input:
51      integer llm
52      character (len=*) :: nom,titre,unite
53      integer imjmax
54      parameter (imjmax=100000)
55      real px(imjmax*llm)
56
57! Local variables:
58
59      real*4 date
60      real*4 zx(imjmax*llm)
61
62
63      integer ierr,ndim,dim_cc(4)
64      integer iq
65      integer i,j,l
66
67      integer zitau
68      character firstnom*20
69      SAVE firstnom
70      SAVE zitau
71      SAVE date
72      data firstnom /'1234567890'/
73      data zitau /0/
74
75! Ajouts
76      integer, save :: ntime=0
77      integer :: idim,varid
78      character (len =50):: fichnom
79      integer, dimension(4) :: id
80      integer, dimension(4) :: edges,corner
81     
82
83!***************************************************************
84! Initialisation of 'firstnom' and create/open the "diagfi.nc" NetCDF file
85! ------------------------------------------------------------------------
86! (Au tout premier appel de la subroutine durant le run.)
87
88
89!--------------------------------------------------------
90! Write the variables to output file if it's time to do so
91!--------------------------------------------------------
92
93
94! Compute/write/extend 'Time' coordinate (date given in days)
95! (done every "first call" (at given time level) to writediagfi)
96! Note: date is incremented as 1 step ahead of physics time
97!--------------------------------------------------------
98
99        zx(1:imax*jmax*llm)=px(1:imax*jmax*llm)
100        if (firstnom =='1234567890') then
101            firstnom=nom
102        endif
103
104!      print*,'nom ',nom,firstnom
105
106!! Quand on tombe sur la premiere variable on ajoute un pas de temps
107        if (nom.eq.firstnom) then
108        ! We have identified a "first call" (at given date)
109
110           ntime=ntime+1 ! increment # of stored time steps
111
112!!          print*,'ntime ',ntime
113           date=ntime
114!          date= float (zitau +1)/float (day_step)
115
116           ! compute corresponding date (in days and fractions thereof)
117           ! Get NetCDF ID of 'Time' variable
118
119           ierr=NF_SYNC(nid)
120
121           ierr= NF_INQ_VARID(nid,"Time",varid)
122           ! Write (append) the new date to the 'Time' array
123
124
125           ierr= NF_PUT_VARA_REAL(nid,varid,ntime,1,date)
126
127        !  print*,'date ',date,ierr,nid
128        !print*,'IOTD Date ,varid,nid,ntime,date',varid,nid,ntime,date
129
130           if (ierr.ne.NF_NOERR) then
131              write(*,*) "***** PUT_VAR matter in writediagfi_nc"
132              write(*,*) "***** with time"
133              write(*,*) 'ierr=', ierr   
134              stop
135           endif
136
137!          write(6,*)'WRITEDIAGFI: date= ', date
138        end if ! of if (nom.eq.firstnom)
139
140
141!Case of a 3D variable
142!---------------------
143        if (llm==lmax) then
144           ndim=4
145           corner(1)=1
146           corner(2)=1
147           corner(3)=1
148           corner(4)=ntime
149           edges(1)=imax
150           edges(2)=jmax
151           edges(3)=llm
152           edges(4)=1
153           dim_cc=dim_coord
154
155
156!Case of a 2D variable
157!---------------------
158
159        else if (llm==1) then
160           ndim=3
161           corner(1)=1
162           corner(2)=1
163           corner(3)=ntime
164           corner(4)=1
165           edges(1)=imax
166           edges(2)=jmax
167           edges(3)=1
168           edges(4)=1
169           dim_cc(1:2)=dim_coord(1:2)
170           dim_cc(3)=dim_coord(4)
171
172        endif ! of if llm=1 ou llm
173
174! AU premier pas de temps, on crée les variables
175!-----------------------------------------------
176
177      if (ntime==1) then
178          ierr = NF_REDEF (nid)
179          ierr = NF_DEF_VAR(nid,nom,NF_FLOAT,ndim,dim_cc,varid)
180          print*,'DEF ',nom,nid,varid
181          ierr = NF_ENDDEF(nid)
182      else
183         ierr= NF_INQ_VARID(nid,nom,varid)
184          print*,'INQ ',nom,nid,varid
185! Commandes pour recuperer automatiquement les coordonnees
186!             ierr= NF_INQ_DIMID(nid,"longitude",id(1))
187      endif
188
189
190      ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,zx)
191
192      if (ierr.ne.NF_NOERR) then
193           write(*,*) "***** PUT_VAR problem in writediagfi"
194           write(*,*) "***** with ",nom
195           write(*,*) 'ierr=', ierr
196      endif
197
198
199      end
Note: See TracBrowser for help on using the repository browser.