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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 56 - (show annotations)
Tue Jan 10 19:02:02 2012 UTC (12 years, 4 months ago) by guez
File size: 4669 byte(s)
Imported "writehist.f" from LMDZ.

Moved module variable "histaveid" from "com_io_dyn" to "initdynav_m".

In "inithist", access directly module variables from "com_io_dyn"
instead of going through the arguments. Copying from LMDZ, write "u"
and scalar variables to separate files. Create a new variable for the
new file in "com_io_dyn". Copying from LMDZ, change the vertical axes
of the three files.

Removed some useless initializations in "dissip".

In "bilan_dyn", removed useless variable "time". Avoiding the
approximate test on "dt_cum" being a multiple of "dt_app", just
compute "ncum" from known usage of "bilan_dyn" and compute "dt_cum"
from "ncum". Change "periodav" from real to integer in
"conf_gcm_m". Since "day_step" is required to be a multiple of
"iperiod", so is "ncum".

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

  ViewVC Help
Powered by ViewVC 1.1.21