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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 30 - (show annotations)
Thu Apr 1 09:07:28 2010 UTC (14 years, 1 month ago) by guez
File size: 4009 byte(s)
Imported Source files of the external library "IOIPSL_Lionel" into
"libf/IOIPSL".

Split "cray.f90" into "scopy.f90" and "ssum.f90".

Rewrote "leapfrog" in order to have a clearer algorithmic structure.

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
14 USE histcom, ONLY: histbeg_totreg, histdef, histend, histvert
15 use calendar, ONLY: ymds2ju
16 USE dimens_m, ONLY : llm
17 USE paramet_m, ONLY : iip1, jjp1
18 USE comconst, ONLY : pi
19 USE comvert, ONLY : nivsigs
20 USE comgeom, ONLY : rlatu, rlonv
21 USE temps, ONLY : itau_dyn
22 USE iniadvtrac_m, ONLY : ttext
23
24 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
25 ! au format IOIPSL. Initialisation du fichier histoire moyenne.
26
27 ! Appels succesifs des routines: histbeg
28 ! histhori
29 ! histver
30 ! histdef
31 ! histend
32
33 ! Entree:
34 ! day0,anne0: date de reference
35 ! tstep : frequence d'ecriture
36 ! t_ops: frequence de l'operation pour IOIPSL
37 ! t_wrt: frequence d'ecriture sur le fichier
38 ! nq: nombre de traceurs
39
40 ! Sortie:
41 ! fileid: ID du fichier netcdf cree
42
43 ! Arguments
44 integer day0, anne0
45 real, intent(in):: tstep, t_ops, t_wrt
46 integer fileid
47 integer nq
48 integer thoriid, zvertiid
49
50 ! Variables locales
51
52 real zjulian
53 integer iq
54 real rlong(iip1,jjp1), rlat(iip1,jjp1)
55 integer ii,jj
56 integer zan, dayref
57
58 !----------------------------------------------------
59
60 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
61
62 zan = anne0
63 dayref = day0
64 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
65
66 do jj = 1, jjp1
67 do ii = 1, iip1
68 rlong(ii,jj) = rlonv(ii) * 180. / pi
69 rlat(ii,jj) = rlatu(jj) * 180. / pi
70 enddo
71 enddo
72
73 call histbeg_totreg('dyn_hist_ave.nc', rlong(:,1), rlat(1,:), &
74 1, iip1, 1, jjp1, &
75 itau_dyn, zjulian, tstep, thoriid, fileid)
76
77
78 ! Appel a histvert pour la grille verticale
79
80 call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', &
81 llm, nivsigs, zvertiid)
82
83 ! Appels a histdef pour la definition des variables a sauvegarder
84
85 ! Vents U
86
87 write(6,*)'inithistave',tstep
88 call histdef(fileid, 'u', 'vents u scalaires moyennes', &
89 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
90 'ave(X)', t_ops, t_wrt)
91
92
93 ! Vents V
94
95 call histdef(fileid, 'v', 'vents v scalaires moyennes', &
96 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
97 'ave(X)', t_ops, t_wrt)
98
99
100 ! Temperature
101
102 call histdef(fileid, 'temp', 'temperature moyennee', 'K', &
103 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
104 'ave(X)', t_ops, t_wrt)
105
106 ! Temperature potentielle
107
108 call histdef(fileid, 'theta', 'temperature potentielle', 'K', &
109 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
110 'ave(X)', t_ops, t_wrt)
111
112
113
114 ! Geopotentiel
115
116 call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', &
117 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
118 'ave(X)', t_ops, t_wrt)
119
120 ! Traceurs
121
122 DO iq=1,nq
123 call histdef(fileid, ttext(iq), ttext(iq), '-', &
124 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
125 'ave(X)', t_ops, t_wrt)
126 enddo
127
128 ! Masse
129
130 call histdef(fileid, 'masse', 'masse', 'kg', &
131 iip1, jjp1, thoriid, 1, 1, 1, -99, &
132 'ave(X)', t_ops, t_wrt)
133
134 ! Pression au sol
135
136 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
137 iip1, jjp1, thoriid, 1, 1, 1, -99, &
138 'ave(X)', t_ops, t_wrt)
139
140 ! Pression au sol
141
142 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
143 iip1, jjp1, thoriid, 1, 1, 1, -99, &
144 'ave(X)', t_ops, t_wrt)
145
146 call histend(fileid)
147
148 end subroutine initdynav
149
150 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21