source: trunk/libf/bibio/initdynav.F @ 6

Last change on this file since 6 was 1, checked in by emillour, 14 years ago

Import initial LMDZ5

File size: 5.8 KB
Line 
1!
2! $Id: initdynav.F 1403 2010-07-01 09:02:53Z fairhead $
3!
4      subroutine initdynav(day0,anne0,tstep,t_ops,t_wrt)
5
6#ifdef CPP_IOIPSL
7       USE IOIPSL
8#endif
9       USE infotrac, ONLY : nqtot, ttext
10      use com_io_dyn_mod, only : histaveid,histvaveid,histuaveid,       &
11     &        dynhistave_file,dynhistvave_file,dynhistuave_file
12      implicit none
13
14C
15C   Routine d'initialisation des ecritures des fichiers histoires LMDZ
16C   au format IOIPSL. Initialisation du fichier histoire moyenne.
17C
18C   Appels succesifs des routines: histbeg
19C                                  histhori
20C                                  histver
21C                                  histdef
22C                                  histend
23C
24C   Entree:
25C
26C      infile: nom du fichier histoire a creer
27C      day0,anne0: date de reference
28C      tstep : frequence d'ecriture
29C      t_ops: frequence de l'operation pour IOIPSL
30C      t_wrt: frequence d'ecriture sur le fichier
31C
32C
33C   L. Fairhead, LMD, 03/99
34C
35C =====================================================================
36C
37C   Declarations
38#include "dimensions.h"
39#include "paramet.h"
40#include "comconst.h"
41#include "comvert.h"
42#include "comgeom.h"
43#include "temps.h"
44#include "ener.h"
45#include "logic.h"
46#include "description.h"
47#include "serre.h"
48#include "iniprint.h"
49
50C   Arguments
51C
52      integer day0, anne0
53      real tstep, t_ops, t_wrt
54
55#ifdef CPP_IOIPSL
56! This routine needs IOIPSL to work
57C   Variables locales
58C
59      integer tau0
60      real zjulian
61      integer iq
62      real rlong(iip1,jjp1), rlat(iip1,jjp1)
63      integer uhoriid, vhoriid, thoriid, zvertiid
64      integer ii,jj
65      integer zan, dayref
66C
67C  Initialisations
68C
69      pi = 4. * atan (1.)
70C
71C  Appel a histbeg: creation du fichier netcdf et initialisations diverses
72C         
73
74      zan = anne0
75      dayref = day0
76      CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
77      tau0 = itau_dyn
78     
79      do jj = 1, jjp1
80        do ii = 1, iip1
81          rlong(ii,jj) = rlonv(ii) * 180. / pi
82          rlat(ii,jj)  = rlatu(jj) * 180. / pi
83        enddo
84      enddo
85       
86! Creation de 3 fichiers pour les differentes grilles horizontales
87! Restriction de IOIPSL: seulement 2 coordonnees dans le meme fichier
88! Grille Scalaire       
89      call histbeg(dynhistave_file, iip1, rlong(:,1), jjp1, rlat(1,:),
90     .             1, iip1, 1, jjp1,
91     .             tau0, zjulian, tstep, thoriid,histaveid)
92
93C  Creation du fichier histoire pour les grilles en V et U (oblige pour l'instant,
94C  IOIPSL ne permet pas de grilles avec des nombres de point differents dans
95C  un meme fichier)
96! Grille V
97      do jj = 1, jjm
98        do ii = 1, iip1
99          rlong(ii,jj) = rlonv(ii) * 180. / pi
100          rlat(ii,jj) = rlatv(jj) * 180. / pi
101        enddo
102      enddo
103
104      call histbeg(dynhistvave_file, iip1, rlong(:,1), jjm, rlat(1,:),
105     .             1, iip1, 1, jjm,
106     .             tau0, zjulian, tstep, vhoriid,histvaveid)
107! Grille U
108      do jj = 1, jjp1
109        do ii = 1, iip1
110          rlong(ii,jj) = rlonu(ii) * 180. / pi
111          rlat(ii,jj) = rlatu(jj) * 180. / pi
112        enddo
113      enddo
114
115      call histbeg(dynhistuave_file, iip1, rlong(:,1),jjp1, rlat(1,:),
116     .             1, iip1, 1, jjp1,
117     .             tau0, zjulian, tstep, uhoriid,histuaveid)
118C
119C  Appel a histvert pour la grille verticale
120C
121      call histvert(histaveid,'presnivs','Niveaux Pression
122     &     approximatifs','mb',llm, presnivs/100., zvertiid,'down')
123      call histvert(histuaveid,'presnivs','Niveaux Pression
124     &     approximatifs','mb',llm, presnivs/100., zvertiid,'down')
125      call histvert(histvaveid,'presnivs','Niveaux Pression
126     &     approximatifs','mb',llm, presnivs/100., zvertiid,'down')
127C
128C  Appels a histdef pour la definition des variables a sauvegarder
129C
130C  Vents U
131C
132!      write(6,*)'inithistave',tstep
133      call histdef(histuaveid, 'u', 'vent u moyen ',
134     .             'm/s', iip1, jjp1, uhoriid, llm, 1, llm, zvertiid,
135     .             32, 'ave(X)', t_ops, t_wrt)
136
137C  Vents V
138C
139      call histdef(histvaveid, 'v', 'vent v moyen',
140     .             'm/s', iip1, jjm, vhoriid, llm, 1, llm, zvertiid,
141     .             32, 'ave(X)', t_ops, t_wrt)
142
143C
144C  Temperature
145C
146      call histdef(histaveid, 'temp', 'temperature moyenne', 'K',
147     .             iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
148     .             32, 'ave(X)', t_ops, t_wrt)
149C
150C  Temperature potentielle
151C
152      call histdef(histaveid, 'theta', 'temperature potentielle', 'K',
153     .             iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
154     .             32, 'ave(X)', t_ops, t_wrt)
155C
156C  Geopotentiel
157C
158      call histdef(histaveid, 'phi', 'geopotentiel moyen', '-',
159     .             iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
160     .             32, 'ave(X)', t_ops, t_wrt)
161C
162C  Traceurs
163C
164!        DO iq=1,nqtot
165!          call histdef(histaveid, ttext(iq), ttext(iq), '-',
166!     .             iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
167!     .             32, 'ave(X)', t_ops, t_wrt)
168!        enddo
169C
170C  Masse
171C
172      call histdef(histaveid, 'masse', 'masse', 'kg',
173     .             iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
174     .             32, 'ave(X)', t_ops, t_wrt)
175C
176C  Pression au sol
177C
178      call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa',
179     .             iip1, jjp1, thoriid, 1, 1, 1, -99,
180     .             32, 'ave(X)', t_ops, t_wrt)
181C
182C  Geopotentiel au sol
183C
184!      call histdef(histaveid, 'phis', 'geopotentiel au sol', '-',
185!     .             iip1, jjp1, thoriid, 1, 1, 1, -99,
186!     .             32, 'ave(X)', t_ops, t_wrt)
187!C
188C  Fin
189C
190      call histend(histaveid)
191      call histend(histuaveid)
192      call histend(histvaveid)
193#else
194! tell the user this routine should be run with ioipsl
195      write(lunout,*)"initdynav: Warning this routine should not be",
196     &               " used without ioipsl"
197#endif
198! of #ifdef CPP_IOIPSL
199      return
200      end
Note: See TracBrowser for help on using the repository browser.