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

Contents of /trunk/dyn3d/initdynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 40 - (show annotations)
Tue Feb 22 13:49:36 2011 UTC (13 years, 3 months ago) by guez
Original Path: trunk/libf/bibio/initdynav.f90
File size: 3332 byte(s)
"alpha" useless, always 0, in "exner_hyb".

1 module initdynav_m
2
3 ! This module is clean: no C preprocessor directive, no include line
4
5 implicit none
6
7 contains
8
9 subroutine initdynav(day0, anne0, tstep, nq, fileid, 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 real, intent(in):: t_ops ! frequence de l'operation pour IOIPSL
29 real, intent(in):: t_wrt ! frequence d'ecriture sur le fichier
30 integer, intent(out):: fileid ! ID du fichier netcdf cree
31 integer, intent(in):: nq ! nombre de traceurs
32
33 ! Variables locales
34 integer thoriid, zvertiid
35 real zjulian
36 integer iq
37 real rlong(iip1, jjp1), rlat(iip1, jjp1)
38 integer ii, jj
39 integer zan, dayref
40
41 !----------------------------------------------------
42
43 print *, "Call sequence information: initdynav"
44
45 zan = anne0
46 dayref = day0
47 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
48
49 do jj = 1, jjp1
50 do ii = 1, iip1
51 rlong(ii, jj) = rlonv(ii) * 180. / pi
52 rlat(ii, jj) = rlatu(jj) * 180. / pi
53 enddo
54 enddo
55
56 call histbeg_totreg('dyn_hist_ave.nc', rlong(:, 1), rlat(1, :), 1, iip1, &
57 1, jjp1, itau_dyn, zjulian, tstep, thoriid, fileid)
58 call histvert(fileid, 'sigss', 'Niveaux sigma', 'Pa', llm, nivsigs, &
59 zvertiid)
60
61 call histdef(fileid, 'u', 'vents u scalaires moyennes', &
62 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
63 'ave(X)', t_ops, t_wrt)
64 call histdef(fileid, 'v', 'vents v scalaires moyennes', &
65 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
66 'ave(X)', t_ops, t_wrt)
67 call histdef(fileid, 'temp', 'temperature moyennee', 'K', &
68 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
69 'ave(X)', t_ops, t_wrt)
70 call histdef(fileid, 'theta', 'temperature potentielle', 'K', &
71 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
72 'ave(X)', t_ops, t_wrt)
73 call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', &
74 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
75 'ave(X)', t_ops, t_wrt)
76
77 ! Traceurs
78 DO iq = 1, nq
79 call histdef(fileid, ttext(iq), ttext(iq), '-', &
80 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
81 'ave(X)', t_ops, t_wrt)
82 enddo
83
84 call histdef(fileid, 'masse', 'masse', 'kg', &
85 iip1, jjp1, thoriid, 1, 1, 1, -99, &
86 'ave(X)', t_ops, t_wrt)
87 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
88 iip1, jjp1, thoriid, 1, 1, 1, -99, &
89 'ave(X)', t_ops, t_wrt)
90 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
91 iip1, jjp1, thoriid, 1, 1, 1, -99, &
92 'ave(X)', t_ops, t_wrt)
93
94 call histend(fileid)
95
96 end subroutine initdynav
97
98 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21