1 |
module inithist_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
integer histid, histvid, histuid |
6 |
|
7 |
contains |
8 |
|
9 |
subroutine inithist(t_ops, t_wrt) |
10 |
|
11 |
! From inithist.F, version 1.1.1.1, 2004/05/19 12:53:05 |
12 |
! L. Fairhead, LMD, 03/99 |
13 |
|
14 |
! Routine d'initialisation des écritures des fichiers histoires au |
15 |
! format IOIPSL. |
16 |
|
17 |
use comconst, only: dtvr |
18 |
USE dimensions, ONLY: jjm, llm, nqmx |
19 |
USE disvert_m, ONLY: presnivs |
20 |
use dynetat0_m, only: day_ref, annee_ref, rlatu, rlatv, rlonu, rlonv |
21 |
use histbeg_totreg_m, only: histbeg_totreg |
22 |
USE histdef_m, ONLY: histdef |
23 |
USE histend_m, ONLY: histend |
24 |
USE histvert_m, ONLY: histvert |
25 |
USE iniadvtrac_m, ONLY: ttext |
26 |
USE nr_util, ONLY: pi |
27 |
USE paramet_m, ONLY: iip1, jjp1 |
28 |
USE temps, ONLY: itau_dyn |
29 |
use ymds2ju_m, ONLY: ymds2ju |
30 |
|
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 |
! Local: |
35 |
real julian |
36 |
integer iq |
37 |
integer uhoriid, vhoriid, thoriid, zvertiid |
38 |
|
39 |
!----------------------------------------------------------------------- |
40 |
|
41 |
print *, "Call sequence information: inithist" |
42 |
CALL ymds2ju(annee_ref, 1, day_ref, 0., julian) |
43 |
|
44 |
call histbeg_totreg("dyn_histu.nc", rlonu * 180. / pi, rlatu * 180. / pi, & |
45 |
1, iip1, 1, jjp1, itau_dyn, julian, dtvr, uhoriid, histuid) |
46 |
|
47 |
! Creation du fichier histoire pour la grille en V (oblige pour l'instant, |
48 |
! IOIPSL ne permet pas de grilles avec des nombres de point differents dans |
49 |
! un meme fichier) |
50 |
|
51 |
call histbeg_totreg('dyn_histv.nc', rlonv * 180. / pi, rlatv * 180. / pi, & |
52 |
1, iip1, 1, jjm, itau_dyn, julian, dtvr, vhoriid, histvid) |
53 |
|
54 |
call histbeg_totreg("dyn_hist.nc", rlonv * 180. / pi, rlatu * 180. / pi, & |
55 |
1, iip1, 1, jjp1, itau_dyn, julian, dtvr, thoriid, histid) |
56 |
|
57 |
call histvert(histid, 'presnivs', 'Niveaux pression', 'mb', presnivs/100., & |
58 |
zvertiid, 'down') |
59 |
call histvert(histvid, 'presnivs', 'Niveaux pression', 'mb', & |
60 |
presnivs/100., zvertiid, 'down') |
61 |
call histvert(histuid, 'presnivs', 'Niveaux pression', 'mb', & |
62 |
presnivs/100., zvertiid, 'down') |
63 |
|
64 |
call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, & |
65 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
66 |
call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, & |
67 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
68 |
call histdef(histid, 'temp', 'temperature', 'K', iip1, jjp1, & |
69 |
thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
70 |
call histdef(histid, 'theta', 'temperature potentielle', 'K', iip1, jjp1, & |
71 |
thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
72 |
call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, & |
73 |
llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
74 |
|
75 |
! Traceurs |
76 |
DO iq = 1, nqmx |
77 |
call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, & |
78 |
llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
79 |
enddo |
80 |
|
81 |
call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, & |
82 |
llm, zvertiid, 'inst(X)', t_ops, t_wrt) |
83 |
call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, & |
84 |
thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt) |
85 |
|
86 |
call histend(histid) |
87 |
call histend(histuid) |
88 |
call histend(histvid) |
89 |
|
90 |
end subroutine inithist |
91 |
|
92 |
end module inithist_m |