/[lmdze]/trunk/Sources/dyn3d/initdynav.f
ViewVC logotype

Annotation of /trunk/Sources/dyn3d/initdynav.f

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21