1 |
guez |
3 |
module inithist_m |
2 |
|
|
|
3 |
|
|
implicit none |
4 |
|
|
|
5 |
|
|
contains |
6 |
|
|
|
7 |
guez |
56 |
subroutine inithist(day0, anne0, tstep, nq, t_ops, t_wrt) |
8 |
guez |
3 |
|
9 |
guez |
39 |
! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05 |
10 |
guez |
3 |
|
11 |
guez |
39 |
! Routine d'initialisation des écritures des fichiers histoires LMDZ |
12 |
|
|
! au format IOIPSL |
13 |
|
|
! Appels successifs des routines : histbeg, histhori, histver, |
14 |
|
|
! histdef, histend |
15 |
guez |
3 |
|
16 |
guez |
39 |
! Entrées : |
17 |
|
|
! day0, anne0: date de référence |
18 |
|
|
! tstep : durée du pas de temps en secondes |
19 |
|
|
! t_ops : fréquence de l'opération pour IOIPSL |
20 |
|
|
! t_wrt : fréquence d'écriture sur le fichier |
21 |
|
|
! nq : nombre de traceurs |
22 |
guez |
3 |
|
23 |
guez |
39 |
! L. Fairhead, LMD, 03/99 |
24 |
guez |
3 |
|
25 |
guez |
30 |
USE calendar |
26 |
guez |
56 |
use com_io_dyn, only: histid, histvid, histuid |
27 |
guez |
30 |
use histcom |
28 |
guez |
3 |
use dimens_m |
29 |
|
|
use paramet_m |
30 |
|
|
use comconst |
31 |
|
|
use comvert |
32 |
|
|
use logic |
33 |
|
|
use comgeom |
34 |
|
|
use serre |
35 |
|
|
use temps |
36 |
|
|
use ener |
37 |
guez |
18 |
use iniadvtrac_m |
38 |
guez |
39 |
use nr_util, only: pi |
39 |
guez |
3 |
|
40 |
guez |
39 |
! Arguments |
41 |
guez |
3 |
integer day0, anne0 |
42 |
|
|
real, intent(in):: tstep, t_ops, t_wrt |
43 |
|
|
integer nq |
44 |
|
|
|
45 |
guez |
39 |
! Variables locales |
46 |
guez |
3 |
real zjulian |
47 |
|
|
integer iq |
48 |
guez |
39 |
real rlong(iip1, jjp1), rlat(iip1, jjp1) |
49 |
guez |
3 |
integer uhoriid, vhoriid, thoriid, zvertiid |
50 |
guez |
39 |
integer ii, jj |
51 |
guez |
3 |
integer zan, dayref |
52 |
|
|
|
53 |
|
|
!----------------------------------------------------------------------- |
54 |
|
|
|
55 |
guez |
39 |
! Appel a histbeg: creation du fichier netcdf et initialisations diverses |
56 |
guez |
3 |
|
57 |
|
|
zan = anne0 |
58 |
|
|
dayref = day0 |
59 |
|
|
CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) |
60 |
|
|
|
61 |
|
|
do jj = 1, jjp1 |
62 |
|
|
do ii = 1, iip1 |
63 |
guez |
39 |
rlong(ii, jj) = rlonu(ii) * 180. / pi |
64 |
|
|
rlat(ii, jj) = rlatu(jj) * 180. / pi |
65 |
guez |
3 |
enddo |
66 |
|
|
enddo |
67 |
|
|
|
68 |
guez |
56 |
call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, & |
69 |
|
|
jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid) |
70 |
guez |
3 |
|
71 |
|
|
do jj = 1, jjm |
72 |
|
|
do ii = 1, iip1 |
73 |
guez |
39 |
rlong(ii, jj) = rlonv(ii) * 180. / pi |
74 |
|
|
rlat(ii, jj) = rlatv(jj) * 180. / pi |
75 |
guez |
3 |
enddo |
76 |
|
|
enddo |
77 |
|
|
|
78 |
guez |
56 |
! Creation du fichier histoire pour la grille en V (oblige pour l'instant, |
79 |
|
|
! IOIPSL ne permet pas de grilles avec des nombres de point differents dans |
80 |
|
|
! un meme fichier) |
81 |
|
|
|
82 |
guez |
39 |
call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), & |
83 |
guez |
3 |
1, iip1, 1, jjm, & |
84 |
guez |
56 |
itau_dyn, zjulian, tstep, vhoriid, histvid) |
85 |
guez |
3 |
! |
86 |
guez |
39 |
! Appel a histhori pour rajouter les autres grilles horizontales |
87 |
guez |
3 |
! |
88 |
|
|
do jj = 1, jjp1 |
89 |
|
|
do ii = 1, iip1 |
90 |
guez |
39 |
rlong(ii, jj) = rlonv(ii) * 180. / pi |
91 |
|
|
rlat(ii, jj) = rlatu(jj) * 180. / pi |
92 |
guez |
3 |
enddo |
93 |
|
|
enddo |
94 |
|
|
|
95 |
guez |
56 |
call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, & |
96 |
|
|
jjp1, itau_dyn, zjulian, tstep, thoriid, histid) |
97 |
guez |
3 |
! |
98 |
guez |
39 |
! Appel a histvert pour la grille verticale |
99 |
guez |
3 |
! |
100 |
guez |
56 |
call histvert(histid, 'presnivs', 'Niveaux pression','mb', llm, & |
101 |
|
|
presnivs/100., zvertiid,'down') |
102 |
|
|
call histvert(histvid, 'presnivs', 'Niveaux pression','mb', llm, & |
103 |
|
|
presnivs/100., zvertiid,'down') |
104 |
|
|
call histvert(histuid, 'presnivs', 'Niveaux pression','mb', llm, & |
105 |
|
|
presnivs/100., zvertiid,'down') |
106 |
guez |
3 |
! |
107 |
guez |
39 |
! Appels a histdef pour la definition des variables a sauvegarder |
108 |
guez |
3 |
! |
109 |
guez |
39 |
! Vents U |
110 |
guez |
3 |
! |
111 |
guez |
56 |
call histdef(histuid, 'u', 'vent u', 'm/s', & |
112 |
guez |
3 |
iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, & |
113 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
114 |
guez |
3 |
! |
115 |
guez |
39 |
! Vents V |
116 |
guez |
3 |
! |
117 |
guez |
56 |
call histdef(histvid, 'v', 'vent v', 'm/s', & |
118 |
guez |
3 |
iip1, jjm, vhoriid, llm, 1, llm, zvertiid, & |
119 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
120 |
guez |
3 |
|
121 |
|
|
! |
122 |
guez |
39 |
! Temperature potentielle |
123 |
guez |
3 |
! |
124 |
guez |
56 |
call histdef(histid, 'teta', 'temperature potentielle', '-', & |
125 |
guez |
3 |
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
126 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
127 |
guez |
3 |
! |
128 |
guez |
39 |
! Geopotentiel |
129 |
guez |
3 |
! |
130 |
guez |
56 |
call histdef(histid, 'phi', 'geopotentiel', '-', & |
131 |
guez |
3 |
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
132 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
133 |
guez |
3 |
! |
134 |
guez |
39 |
! Traceurs |
135 |
guez |
3 |
! |
136 |
guez |
39 |
DO iq=1, nq |
137 |
guez |
56 |
call histdef(histid, ttext(iq), ttext(iq), '-', & |
138 |
guez |
3 |
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
139 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
140 |
guez |
3 |
enddo |
141 |
|
|
! |
142 |
guez |
39 |
! Masse |
143 |
guez |
3 |
! |
144 |
guez |
56 |
call histdef(histid, 'masse', 'masse', 'kg', & |
145 |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
146 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
147 |
guez |
3 |
! |
148 |
guez |
39 |
! Pression au sol |
149 |
guez |
3 |
! |
150 |
guez |
56 |
call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', & |
151 |
guez |
3 |
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
152 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
153 |
guez |
3 |
! |
154 |
guez |
56 |
! Geopotentiel au sol |
155 |
guez |
3 |
! |
156 |
guez |
56 |
call histdef(histid, 'phis', 'geopotentiel au sol', '-', & |
157 |
guez |
3 |
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
158 |
guez |
15 |
'inst(X)', t_ops, t_wrt) |
159 |
guez |
3 |
! |
160 |
guez |
39 |
! Fin |
161 |
guez |
3 |
! |
162 |
guez |
56 |
call histend(histid) |
163 |
|
|
call histend(histuid) |
164 |
|
|
call histend(histvid) |
165 |
guez |
3 |
|
166 |
|
|
end subroutine inithist |
167 |
|
|
|
168 |
|
|
end module inithist_m |