1 | ! |
---|
2 | ! $Id: initdynav_p.F 1823 2013-07-31 10:38:37Z fhourdin $ |
---|
3 | ! |
---|
4 | subroutine initdynav_p(infile,day0,anne0,tstep,t_ops,t_wrt,fileid) |
---|
5 | |
---|
6 | #ifdef CPP_IOIPSL |
---|
7 | ! This routine needs IOIPSL |
---|
8 | USE IOIPSL |
---|
9 | #endif |
---|
10 | USE parallel_lmdz |
---|
11 | use Write_field |
---|
12 | use misc_mod |
---|
13 | USE infotrac |
---|
14 | |
---|
15 | implicit none |
---|
16 | |
---|
17 | C |
---|
18 | C Routine d'initialisation des ecritures des fichiers histoires LMDZ |
---|
19 | C au format IOIPSL. Initialisation du fichier histoire moyenne. |
---|
20 | C |
---|
21 | C Appels succesifs des routines: histbeg |
---|
22 | C histhori |
---|
23 | C histver |
---|
24 | C histdef |
---|
25 | C histend |
---|
26 | C |
---|
27 | C Entree: |
---|
28 | C |
---|
29 | C infile: nom du fichier histoire a creer |
---|
30 | C day0,anne0: date de reference |
---|
31 | C tstep : frequence d'ecriture |
---|
32 | C t_ops: frequence de l'operation pour IOIPSL |
---|
33 | C t_wrt: frequence d'ecriture sur le fichier |
---|
34 | C |
---|
35 | C Sortie: |
---|
36 | C fileid: ID du fichier netcdf cree |
---|
37 | C |
---|
38 | C L. Fairhead, LMD, 03/99 |
---|
39 | C |
---|
40 | C ===================================================================== |
---|
41 | C |
---|
42 | C Declarations |
---|
43 | #include "dimensions.h" |
---|
44 | #include "paramet.h" |
---|
45 | #include "comconst.h" |
---|
46 | #include "comvert.h" |
---|
47 | #include "comgeom.h" |
---|
48 | #include "temps.h" |
---|
49 | #include "ener.h" |
---|
50 | #include "logic.h" |
---|
51 | #include "description.h" |
---|
52 | #include "serre.h" |
---|
53 | #include "iniprint.h" |
---|
54 | |
---|
55 | C Arguments |
---|
56 | C |
---|
57 | character*(*) infile |
---|
58 | integer*4 day0, anne0 |
---|
59 | real tstep, t_ops, t_wrt |
---|
60 | integer fileid |
---|
61 | |
---|
62 | #ifdef CPP_IOIPSL |
---|
63 | ! This routine needs IOIPSL |
---|
64 | C Variables locales |
---|
65 | C |
---|
66 | integer thoriid, zvertiid |
---|
67 | integer tau0 |
---|
68 | real zjulian |
---|
69 | integer iq |
---|
70 | real rlong(iip1,jjp1), rlat(iip1,jjp1) |
---|
71 | integer ii,jj |
---|
72 | integer zan, dayref |
---|
73 | integer :: jjb,jje,jjn |
---|
74 | |
---|
75 | ! definition du domaine d'ecriture pour le rebuild |
---|
76 | |
---|
77 | INTEGER,DIMENSION(2) :: ddid |
---|
78 | INTEGER,DIMENSION(2) :: dsg |
---|
79 | INTEGER,DIMENSION(2) :: dsl |
---|
80 | INTEGER,DIMENSION(2) :: dpf |
---|
81 | INTEGER,DIMENSION(2) :: dpl |
---|
82 | INTEGER,DIMENSION(2) :: dhs |
---|
83 | INTEGER,DIMENSION(2) :: dhe |
---|
84 | |
---|
85 | INTEGER :: dynave_domain_id |
---|
86 | |
---|
87 | if (adjust) return |
---|
88 | C |
---|
89 | C Initialisations |
---|
90 | C |
---|
91 | pi = 4. * atan (1.) |
---|
92 | C |
---|
93 | C Appel a histbeg: creation du fichier netcdf et initialisations diverses |
---|
94 | C |
---|
95 | |
---|
96 | zan = anne0 |
---|
97 | dayref = day0 |
---|
98 | CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) |
---|
99 | tau0 = itau_dyn |
---|
100 | |
---|
101 | do jj = 1, jjp1 |
---|
102 | do ii = 1, iip1 |
---|
103 | rlong(ii,jj) = rlonv(ii) * 180. / pi |
---|
104 | rlat(ii,jj) = rlatu(jj) * 180. / pi |
---|
105 | enddo |
---|
106 | enddo |
---|
107 | |
---|
108 | jjb=jj_begin |
---|
109 | jje=jj_end |
---|
110 | jjn=jj_nb |
---|
111 | |
---|
112 | ddid=(/ 1,2 /) |
---|
113 | dsg=(/ iip1,jjp1 /) |
---|
114 | dsl=(/ iip1,jjn /) |
---|
115 | dpf=(/ 1,jjb /) |
---|
116 | dpl=(/ iip1,jje /) |
---|
117 | dhs=(/ 0,0 /) |
---|
118 | dhe=(/ 0,0 /) |
---|
119 | |
---|
120 | call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, |
---|
121 | . 'box',dynave_domain_id) |
---|
122 | |
---|
123 | call histbeg(trim(infile),iip1, rlong(:,1), jjn, rlat(1,jjb:jje), |
---|
124 | . 1, iip1, 1, jjn,tau0, zjulian, tstep, thoriid, |
---|
125 | . fileid,dynave_domain_id) |
---|
126 | |
---|
127 | C |
---|
128 | C Appel a histvert pour la grille verticale |
---|
129 | C |
---|
130 | call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', |
---|
131 | . llm, nivsigs, zvertiid) |
---|
132 | C |
---|
133 | C Appels a histdef pour la definition des variables a sauvegarder |
---|
134 | C |
---|
135 | C Vents U |
---|
136 | C |
---|
137 | write(6,*)'inithistave',tstep |
---|
138 | call histdef(fileid, 'u', 'vents u scalaires moyennes', |
---|
139 | . 'm/s', iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
140 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
141 | |
---|
142 | C |
---|
143 | C Vents V |
---|
144 | C |
---|
145 | call histdef(fileid, 'v', 'vents v scalaires moyennes', |
---|
146 | . 'm/s', iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
147 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
148 | |
---|
149 | C |
---|
150 | C Temperature |
---|
151 | C |
---|
152 | call histdef(fileid, 'temp', 'temperature moyennee', 'K', |
---|
153 | . iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
154 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
155 | C |
---|
156 | C Temperature potentielle |
---|
157 | C |
---|
158 | call histdef(fileid, 'theta', 'temperature potentielle', 'K', |
---|
159 | . iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
160 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
161 | |
---|
162 | |
---|
163 | C |
---|
164 | C Geopotentiel |
---|
165 | C |
---|
166 | call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', |
---|
167 | . iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
168 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
169 | C |
---|
170 | C Traceurs |
---|
171 | C |
---|
172 | DO iq=1,nqtot |
---|
173 | call histdef(fileid, ttext(iq), ttext(iq), '-', |
---|
174 | . iip1, jjn, thoriid, llm, 1, llm, zvertiid, |
---|
175 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
176 | enddo |
---|
177 | C |
---|
178 | C Masse |
---|
179 | C |
---|
180 | call histdef(fileid, 'masse', 'masse', 'kg', |
---|
181 | . iip1, jjn, thoriid, 1, 1, 1, -99, |
---|
182 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
183 | C |
---|
184 | C Pression au sol |
---|
185 | C |
---|
186 | call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', |
---|
187 | . iip1, jjn, thoriid, 1, 1, 1, -99, |
---|
188 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
189 | C |
---|
190 | C Pression au sol |
---|
191 | C |
---|
192 | call histdef(fileid, 'phis', 'geopotentiel au sol', '-', |
---|
193 | . iip1, jjn, thoriid, 1, 1, 1, -99, |
---|
194 | . 32, 'ave(X)', t_ops, t_wrt) |
---|
195 | C |
---|
196 | C Fin |
---|
197 | C |
---|
198 | call histend(fileid) |
---|
199 | #else |
---|
200 | write(lunout,*)'initdynav_p: Needs IOIPSL to function' |
---|
201 | #endif |
---|
202 | ! #endif of #ifdef CPP_IOIPSL |
---|
203 | return |
---|
204 | end |
---|