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

Contents of /trunk/libf/bibio/initdynav.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: 3247 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 initdynav_m
2
3 implicit none
4
5 integer histaveid
6
7 contains
8
9 subroutine initdynav(day0, anne0, tstep, nq, 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 ! Routine d'initialisation des écritures des fichiers histoires LMDZ
14 ! au format IOIPSL. Initialisation du fichier histoire moyenne.
15
16 use calendar, ONLY: ymds2ju
17 USE comvert, ONLY : nivsigs
18 USE comgeom, ONLY : rlatu, rlonv
19 USE dimens_m, ONLY : llm
20 USE histcom, ONLY: histbeg_totreg, histdef, histend, histvert
21 USE iniadvtrac_m, ONLY : ttext
22 USE nr_util, ONLY : pi
23 USE paramet_m, ONLY : iip1, jjp1
24 USE temps, ONLY : itau_dyn
25
26 integer, intent(in):: day0, anne0 ! date de reference
27 real, intent(in):: tstep ! frequence d'ecriture
28 integer, intent(in):: nq ! nombre de traceurs
29 real, intent(in):: t_ops ! frequence de l'operation pour IOIPSL
30 real, intent(in):: t_wrt ! frequence d'ecriture sur le fichier
31
32 ! Variables locales
33 integer thoriid, zvertiid
34 real zjulian
35 integer iq
36 real rlong(iip1, jjp1), rlat(iip1, jjp1)
37 integer ii, jj
38 integer zan, dayref
39
40 !----------------------------------------------------
41
42 print *, "Call sequence information: initdynav"
43
44 zan = anne0
45 dayref = day0
46 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
47
48 do jj = 1, jjp1
49 do ii = 1, iip1
50 rlong(ii, jj) = rlonv(ii) * 180. / pi
51 rlat(ii, jj) = rlatu(jj) * 180. / pi
52 enddo
53 enddo
54
55 call histbeg_totreg('dyn_hist_ave.nc', rlong(:, 1), rlat(1, :), 1, iip1, &
56 1, jjp1, itau_dyn, zjulian, tstep, thoriid, histaveid)
57 call histvert(histaveid, 'sigss', 'Niveaux sigma', 'Pa', llm, nivsigs, &
58 zvertiid)
59
60 call histdef(histaveid, 'u', 'vents u scalaires moyennes', &
61 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
62 'ave(X)', t_ops, t_wrt)
63 call histdef(histaveid, 'v', 'vents v scalaires moyennes', &
64 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
65 'ave(X)', t_ops, t_wrt)
66 call histdef(histaveid, 'temp', 'temperature moyennee', 'K', &
67 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
68 'ave(X)', t_ops, t_wrt)
69 call histdef(histaveid, 'theta', 'temperature potentielle', 'K', &
70 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
71 'ave(X)', t_ops, t_wrt)
72 call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', &
73 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
74 'ave(X)', t_ops, t_wrt)
75
76 ! Traceurs
77 DO iq = 1, nq
78 call histdef(histaveid, ttext(iq), ttext(iq), '-', &
79 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
80 'ave(X)', t_ops, t_wrt)
81 enddo
82
83 call histdef(histaveid, 'masse', 'masse', 'kg', &
84 iip1, jjp1, thoriid, 1, 1, 1, -99, &
85 'ave(X)', t_ops, t_wrt)
86 call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', &
87 iip1, jjp1, thoriid, 1, 1, 1, -99, &
88 'ave(X)', t_ops, t_wrt)
89 call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', &
90 iip1, jjp1, thoriid, 1, 1, 1, -99, &
91 'ave(X)', t_ops, t_wrt)
92
93 call histend(histaveid)
94
95 end subroutine initdynav
96
97 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21