/[lmdze]/trunk/Sources/dyn3d/initdynav.f
ViewVC logotype

Contents of /trunk/Sources/dyn3d/initdynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 27 - (show annotations)
Thu Mar 25 14:29:07 2010 UTC (14 years, 2 months ago) by guez
Original Path: trunk/libf/bibio/initdynav.f90
File size: 3986 byte(s)
"dyn3d" and "filtrez" do not contain any included file so make rules
have been updated.

"comdissip.f90" was useless, removed it.

"dynredem0" wrote undefined value in "controle(31)", that was
overwritten by "dynredem1". Now "dynredem0" just writes 0 to
"controle(31)".

Removed arguments of "inidissip". "inidissip" now accesses the
variables by use association.

In program "etat0_lim", "itaufin" is not defined so "dynredem1" wrote
undefined value to "controle(31)". Added argument "itau" of
"dynredem1" to correct that.

"itaufin" does not need to be a module variable (of "temps"), made it
a local variable of "leapfrog".

Removed calls to "diagedyn" from "leapfrog".

1 module initdynav_m
2
3 ! This module is clean: no C preprocessor directive, no include line
4
5 implicit none
6
7 contains
8
9 subroutine initdynav(day0, anne0, tstep, nq, fileid, 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 USE ioipsl, ONLY : histbeg_totreg, histdef, histend, histvert, ymds2ju
15 USE dimens_m, ONLY : llm
16 USE paramet_m, ONLY : iip1, jjp1
17 USE comconst, ONLY : pi
18 USE comvert, ONLY : nivsigs
19 USE comgeom, ONLY : rlatu, rlonv
20 USE temps, ONLY : itau_dyn
21 USE iniadvtrac_m, ONLY : ttext
22
23 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
24 ! au format IOIPSL. Initialisation du fichier histoire moyenne.
25
26 ! Appels succesifs des routines: histbeg
27 ! histhori
28 ! histver
29 ! histdef
30 ! histend
31
32 ! Entree:
33 ! day0,anne0: date de reference
34 ! tstep : frequence d'ecriture
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 ! Sortie:
40 ! fileid: ID du fichier netcdf cree
41
42 ! Arguments
43 integer day0, anne0
44 real, intent(in):: tstep, t_ops, t_wrt
45 integer fileid
46 integer nq
47 integer thoriid, zvertiid
48
49 ! Variables locales
50
51 real zjulian
52 integer iq
53 real rlong(iip1,jjp1), rlat(iip1,jjp1)
54 integer ii,jj
55 integer zan, dayref
56
57 !----------------------------------------------------
58
59 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
60
61 zan = anne0
62 dayref = day0
63 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
64
65 do jj = 1, jjp1
66 do ii = 1, iip1
67 rlong(ii,jj) = rlonv(ii) * 180. / pi
68 rlat(ii,jj) = rlatu(jj) * 180. / pi
69 enddo
70 enddo
71
72 call histbeg_totreg('dyn_hist_ave.nc', rlong(:,1), rlat(1,:), &
73 1, iip1, 1, jjp1, &
74 itau_dyn, zjulian, tstep, thoriid, fileid)
75
76
77 ! Appel a histvert pour la grille verticale
78
79 call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', &
80 llm, nivsigs, zvertiid)
81
82 ! Appels a histdef pour la definition des variables a sauvegarder
83
84 ! Vents U
85
86 write(6,*)'inithistave',tstep
87 call histdef(fileid, 'u', 'vents u scalaires moyennes', &
88 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
89 'ave(X)', t_ops, t_wrt)
90
91
92 ! Vents V
93
94 call histdef(fileid, 'v', 'vents v scalaires moyennes', &
95 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
96 'ave(X)', t_ops, t_wrt)
97
98
99 ! Temperature
100
101 call histdef(fileid, 'temp', 'temperature moyennee', 'K', &
102 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
103 'ave(X)', t_ops, t_wrt)
104
105 ! Temperature potentielle
106
107 call histdef(fileid, 'theta', 'temperature potentielle', 'K', &
108 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
109 'ave(X)', t_ops, t_wrt)
110
111
112
113 ! Geopotentiel
114
115 call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', &
116 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
117 'ave(X)', t_ops, t_wrt)
118
119 ! Traceurs
120
121 DO iq=1,nq
122 call histdef(fileid, ttext(iq), ttext(iq), '-', &
123 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
124 'ave(X)', t_ops, t_wrt)
125 enddo
126
127 ! Masse
128
129 call histdef(fileid, 'masse', 'masse', 'kg', &
130 iip1, jjp1, thoriid, 1, 1, 1, -99, &
131 'ave(X)', t_ops, t_wrt)
132
133 ! Pression au sol
134
135 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
136 iip1, jjp1, thoriid, 1, 1, 1, -99, &
137 'ave(X)', t_ops, t_wrt)
138
139 ! Pression au sol
140
141 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
142 iip1, jjp1, thoriid, 1, 1, 1, -99, &
143 'ave(X)', t_ops, t_wrt)
144
145 call histend(fileid)
146
147 end subroutine initdynav
148
149 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21