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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 39 - (show annotations)
Tue Jan 25 15:11:05 2011 UTC (13 years, 3 months ago) by guez
File size: 4708 byte(s)
"pi" comes from "nr_util". Removed subroutine "initialize" in module
"comconst".

Copied the content of "fxy_sin.h" into "fxysinus", instead of getting
it from an "include" line. Removed file "fxy_sin.h".

"ps" has rank 2 in "gcm" and "dynetat0".

Assumed-shape for argument "q" of "integrd".

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

  ViewVC Help
Powered by ViewVC 1.1.21