/[lmdze]/trunk/libf/bibio/inithist.f90
ViewVC logotype

Contents of /trunk/libf/bibio/inithist.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 27 - (show annotations)
Thu Mar 25 14:29:07 2010 UTC (14 years, 1 month ago) by guez
File size: 4888 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 inithist_m
2
3 ! This module is clean: no C preprocessor directive, no include line
4
5 implicit none
6
7 contains
8
9 subroutine inithist(day0, anne0, tstep, nq, fileid, filevid, t_ops, &
10 t_wrt)
11
12 ! From inithist.F,v 1.1.1.1 2004/05/19 12:53:05
13
14 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
15 ! au format IOIPSL
16
17 ! Appels succesifs des routines: histbeg
18 ! histhori
19 ! histver
20 ! histdef
21 ! histend
22
23 ! Entree:
24 ! day0,anne0: date de reference
25 ! tstep: duree du pas de temps en seconde
26 ! t_ops: frequence de l'operation pour IOIPSL
27 ! t_wrt: frequence d'ecriture sur le fichier
28 ! nq: nombre de traceurs
29
30 ! Sortie:
31 ! fileid: ID du fichier netcdf cree
32 ! filevid:ID du fichier netcdf pour la grille v
33
34 ! L. Fairhead, LMD, 03/99
35
36 USE IOIPSL
37 use dimens_m
38 use paramet_m
39 use comconst
40 use comvert
41 use logic
42 use comgeom
43 use serre
44 use temps
45 use ener
46 use iniadvtrac_m
47
48 ! Arguments
49 integer day0, anne0
50 real, intent(in):: tstep, t_ops, t_wrt
51 integer fileid, filevid
52 integer nq
53
54 ! Variables locales
55 real zjulian
56 integer iq
57 real rlong(iip1,jjp1), rlat(iip1,jjp1)
58 integer uhoriid, vhoriid, thoriid, zvertiid
59 integer ii,jj
60 integer zan, dayref
61
62 !-----------------------------------------------------------------------
63
64 ! Initialisations
65
66 pi = 4. * atan (1.)
67
68 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
69
70 zan = anne0
71 dayref = day0
72 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
73
74 do jj = 1, jjp1
75 do ii = 1, iip1
76 rlong(ii,jj) = rlonu(ii) * 180. / pi
77 rlat(ii,jj) = rlatu(jj) * 180. / pi
78 enddo
79 enddo
80
81 call histbeg_totreg("dyn_hist.nc", rlong(:,1), rlat(1,:), &
82 1, iip1, 1, jjp1, &
83 itau_dyn, zjulian, tstep, uhoriid, fileid)
84 !
85 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
86 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
87 ! un meme fichier)
88
89 do jj = 1, jjm
90 do ii = 1, iip1
91 rlong(ii,jj) = rlonv(ii) * 180. / pi
92 rlat(ii,jj) = rlatv(jj) * 180. / pi
93 enddo
94 enddo
95
96 call histbeg_totreg('dyn_histv.nc', rlong(:,1), rlat(1,:jjm), &
97 1, iip1, 1, jjm, &
98 itau_dyn, zjulian, tstep, vhoriid, filevid)
99 !
100 ! Appel a histhori pour rajouter les autres grilles horizontales
101 !
102 do jj = 1, jjp1
103 do ii = 1, iip1
104 rlong(ii,jj) = rlonv(ii) * 180. / pi
105 rlat(ii,jj) = rlatu(jj) * 180. / pi
106 enddo
107 enddo
108
109 call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &
110 'Grille points scalaires', thoriid)
111 !
112 ! Appel a histvert pour la grille verticale
113 !
114 call histvert(fileid, 'sig_s', 'Niveaux sigma','-', &
115 llm, nivsigs, zvertiid)
116 ! Pour le fichier V
117 call histvert(filevid, 'sig_s', 'Niveaux sigma','-', &
118 llm, nivsigs, zvertiid)
119 !
120 ! Appels a histdef pour la definition des variables a sauvegarder
121 !
122 ! Vents U
123 !
124 call histdef(fileid, 'ucov', 'vents u covariants', 'm/s', &
125 iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
126 'inst(X)', t_ops, t_wrt)
127 !
128 ! Vents V
129 !
130 call histdef(filevid, 'vcov', 'vents v covariants', 'm/s', &
131 iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
132 'inst(X)', t_ops, t_wrt)
133
134 !
135 ! Temperature potentielle
136 !
137 call histdef(fileid, 'teta', 'temperature potentielle', '-', &
138 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
139 'inst(X)', t_ops, t_wrt)
140 !
141 ! Geopotentiel
142 !
143 call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &
144 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
145 'inst(X)', t_ops, t_wrt)
146 !
147 ! Traceurs
148 !
149 DO iq=1,nq
150 call histdef(fileid, ttext(iq), ttext(iq), '-', &
151 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
152 'inst(X)', t_ops, t_wrt)
153 enddo
154 !
155 ! Masse
156 !
157 call histdef(fileid, 'masse', 'masse', 'kg', &
158 iip1, jjp1, thoriid, 1, 1, 1, -99, &
159 'inst(X)', t_ops, t_wrt)
160 !
161 ! Pression au sol
162 !
163 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
164 iip1, jjp1, thoriid, 1, 1, 1, -99, &
165 'inst(X)', t_ops, t_wrt)
166 !
167 ! Pression au sol
168 !
169 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
170 iip1, jjp1, thoriid, 1, 1, 1, -99, &
171 'inst(X)', t_ops, t_wrt)
172 !
173 ! Fin
174 !
175 call histend(fileid)
176 call histend(filevid)
177
178 end subroutine inithist
179
180 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21