1 |
module initdynav_m |
module initdynav_m |
2 |
|
|
|
! This module is clean: no C preprocessor directive, no include line |
|
|
|
|
3 |
implicit none |
implicit none |
4 |
|
|
5 |
|
integer histaveid |
6 |
|
|
7 |
contains |
contains |
8 |
|
|
9 |
subroutine initdynav(day0, anne0, tstep, nq, fileid, infile, t_ops, t_wrt) |
subroutine initdynav(day0, anne0, tstep, nq, t_ops, t_wrt) |
10 |
|
|
11 |
|
! From initdynav.F, version 1.1.1.1, 2004/05/19 12:53:05 |
12 |
|
! L. Fairhead, LMD |
13 |
|
|
14 |
! From initdynav.F,v 1.1.1.1 2004/05/19 12:53:05 |
! Routine d'initialisation des écritures des fichiers histoires au |
15 |
|
! format IOIPSL. Initialisation du fichier histoire moyenne. |
16 |
|
|
17 |
USE IOIPSL, only: ymds2ju, histbeg_totreg, histvert, histdef, histend |
use calendar, ONLY: ymds2ju |
18 |
use dimens_m |
USE comgeom, ONLY: rlatu, rlonv |
19 |
use paramet_m |
USE dimens_m, ONLY: llm |
20 |
use comconst, only: pi |
USE disvert_m, ONLY: nivsigs |
21 |
use comvert, only: nivsigs |
USE histbeg_totreg_m, ONLY: histbeg_totreg |
22 |
use logic |
USE histdef_m, ONLY: histdef |
23 |
use comgeom |
USE histend_m, ONLY: histend |
24 |
use serre |
USE histvert_m, ONLY: histvert |
25 |
use temps |
USE iniadvtrac_m, ONLY: ttext |
26 |
use ener |
USE nr_util, ONLY: pi |
27 |
use advtrac_m, only: ttext |
USE paramet_m, ONLY: iip1, jjp1 |
28 |
|
USE temps, ONLY: itau_dyn |
29 |
! Routine d'initialisation des ecritures des fichiers histoires LMDZ |
|
30 |
! au format IOIPSL. Initialisation du fichier histoire moyenne. |
integer, intent(in):: day0, anne0 ! date de référence |
31 |
|
real, intent(in):: tstep ! fréquence d'écriture |
32 |
! Appels succesifs des routines: histbeg |
integer, intent(in):: nq ! nombre de traceurs |
33 |
! histhori |
real, intent(in):: t_ops ! fréquence de l'opération pour IOIPSL |
34 |
! histver |
real, intent(in):: t_wrt ! fréquence d'écriture sur le fichier |
35 |
! histdef |
|
36 |
! histend |
! Variables locales |
37 |
|
integer horiid, zvertiid |
38 |
! Entree: |
real julian |
39 |
! infile: nom du fichier histoire a creer |
integer iq, ii, jj |
|
! day0,anne0: date de reference |
|
|
! tstep : frequence d'ecriture |
|
|
! t_ops: frequence de l'operation pour IOIPSL |
|
|
! t_wrt: frequence d'ecriture sur le fichier |
|
|
! nq: nombre de traceurs |
|
|
|
|
|
! Sortie: |
|
|
! fileid: ID du fichier netcdf cree |
|
|
|
|
|
! L. Fairhead, LMD, 03/99 |
|
|
|
|
|
! Arguments |
|
|
character(len=*) infile |
|
|
integer day0, anne0 |
|
|
real, intent(in):: tstep, t_ops, t_wrt |
|
|
integer fileid |
|
|
integer nq |
|
|
integer thoriid, zvertiid |
|
|
|
|
|
! Variables locales |
|
|
|
|
|
integer tau0 |
|
|
real zjulian |
|
|
integer iq |
|
|
real rlong(iip1,jjp1), rlat(iip1,jjp1) |
|
|
integer ii,jj |
|
|
integer zan, dayref |
|
40 |
|
|
41 |
!---------------------------------------------------- |
!---------------------------------------------------- |
42 |
|
|
43 |
! Initialisations |
print *, "Call sequence information: initdynav" |
44 |
|
|
45 |
pi = 4. * atan (1.) |
CALL ymds2ju(anne0, 1, day0, 0., julian) |
46 |
! |
call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, & |
47 |
! Appel a histbeg: creation du fichier netcdf et initialisations diverses |
rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, & |
48 |
! |
horiid, histaveid) |
49 |
|
call histvert(histaveid, 'sigss', 'Niveaux sigma', 'Pa', nivsigs, zvertiid) |
50 |
zan = anne0 |
|
51 |
dayref = day0 |
call histdef(histaveid, 'u', 'vents u scalaires moyennes', 'm/s', iip1, & |
52 |
CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) |
jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
53 |
tau0 = itau_dyn |
call histdef(histaveid, 'v', 'vents v scalaires moyennes', 'm/s', iip1, & |
54 |
|
jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
55 |
do jj = 1, jjp1 |
call histdef(histaveid, 'temp', 'temperature moyennee', 'K', iip1, jjp1, & |
56 |
do ii = 1, iip1 |
horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
57 |
rlong(ii,jj) = rlonv(ii) * 180. / pi |
call histdef(histaveid, 'theta', 'temperature potentielle', 'K', iip1, & |
58 |
rlat(ii,jj) = rlatu(jj) * 180. / pi |
jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
59 |
enddo |
call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', iip1, jjp1, & |
60 |
|
horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
61 |
|
|
62 |
|
! Traceurs |
63 |
|
DO iq = 1, nq |
64 |
|
call histdef(histaveid, ttext(iq), ttext(iq), '-', iip1, jjp1, & |
65 |
|
horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt) |
66 |
enddo |
enddo |
67 |
|
|
68 |
call histbeg_totreg(infile, iip1, rlong(:,1), jjp1, rlat(1,:), & |
call histdef(histaveid, 'masse', 'masse', 'kg', iip1, jjp1, horiid, 1, & |
69 |
1, iip1, 1, jjp1, & |
1, 1, -99, 'ave(X)', t_ops, t_wrt) |
70 |
tau0, zjulian, tstep, thoriid, fileid) |
call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', iip1, & |
71 |
|
jjp1, horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt) |
72 |
! |
call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, & |
73 |
! Appel a histvert pour la grille verticale |
horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt) |
74 |
! |
|
75 |
call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', & |
call histend(histaveid) |
|
llm, nivsigs, zvertiid) |
|
|
! |
|
|
! Appels a histdef pour la definition des variables a sauvegarder |
|
|
! |
|
|
! Vents U |
|
|
! |
|
|
write(6,*)'inithistave',tstep |
|
|
call histdef(fileid, 'u', 'vents u scalaires moyennes', & |
|
|
'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
|
|
|
! |
|
|
! Vents V |
|
|
! |
|
|
call histdef(fileid, 'v', 'vents v scalaires moyennes', & |
|
|
'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
|
|
|
! |
|
|
! Temperature |
|
|
! |
|
|
call histdef(fileid, 'temp', 'temperature moyennee', 'K', & |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
! |
|
|
! Temperature potentielle |
|
|
! |
|
|
call histdef(fileid, 'theta', 'temperature potentielle', 'K', & |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
|
|
|
|
|
|
! |
|
|
! Geopotentiel |
|
|
! |
|
|
call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', & |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
! |
|
|
! Traceurs |
|
|
! |
|
|
DO iq=1,nq |
|
|
call histdef(fileid, ttext(iq), ttext(iq), '-', & |
|
|
iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
enddo |
|
|
! |
|
|
! Masse |
|
|
! |
|
|
call histdef(fileid, 'masse', 'masse', 'kg', & |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
! |
|
|
! Pression au sol |
|
|
! |
|
|
call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', & |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
! |
|
|
! Pression au sol |
|
|
! |
|
|
call histdef(fileid, 'phis', 'geopotentiel au sol', '-', & |
|
|
iip1, jjp1, thoriid, 1, 1, 1, -99, & |
|
|
32, 'ave(X)', t_ops, t_wrt) |
|
|
! |
|
|
! Fin |
|
|
! |
|
|
call histend(fileid) |
|
76 |
|
|
77 |
end subroutine initdynav |
end subroutine initdynav |
78 |
|
|