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

Contents of /trunk/libf/bibio/inithist.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: 4906 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 inithist_m
2
3 ! This module is clean: no C preprocessor directive, no include line
4
5 implicit none
6
7 contains
8
9 subroutine inithist(day0, anne0, tstep, nq, fileid, filevid, t_ops, &
10 t_wrt)
11
12 ! From inithist.F,v 1.1.1.1 2004/05/19 12:53:05
13
14 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
15 ! au format IOIPSL
16
17 ! Appels succesifs des routines: histbeg
18 ! histhori
19 ! histver
20 ! histdef
21 ! histend
22
23 ! Entree:
24 ! day0,anne0: date de reference
25 ! tstep: duree du pas de temps en seconde
26 ! t_ops: frequence de l'operation pour IOIPSL
27 ! t_wrt: frequence d'ecriture sur le fichier
28 ! nq: nombre de traceurs
29
30 ! Sortie:
31 ! fileid: ID du fichier netcdf cree
32 ! filevid:ID du fichier netcdf pour la grille v
33
34 ! L. Fairhead, LMD, 03/99
35
36 USE calendar
37 use histcom
38 use dimens_m
39 use paramet_m
40 use comconst
41 use comvert
42 use logic
43 use comgeom
44 use serre
45 use temps
46 use ener
47 use iniadvtrac_m
48
49 ! Arguments
50 integer day0, anne0
51 real, intent(in):: tstep, t_ops, t_wrt
52 integer fileid, filevid
53 integer nq
54
55 ! Variables locales
56 real zjulian
57 integer iq
58 real rlong(iip1,jjp1), rlat(iip1,jjp1)
59 integer uhoriid, vhoriid, thoriid, zvertiid
60 integer ii,jj
61 integer zan, dayref
62
63 !-----------------------------------------------------------------------
64
65 ! Initialisations
66
67 pi = 4. * atan (1.)
68
69 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
70
71 zan = anne0
72 dayref = day0
73 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
74
75 do jj = 1, jjp1
76 do ii = 1, iip1
77 rlong(ii,jj) = rlonu(ii) * 180. / pi
78 rlat(ii,jj) = rlatu(jj) * 180. / pi
79 enddo
80 enddo
81
82 call histbeg_totreg("dyn_hist.nc", rlong(:,1), rlat(1,:), &
83 1, iip1, 1, jjp1, &
84 itau_dyn, zjulian, tstep, uhoriid, fileid)
85 !
86 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
87 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
88 ! un meme fichier)
89
90 do jj = 1, jjm
91 do ii = 1, iip1
92 rlong(ii,jj) = rlonv(ii) * 180. / pi
93 rlat(ii,jj) = rlatv(jj) * 180. / pi
94 enddo
95 enddo
96
97 call histbeg_totreg('dyn_histv.nc', rlong(:,1), rlat(1,:jjm), &
98 1, iip1, 1, jjm, &
99 itau_dyn, zjulian, tstep, vhoriid, filevid)
100 !
101 ! Appel a histhori pour rajouter les autres grilles horizontales
102 !
103 do jj = 1, jjp1
104 do ii = 1, iip1
105 rlong(ii,jj) = rlonv(ii) * 180. / pi
106 rlat(ii,jj) = rlatu(jj) * 180. / pi
107 enddo
108 enddo
109
110 call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &
111 'Grille points scalaires', thoriid)
112 !
113 ! Appel a histvert pour la grille verticale
114 !
115 call histvert(fileid, 'sig_s', 'Niveaux sigma','-', &
116 llm, nivsigs, zvertiid)
117 ! Pour le fichier V
118 call histvert(filevid, 'sig_s', 'Niveaux sigma','-', &
119 llm, nivsigs, zvertiid)
120 !
121 ! Appels a histdef pour la definition des variables a sauvegarder
122 !
123 ! Vents U
124 !
125 call histdef(fileid, 'ucov', 'vents u covariants', 'm/s', &
126 iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
127 'inst(X)', t_ops, t_wrt)
128 !
129 ! Vents V
130 !
131 call histdef(filevid, 'vcov', 'vents v covariants', 'm/s', &
132 iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
133 'inst(X)', t_ops, t_wrt)
134
135 !
136 ! Temperature potentielle
137 !
138 call histdef(fileid, 'teta', 'temperature potentielle', '-', &
139 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
140 'inst(X)', t_ops, t_wrt)
141 !
142 ! Geopotentiel
143 !
144 call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &
145 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
146 'inst(X)', t_ops, t_wrt)
147 !
148 ! Traceurs
149 !
150 DO iq=1,nq
151 call histdef(fileid, ttext(iq), ttext(iq), '-', &
152 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
153 'inst(X)', t_ops, t_wrt)
154 enddo
155 !
156 ! Masse
157 !
158 call histdef(fileid, 'masse', 'masse', 'kg', &
159 iip1, jjp1, thoriid, 1, 1, 1, -99, &
160 'inst(X)', t_ops, t_wrt)
161 !
162 ! Pression au sol
163 !
164 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
165 iip1, jjp1, thoriid, 1, 1, 1, -99, &
166 'inst(X)', t_ops, t_wrt)
167 !
168 ! Pression au sol
169 !
170 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
171 iip1, jjp1, thoriid, 1, 1, 1, -99, &
172 'inst(X)', t_ops, t_wrt)
173 !
174 ! Fin
175 !
176 call histend(fileid)
177 call histend(filevid)
178
179 end subroutine inithist
180
181 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21