1 |
guez |
3 |
module initdynav_m |
2 |
|
|
|
3 |
|
|
! This module is clean: no C preprocessor directive, no include line |
4 |
|
|
|
5 |
|
|
implicit none |
6 |
|
|
|
7 |
|
|
contains |
8 |
|
|
|
9 |
|
|
subroutine initdynav(day0, anne0, tstep, nq, fileid, infile, t_ops, t_wrt) |
10 |
|
|
|
11 |
|
|
! From initdynav.F,v 1.1.1.1 2004/05/19 12:53:05 |
12 |
|
|
|
13 |
|
|
USE IOIPSL, only: ymds2ju, histbeg_totreg, histvert, histdef, histend |
14 |
|
|
use dimens_m |
15 |
|
|
use paramet_m |
16 |
|
|
use comconst, only: pi |
17 |
|
|
use comvert, only: nivsigs |
18 |
|
|
use logic |
19 |
|
|
use comgeom |
20 |
|
|
use serre |
21 |
|
|
use temps |
22 |
|
|
use ener |
23 |
|
|
use advtrac_m, only: ttext |
24 |
|
|
|
25 |
|
|
! Routine d'initialisation des ecritures des fichiers histoires LMDZ |
26 |
|
|
! au format IOIPSL. Initialisation du fichier histoire moyenne. |
27 |
|
|
|
28 |
|
|
! Appels succesifs des routines: histbeg |
29 |
|
|
! histhori |
30 |
|
|
! histver |
31 |
|
|
! histdef |
32 |
|
|
! histend |
33 |
|
|
|
34 |
|
|
! Entree: |
35 |
|
|
! infile: nom du fichier histoire a creer |
36 |
|
|
! day0,anne0: date de reference |
37 |
|
|
! tstep : frequence d'ecriture |
38 |
|
|
! t_ops: frequence de l'operation pour IOIPSL |
39 |
|
|
! t_wrt: frequence d'ecriture sur le fichier |
40 |
|
|
! nq: nombre de traceurs |
41 |
|
|
|
42 |
|
|
! Sortie: |
43 |
|
|
! fileid: ID du fichier netcdf cree |
44 |
|
|
|
45 |
|
|
! L. Fairhead, LMD, 03/99 |
46 |
|
|
|
47 |
|
|
! Arguments |
48 |
|
|
character(len=*) infile |
49 |
|
|
integer day0, anne0 |
50 |
|
|
real, intent(in):: tstep, t_ops, t_wrt |
51 |
|
|
integer fileid |
52 |
|
|
integer nq |
53 |
|
|
integer thoriid, zvertiid |
54 |
|
|
|
55 |
|
|
! Variables locales |
56 |
|
|
|
57 |
|
|
integer tau0 |
58 |
|
|
real zjulian |
59 |
|
|
integer iq |
60 |
|
|
real rlong(iip1,jjp1), rlat(iip1,jjp1) |
61 |
|
|
integer ii,jj |
62 |
|
|
integer zan, dayref |
63 |
|
|
|
64 |
|
|
!---------------------------------------------------- |
65 |
|
|
|
66 |
|
|
! Initialisations |
67 |
|
|
|
68 |
|
|
pi = 4. * atan (1.) |
69 |
|
|
! |
70 |
|
|
! Appel a histbeg: creation du fichier netcdf et initialisations diverses |
71 |
|
|
! |
72 |
|
|
|
73 |
|
|
zan = anne0 |
74 |
|
|
dayref = day0 |
75 |
|
|
CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) |
76 |
|
|
tau0 = itau_dyn |
77 |
|
|
|
78 |
|
|
do jj = 1, jjp1 |
79 |
|
|
do ii = 1, iip1 |
80 |
|
|
rlong(ii,jj) = rlonv(ii) * 180. / pi |
81 |
|
|
rlat(ii,jj) = rlatu(jj) * 180. / pi |
82 |
|
|
enddo |
83 |
|
|
enddo |
84 |
|
|
|
85 |
|
|
call histbeg_totreg(infile, iip1, rlong(:,1), jjp1, rlat(1,:), & |
86 |
|
|
1, iip1, 1, jjp1, & |
87 |
|
|
tau0, zjulian, tstep, thoriid, fileid) |
88 |
|
|
|
89 |
|
|
! |
90 |
|
|
! Appel a histvert pour la grille verticale |
91 |
|
|
! |
92 |
|
|
call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', & |
93 |
|
|
llm, nivsigs, zvertiid) |
94 |
|
|
! |
95 |
|
|
! Appels a histdef pour la definition des variables a sauvegarder |
96 |
|
|
! |
97 |
|
|
! Vents U |
98 |
|
|
! |
99 |
|
|
write(6,*)'inithistave',tstep |
100 |
|
|
call histdef(fileid, 'u', 'vents u scalaires moyennes', & |
101 |
|
|
'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
102 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
103 |
|
|
|
104 |
|
|
! |
105 |
|
|
! Vents V |
106 |
|
|
! |
107 |
|
|
call histdef(fileid, 'v', 'vents v scalaires moyennes', & |
108 |
|
|
'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
109 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
110 |
|
|
|
111 |
|
|
! |
112 |
|
|
! Temperature |
113 |
|
|
! |
114 |
|
|
call histdef(fileid, 'temp', 'temperature moyennee', 'K', & |
115 |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
116 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
117 |
|
|
! |
118 |
|
|
! Temperature potentielle |
119 |
|
|
! |
120 |
|
|
call histdef(fileid, 'theta', 'temperature potentielle', 'K', & |
121 |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
122 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
123 |
|
|
|
124 |
|
|
|
125 |
|
|
! |
126 |
|
|
! Geopotentiel |
127 |
|
|
! |
128 |
|
|
call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', & |
129 |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
130 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
131 |
|
|
! |
132 |
|
|
! Traceurs |
133 |
|
|
! |
134 |
|
|
DO iq=1,nq |
135 |
|
|
call histdef(fileid, ttext(iq), ttext(iq), '-', & |
136 |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
137 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
138 |
|
|
enddo |
139 |
|
|
! |
140 |
|
|
! Masse |
141 |
|
|
! |
142 |
|
|
call histdef(fileid, 'masse', 'masse', 'kg', & |
143 |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
144 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
145 |
|
|
! |
146 |
|
|
! Pression au sol |
147 |
|
|
! |
148 |
|
|
call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', & |
149 |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
150 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
151 |
|
|
! |
152 |
|
|
! Pression au sol |
153 |
|
|
! |
154 |
|
|
call histdef(fileid, 'phis', 'geopotentiel au sol', '-', & |
155 |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
156 |
|
|
32, 'ave(X)', t_ops, t_wrt) |
157 |
|
|
! |
158 |
|
|
! Fin |
159 |
|
|
! |
160 |
|
|
call histend(fileid) |
161 |
|
|
|
162 |
|
|
end subroutine initdynav |
163 |
|
|
|
164 |
|
|
end module initdynav_m |