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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 57 - (show annotations)
Mon Jan 30 12:54:02 2012 UTC (12 years, 3 months ago) by guez
File size: 4655 byte(s)
Write used namelists to file "" instead of standard output.

Avoid aliasing in "inidissip" in calls to "divgrad2", "divgrad",
"gradiv2", "gradiv", "nxgraro2" and "nxgrarot". Add a degenerate
dimension to arrays so they have rank 3, like the dummy arguments in
"divgrad2", "divgrad", "gradiv2", "gradiv", "nxgraro2" and "nxgrarot".

Extract the initialization part from "bilan_dyn" and make a separate
procedure, "init_dynzon", from it.

Move variables from modules "iniprint" and "logic" to module
"conf_gcm_m".

Promote internal procedures of "fxy" to private procedures of module
"fxy_m".

Extracted documentation from "inigeom". Removed useless "save"
attributes. Removed useless intermediate variables. Extracted
processing of poles from loop on latitudes. Write coordinates to file
"longitude_latitude.txt" instead of standard output.

Do not use ozone tracer for radiative transfer.

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 comgeom
33 use serre
34 use temps
35 use ener
36 use iniadvtrac_m
37 use nr_util, only: pi
38
39 ! Arguments
40 integer day0, anne0
41 real, intent(in):: tstep, t_ops, t_wrt
42 integer nq
43
44 ! Variables locales
45 real zjulian
46 integer iq
47 real rlong(iip1, jjp1), rlat(iip1, jjp1)
48 integer uhoriid, vhoriid, thoriid, zvertiid
49 integer ii, jj
50 integer zan, dayref
51
52 !-----------------------------------------------------------------------
53
54 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
55
56 zan = anne0
57 dayref = day0
58 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
59
60 do jj = 1, jjp1
61 do ii = 1, iip1
62 rlong(ii, jj) = rlonu(ii) * 180. / pi
63 rlat(ii, jj) = rlatu(jj) * 180. / pi
64 enddo
65 enddo
66
67 call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, &
68 jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid)
69
70 do jj = 1, jjm
71 do ii = 1, iip1
72 rlong(ii, jj) = rlonv(ii) * 180. / pi
73 rlat(ii, jj) = rlatv(jj) * 180. / pi
74 enddo
75 enddo
76
77 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
78 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
79 ! un meme fichier)
80
81 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), &
82 1, iip1, 1, jjm, &
83 itau_dyn, zjulian, tstep, vhoriid, histvid)
84 !
85 ! Appel a histhori pour rajouter les autres grilles horizontales
86 !
87 do jj = 1, jjp1
88 do ii = 1, iip1
89 rlong(ii, jj) = rlonv(ii) * 180. / pi
90 rlat(ii, jj) = rlatu(jj) * 180. / pi
91 enddo
92 enddo
93
94 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, &
95 jjp1, itau_dyn, zjulian, tstep, thoriid, histid)
96 !
97 ! Appel a histvert pour la grille verticale
98 !
99 call histvert(histid, 'presnivs', 'Niveaux pression','mb', llm, &
100 presnivs/100., zvertiid,'down')
101 call histvert(histvid, 'presnivs', 'Niveaux pression','mb', llm, &
102 presnivs/100., zvertiid,'down')
103 call histvert(histuid, 'presnivs', 'Niveaux pression','mb', llm, &
104 presnivs/100., zvertiid,'down')
105 !
106 ! Appels a histdef pour la definition des variables a sauvegarder
107 !
108 ! Vents U
109 !
110 call histdef(histuid, 'u', 'vent u', 'm/s', &
111 iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
112 'inst(X)', t_ops, t_wrt)
113 !
114 ! Vents V
115 !
116 call histdef(histvid, 'v', 'vent v', 'm/s', &
117 iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
118 'inst(X)', t_ops, t_wrt)
119
120 !
121 ! Temperature potentielle
122 !
123 call histdef(histid, 'teta', 'temperature potentielle', '-', &
124 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
125 'inst(X)', t_ops, t_wrt)
126 !
127 ! Geopotentiel
128 !
129 call histdef(histid, 'phi', 'geopotentiel', '-', &
130 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
131 'inst(X)', t_ops, t_wrt)
132 !
133 ! Traceurs
134 !
135 DO iq=1, nq
136 call histdef(histid, ttext(iq), ttext(iq), '-', &
137 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
138 'inst(X)', t_ops, t_wrt)
139 enddo
140 !
141 ! Masse
142 !
143 call histdef(histid, 'masse', 'masse', 'kg', &
144 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
145 'inst(X)', t_ops, t_wrt)
146 !
147 ! Pression au sol
148 !
149 call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', &
150 iip1, jjp1, thoriid, 1, 1, 1, -99, &
151 'inst(X)', t_ops, t_wrt)
152 !
153 ! Geopotentiel au sol
154 !
155 call histdef(histid, 'phis', 'geopotentiel au sol', '-', &
156 iip1, jjp1, thoriid, 1, 1, 1, -99, &
157 'inst(X)', t_ops, t_wrt)
158 !
159 ! Fin
160 !
161 call histend(histid)
162 call histend(histuid)
163 call histend(histvid)
164
165 end subroutine inithist
166
167 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21