1 | ! |
---|
2 | ! $Id: inithist.f90 5285 2024-10-28 13:33:29Z fairhead $ |
---|
3 | ! |
---|
4 | subroutine inithist(day0,anne0,tstep,t_ops,t_wrt) |
---|
5 | |
---|
6 | USE iniprint_mod_h |
---|
7 | USE comgeom_mod_h |
---|
8 | USE IOIPSL |
---|
9 | USE infotrac, ONLY : nqtot |
---|
10 | use com_io_dyn_mod, only : histid,histvid,histuid, & |
---|
11 | dynhist_file,dynhistv_file,dynhistu_file |
---|
12 | USE comconst_mod, ONLY: pi |
---|
13 | USE comvert_mod, ONLY: presnivs |
---|
14 | USE temps_mod, ONLY: itau_dyn |
---|
15 | |
---|
16 | USE dimensions_mod, ONLY: iim, jjm, llm, ndm |
---|
17 | USE paramet_mod_h |
---|
18 | implicit none |
---|
19 | |
---|
20 | ! |
---|
21 | ! Routine d'initialisation des ecritures des fichiers histoires LMDZ |
---|
22 | ! au format IOIPSL |
---|
23 | ! |
---|
24 | ! Appels succesifs des routines: histbeg |
---|
25 | ! histhori |
---|
26 | ! histver |
---|
27 | ! histdef |
---|
28 | ! histend |
---|
29 | ! |
---|
30 | ! Entree: |
---|
31 | ! |
---|
32 | ! infile: nom du fichier histoire a creer |
---|
33 | ! day0,anne0: date de reference |
---|
34 | ! tstep: duree du pas de temps en seconde |
---|
35 | ! t_ops: frequence de l'operation pour IOIPSL |
---|
36 | ! t_wrt: frequence d'ecriture sur le fichier |
---|
37 | ! nq: nombre de traceurs |
---|
38 | ! |
---|
39 | ! |
---|
40 | ! L. Fairhead, LMD, 03/99 |
---|
41 | ! |
---|
42 | ! ===================================================================== |
---|
43 | ! |
---|
44 | ! Declarations |
---|
45 | |
---|
46 | |
---|
47 | ! Arguments |
---|
48 | ! |
---|
49 | integer :: day0, anne0 |
---|
50 | real :: tstep, t_ops, t_wrt |
---|
51 | |
---|
52 | ! This routine needs IOIPSL to work |
---|
53 | ! Variables locales |
---|
54 | ! |
---|
55 | integer :: tau0 |
---|
56 | real :: zjulian |
---|
57 | integer :: iq |
---|
58 | real :: rlong(iip1,jjp1), rlat(iip1,jjp1) |
---|
59 | integer :: uhoriid, vhoriid, thoriid, zvertiid |
---|
60 | integer :: ii,jj |
---|
61 | integer :: zan, dayref |
---|
62 | ! |
---|
63 | ! Initialisations |
---|
64 | ! |
---|
65 | pi = 4. * atan (1.) |
---|
66 | ! |
---|
67 | ! Appel a histbeg: creation du fichier netcdf et initialisations diverses |
---|
68 | ! |
---|
69 | |
---|
70 | zan = anne0 |
---|
71 | dayref = day0 |
---|
72 | CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) |
---|
73 | tau0 = itau_dyn |
---|
74 | |
---|
75 | ! ------------------------------------------------------------- |
---|
76 | ! Creation des 3 fichiers pour les grilles horizontales U,V,Scal |
---|
77 | ! ------------------------------------------------------------- |
---|
78 | !Grille U |
---|
79 | do jj = 1, jjp1 |
---|
80 | do ii = 1, iip1 |
---|
81 | rlong(ii,jj) = rlonu(ii) * 180. / pi |
---|
82 | rlat(ii,jj) = rlatu(jj) * 180. / pi |
---|
83 | enddo |
---|
84 | enddo |
---|
85 | |
---|
86 | call histbeg(dynhistu_file, iip1, rlong(:,1), jjp1, rlat(1,:), & |
---|
87 | 1, iip1, 1, jjp1, & |
---|
88 | tau0, zjulian, tstep, uhoriid, histuid) |
---|
89 | |
---|
90 | ! Grille V |
---|
91 | do jj = 1, jjm |
---|
92 | do ii = 1, iip1 |
---|
93 | rlong(ii,jj) = rlonv(ii) * 180. / pi |
---|
94 | rlat(ii,jj) = rlatv(jj) * 180. / pi |
---|
95 | enddo |
---|
96 | enddo |
---|
97 | |
---|
98 | call histbeg(dynhistv_file, iip1, rlong(:,1), jjm, rlat(1,:), & |
---|
99 | 1, iip1, 1, jjm, & |
---|
100 | tau0, zjulian, tstep, vhoriid, histvid) |
---|
101 | |
---|
102 | !Grille Scalaire |
---|
103 | do jj = 1, jjp1 |
---|
104 | do ii = 1, iip1 |
---|
105 | rlong(ii,jj) = rlonv(ii) * 180. / pi |
---|
106 | rlat(ii,jj) = rlatu(jj) * 180. / pi |
---|
107 | enddo |
---|
108 | enddo |
---|
109 | |
---|
110 | call histbeg(dynhist_file, iip1, rlong(:,1), jjp1, rlat(1,:), & |
---|
111 | 1, iip1, 1, jjp1, & |
---|
112 | tau0, zjulian, tstep, thoriid, histid) |
---|
113 | ! ------------------------------------------------------------- |
---|
114 | ! Appel a histvert pour la grille verticale |
---|
115 | ! ------------------------------------------------------------- |
---|
116 | call histvert(histid, 'presnivs', 'Niveaux pression','mb', & |
---|
117 | llm, presnivs/100., zvertiid,'down') |
---|
118 | call histvert(histvid, 'presnivs', 'Niveaux pression','mb', & |
---|
119 | llm, presnivs/100., zvertiid,'down') |
---|
120 | call histvert(histuid, 'presnivs', 'Niveaux pression','mb', & |
---|
121 | llm, presnivs/100., zvertiid,'down') |
---|
122 | ! |
---|
123 | ! ------------------------------------------------------------- |
---|
124 | ! Appels a histdef pour la definition des variables a sauvegarder |
---|
125 | ! ------------------------------------------------------------- |
---|
126 | ! |
---|
127 | ! Vents U |
---|
128 | ! |
---|
129 | call histdef(histuid, 'u', 'vent u', 'm/s', & |
---|
130 | iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, & |
---|
131 | 32, 'inst(X)', t_ops, t_wrt) |
---|
132 | ! |
---|
133 | ! Vents V |
---|
134 | ! |
---|
135 | call histdef(histvid, 'v', 'vent v', 'm/s', & |
---|
136 | iip1, jjm, vhoriid, llm, 1, llm, zvertiid, & |
---|
137 | 32, 'inst(X)', t_ops, t_wrt) |
---|
138 | |
---|
139 | ! |
---|
140 | ! Temperature potentielle |
---|
141 | ! |
---|
142 | call histdef(histid, 'teta', 'temperature potentielle', '-', & |
---|
143 | iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
---|
144 | 32, 'inst(X)', t_ops, t_wrt) |
---|
145 | ! |
---|
146 | ! Geopotentiel |
---|
147 | ! |
---|
148 | call histdef(histid, 'phi', 'geopotentiel', '-', & |
---|
149 | iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
---|
150 | 32, 'inst(X)', t_ops, t_wrt) |
---|
151 | ! |
---|
152 | ! Traceurs |
---|
153 | ! |
---|
154 | ! |
---|
155 | ! DO iq=1,nqtot |
---|
156 | ! call histdef(histid, tracers(iq)%name, |
---|
157 | ! tracers(iq)%longName, '-', |
---|
158 | ! . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, |
---|
159 | ! . 32, 'inst(X)', t_ops, t_wrt) |
---|
160 | ! enddo |
---|
161 | !C |
---|
162 | ! Masse |
---|
163 | ! |
---|
164 | call histdef(histid, 'masse', 'masse', 'kg', & |
---|
165 | iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
---|
166 | 32, 'inst(X)', t_ops, t_wrt) |
---|
167 | ! |
---|
168 | ! Pression au sol |
---|
169 | ! |
---|
170 | call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', & |
---|
171 | iip1, jjp1, thoriid, 1, 1, 1, -99, & |
---|
172 | 32, 'inst(X)', t_ops, t_wrt) |
---|
173 | ! |
---|
174 | ! Geopotentiel au sol |
---|
175 | !C |
---|
176 | ! call histdef(histid, 'phis', 'geopotentiel au sol', '-', |
---|
177 | ! . iip1, jjp1, thoriid, 1, 1, 1, -99, |
---|
178 | ! . 32, 'inst(X)', t_ops, t_wrt) |
---|
179 | !C |
---|
180 | ! Fin |
---|
181 | ! |
---|
182 | call histend(histid) |
---|
183 | call histend(histuid) |
---|
184 | call histend(histvid) |
---|
185 | |
---|
186 | |
---|
187 | return |
---|
188 | end subroutine inithist |
---|