1 |
module inithist_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
subroutine inithist(tstep, nq, t_ops, t_wrt) |
8 |
|
9 |
! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05 |
10 |
! L. Fairhead, LMD, 03/99 |
11 |
|
12 |
! Routine d'initialisation des écritures des fichiers histoires |
13 |
! LMDZ au format IOIPSL. |
14 |
|
15 |
USE com_io_dyn, ONLY: histid, histuid, histvid |
16 |
USE dimens_m, ONLY: jjm, llm |
17 |
USE disvert_m, ONLY: presnivs |
18 |
use dynetat0_m, only: day_ref, annee_ref, rlatu, rlatv, rlonu, rlonv |
19 |
USE histbeg_totreg_m, ONLY : histbeg_totreg |
20 |
USE histdef_m, ONLY : histdef |
21 |
USE histend_m, ONLY : histend |
22 |
USE histvert_m, ONLY : histvert |
23 |
USE iniadvtrac_m, ONLY: ttext |
24 |
USE nr_util, ONLY: pi |
25 |
USE paramet_m, ONLY: iip1, jjp1 |
26 |
USE temps, ONLY: itau_dyn |
27 |
USE ymds2ju_m, ONLY: ymds2ju |
28 |
|
29 |
real, intent(in):: tstep ! durée du pas de temps en secondes |
30 |
integer, intent(in):: nq ! nombre de traceurs |
31 |
real, intent(in):: t_ops ! fréquence de l'opération pour IOIPSL |
32 |
real, intent(in):: t_wrt ! fréquence d'écriture sur le fichier |
33 |
|
34 |
! Variables locales: |
35 |
real zjulian |
36 |
integer iq |
37 |
real rlong(iip1, jjp1), rlat(iip1, jjp1) |
38 |
integer uhoriid, vhoriid, thoriid, zvertiid |
39 |
integer ii, jj |
40 |
|
41 |
!----------------------------------------------------------------------- |
42 |
|
43 |
CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian) |
44 |
|
45 |
do jj = 1, jjp1 |
46 |
do ii = 1, iip1 |
47 |
rlong(ii, jj) = rlonu(ii) * 180. / pi |
48 |
rlat(ii, jj) = rlatu(jj) * 180. / pi |
49 |
enddo |
50 |
enddo |
51 |
|
52 |
call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, & |
53 |
jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid) |
54 |
|
55 |
do jj = 1, jjm |
56 |
do ii = 1, iip1 |
57 |
rlong(ii, jj) = rlonv(ii) * 180. / pi |
58 |
rlat(ii, jj) = rlatv(jj) * 180. / pi |
59 |
enddo |
60 |
enddo |
61 |
|
62 |
! Creation du fichier histoire pour la grille en V (oblige pour l'instant, |
63 |
! IOIPSL ne permet pas de grilles avec des nombres de point differents dans |
64 |
! un meme fichier) |
65 |
|
66 |
call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, & |
67 |
1, jjm, itau_dyn, zjulian, tstep, vhoriid, histvid) |
68 |
|
69 |
do jj = 1, jjp1 |
70 |
do ii = 1, iip1 |
71 |
rlong(ii, jj) = rlonv(ii) * 180. / pi |
72 |
rlat(ii, jj) = rlatu(jj) * 180. / pi |
73 |
enddo |
74 |
enddo |
75 |
|
76 |
call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, & |
77 |
jjp1, itau_dyn, zjulian, tstep, thoriid, histid) |
78 |
|
79 |
! Appel a histvert pour la grille verticale |
80 |
|
81 |
call histvert(histid, 'presnivs', 'Niveaux pression','mb', presnivs/100., & |
82 |
zvertiid,'down') |
83 |
call histvert(histvid, 'presnivs', 'Niveaux pression','mb', & |
84 |
presnivs/100., zvertiid,'down') |
85 |
call histvert(histuid, 'presnivs', 'Niveaux pression','mb', & |
86 |
presnivs/100., zvertiid,'down') |
87 |
|
88 |
! Appels a histdef pour la definition des variables a sauvegarder |
89 |
|
90 |
call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, & |
91 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
92 |
call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, & |
93 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
94 |
call histdef(histid, 'teta', 'temperature potentielle', '-', iip1, jjp1, & |
95 |
thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
96 |
call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, & |
97 |
llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
98 |
|
99 |
! Traceurs |
100 |
DO iq=1, nq |
101 |
call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, & |
102 |
llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
103 |
enddo |
104 |
|
105 |
call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, & |
106 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
107 |
call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, & |
108 |
thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt) |
109 |
call histdef(histid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, & |
110 |
thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt) |
111 |
|
112 |
call histend(histid) |
113 |
call histend(histuid) |
114 |
call histend(histvid) |
115 |
|
116 |
end subroutine inithist |
117 |
|
118 |
end module inithist_m |