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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (show annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 7 months ago) by guez
File size: 3260 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

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

  ViewVC Help
Powered by ViewVC 1.1.21