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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (show annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years ago) by guez
File size: 3333 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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

  ViewVC Help
Powered by ViewVC 1.1.21