source: LMDZ5/branches/LMDZ5_SPLA/libf/dyn3dpar/initfluxsto_p.F @ 5360

Last change on this file since 5360 was 1907, checked in by lguez, 11 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
Name of program: LMDZ
Creation date: 1984
Version: LMDZ5
License: CeCILL version 2
Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
See the license file in the root directory

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 7.6 KB
Line 
1!
2! $Id: initfluxsto_p.F 1907 2013-11-26 13:10:46Z abarral $
3!
4      subroutine initfluxsto_p
5     .  (infile,tstep,t_ops,t_wrt,
6     .                    fileid,filevid,filedid)
7
8#ifdef CPP_IOIPSL
9! This routine needs IOIPSL
10       USE IOIPSL
11#endif
12       USE parallel_lmdz
13       use Write_field
14       use misc_mod
15       
16      implicit none
17
18C
19C   Routine d'initialisation des ecritures des fichiers histoires LMDZ
20C   au format IOIPSL
21C
22C   Appels succesifs des routines: histbeg
23C                                  histhori
24C                                  histver
25C                                  histdef
26C                                  histend
27C
28C   Entree:
29C
30C      infile: nom du fichier histoire a creer
31C      day0,anne0: date de reference
32C      tstep: duree du pas de temps en seconde
33C      t_ops: frequence de l'operation pour IOIPSL
34C      t_wrt: frequence d'ecriture sur le fichier
35C
36C   Sortie:
37C      fileid: ID du fichier netcdf cree
38C      filevid:ID du fichier netcdf pour la grille v
39C
40C   L. Fairhead, LMD, 03/99
41C
42C =====================================================================
43C
44C   Declarations
45#include "dimensions.h"
46#include "paramet.h"
47#include "comconst.h"
48#include "comvert.h"
49#include "comgeom.h"
50#include "temps.h"
51#include "ener.h"
52#include "logic.h"
53#include "description.h"
54#include "serre.h"
55#include "iniprint.h"
56
57C   Arguments
58C
59      character*(*) infile
60      real tstep, t_ops, t_wrt
61      integer fileid, filevid,filedid
62
63#ifdef CPP_IOIPSL
64! This routine needs IOIPSL
65C   Variables locales
66C
67      real nivd(1)
68      integer tau0
69      real zjulian
70      character*3 str
71      character*10 ctrac
72      integer iq
73      real rlong(iip1,jjp1), rlat(iip1,jjp1),rl(1,1)
74      integer uhoriid, vhoriid, thoriid, zvertiid,dhoriid,dvertiid
75      integer ii,jj
76      integer zan, idayref
77      logical ok_sync
78      integer :: jjb,jje,jjn
79
80! definition du domaine d'ecriture pour le rebuild
81
82      INTEGER,DIMENSION(2) :: ddid
83      INTEGER,DIMENSION(2) :: dsg
84      INTEGER,DIMENSION(2) :: dsl
85      INTEGER,DIMENSION(2) :: dpf
86      INTEGER,DIMENSION(2) :: dpl
87      INTEGER,DIMENSION(2) :: dhs
88      INTEGER,DIMENSION(2) :: dhe
89     
90      INTEGER :: dynu_domain_id
91      INTEGER :: dynv_domain_id
92
93C
94C  Initialisations
95C
96      pi = 4. * atan (1.)
97      str='q  '
98      ctrac = 'traceur   '
99      ok_sync = .true.
100C
101C  Appel a histbeg: creation du fichier netcdf et initialisations diverses
102C         
103
104      zan = annee_ref
105      idayref = day_ref
106      CALL ymds2ju(zan, 1, idayref, 0.0, zjulian)
107      tau0 = itau_dyn
108       
109        do jj = 1, jjp1
110        do ii = 1, iip1
111          rlong(ii,jj) = rlonu(ii) * 180. / pi
112          rlat(ii,jj) = rlatu(jj) * 180. / pi
113        enddo
114      enddo
115
116      jjb=jj_begin
117      jje=jj_end
118      jjn=jj_nb
119
120      ddid=(/ 1,2 /)
121      dsg=(/ iip1,jjp1 /)
122      dsl=(/ iip1,jjn /)
123      dpf=(/ 1,jjb /)
124      dpl=(/ iip1,jje /)
125      dhs=(/ 0,0 /)
126      dhe=(/ 0,0 /)
127
128      call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe,
129     .                 'box',dynu_domain_id)
130       
131      call histbeg(trim(infile),iip1, rlong(:,1), jjn, rlat(1,jjb:jje),
132     .             1, iip1, 1, jjn, tau0, zjulian, tstep, uhoriid,
133     .             fileid,dynu_domain_id)
134C
135C  Creation du fichier histoire pour la grille en V (oblige pour l'instant,
136C  IOIPSL ne permet pas de grilles avec des nombres de point differents dans
137C  un meme fichier)
138
139
140      do jj = 1, jjm
141        do ii = 1, iip1
142          rlong(ii,jj) = rlonv(ii) * 180. / pi
143          rlat(ii,jj) = rlatv(jj) * 180. / pi
144        enddo
145      enddo
146
147      jjb=jj_begin
148      jje=jj_end
149      jjn=jj_nb
150      if (pole_sud) jje=jj_end-1
151      if (pole_sud) jjn=jj_nb-1
152
153      ddid=(/ 1,2 /)
154      dsg=(/ iip1,jjm /)
155      dsl=(/ iip1,jjn /)
156      dpf=(/ 1,jjb /)
157      dpl=(/ iip1,jje /)
158      dhs=(/ 0,0 /)
159      dhe=(/ 0,0 /)
160
161      call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe,
162     .                 'box',dynv_domain_id)
163     
164      call histbeg('fluxstokev',iip1, rlong(:,1), jjn, rlat(1,jjb:jje),
165     .             1, iip1, 1, jjn,tau0, zjulian, tstep, vhoriid,
166     .             filevid,dynv_domain_id)
167       
168      rl(1,1) = 1.     
169     
170      if (mpi_rank==0) then
171         
172        call histbeg('defstoke.nc', 1, rl, 1, rl,
173     .               1, 1, 1, 1,
174     .               tau0, zjulian, tstep, dhoriid, filedid)
175     
176      endif
177C
178C  Appel a histhori pour rajouter les autres grilles horizontales
179C
180      do jj = 1, jjp1
181        do ii = 1, iip1
182          rlong(ii,jj) = rlonv(ii) * 180. / pi
183          rlat(ii,jj) = rlatu(jj) * 180. / pi
184        enddo
185      enddo
186
187      jjb=jj_begin
188      jje=jj_end
189      jjn=jj_nb
190
191      call histhori(fileid, iip1, rlong(:,jjb:jje),jjn,rlat(:,jjb:jje),
192     .             'scalar','Grille points scalaires', thoriid)
193       
194C
195C  Appel a histvert pour la grille verticale
196C
197      call histvert(fileid, 'sig_s', 'Niveaux sigma',
198     . 'sigma_level',
199     .              llm, nivsigs, zvertiid)
200C Pour le fichier V
201      call histvert(filevid, 'sig_s', 'Niveaux sigma',
202     .  'sigma_level',
203     .              llm, nivsigs, zvertiid)
204c pour le fichier def
205      if (mpi_rank==0) then
206         nivd(1) = 1
207         call histvert(filedid, 'sig_s', 'Niveaux sigma',
208     .        'sigma_level',
209     .        1, nivd, dvertiid)
210      endif
211C
212C  Appels a histdef pour la definition des variables a sauvegarder
213       
214        CALL histdef(fileid, "phis", "Surface geop. height", "-",
215     .                iip1,jjn,thoriid, 1,1,1, -99, 32,
216     .                "once", t_ops, t_wrt)
217
218         CALL histdef(fileid, "aire", "Grid area", "-",
219     .                iip1,jjn,thoriid, 1,1,1, -99, 32,
220     .                "once", t_ops, t_wrt)
221       
222        if (mpi_rank==0) then
223       
224        CALL histdef(filedid, "dtvr", "tps dyn", "s",
225     .                1,1,dhoriid, 1,1,1, -99, 32,
226     .                "once", t_ops, t_wrt)
227       
228         CALL histdef(filedid, "istdyn", "tps stock", "s",
229     .                1,1,dhoriid, 1,1,1, -99, 32,
230     .                "once", t_ops, t_wrt)
231         
232         CALL histdef(filedid, "istphy", "tps stock phy", "s",
233     .                1,1,dhoriid, 1,1,1, -99, 32,
234     .                "once", t_ops, t_wrt)
235
236        endif
237C
238C Masse
239C
240      call histdef(fileid, 'masse', 'Masse', 'kg',
241     .             iip1, jjn, thoriid, llm, 1, llm, zvertiid,
242     .             32, 'inst(X)', t_ops, t_wrt)
243C
244C  Pbaru
245C
246      call histdef(fileid, 'pbaru', 'flx de masse zonal', 'kg m/s',
247     .             iip1, jjn, uhoriid, llm, 1, llm, zvertiid,
248     .             32, 'inst(X)', t_ops, t_wrt)
249
250C
251C  Pbarv
252C
253      if (pole_sud) jjn=jj_nb-1
254     
255      call histdef(filevid, 'pbarv', 'flx de masse mer', 'kg m/s',
256     .             iip1, jjn, vhoriid, llm, 1, llm, zvertiid,
257     .             32, 'inst(X)', t_ops, t_wrt)
258C
259C  w
260C
261      if (pole_sud) jjn=jj_nb
262      call histdef(fileid, 'w', 'flx de masse vert', 'kg m/s',
263     .             iip1, jjn, thoriid, llm, 1, llm, zvertiid,
264     .             32, 'inst(X)', t_ops, t_wrt)
265
266C
267C  Temperature potentielle
268C
269      call histdef(fileid, 'teta', 'temperature potentielle', '-',
270     .             iip1, jjn, thoriid, llm, 1, llm, zvertiid,
271     .             32, 'inst(X)', t_ops, t_wrt)
272C
273
274C
275C Geopotentiel
276C
277      call histdef(fileid, 'phi', 'geopotentiel instantane', '-',
278     .             iip1, jjn, thoriid, llm, 1, llm, zvertiid,
279     .             32, 'inst(X)', t_ops, t_wrt)
280C
281C  Fin
282C
283      call histend(fileid)
284      call histend(filevid)
285      if (mpi_rank==0) call histend(filedid)
286      if (ok_sync) then
287        call histsync(fileid)
288        call histsync(filevid)
289        if (mpi_rank==0) call histsync(filedid)
290      endif
291       
292#else
293      write(lunout,*)'initfluxsto_p: Needs IOIPSL to function'
294#endif
295! #endif of #ifdef CPP_IOIPSL
296      return
297      end
Note: See TracBrowser for help on using the repository browser.