/[lmdze]/trunk/bibio/initfluxsto.f
ViewVC logotype

Annotation of /trunk/bibio/initfluxsto.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 15 - (hide annotations)
Fri Aug 1 15:24:12 2008 UTC (15 years, 9 months ago) by guez
Original Path: trunk/libf/bibio/initfluxsto.f
File size: 6115 byte(s)
-- Minor modification of input/output:

Added variable "Sigma_O3_Royer" to "histday.nc". "ecrit_day" is not
modified in "physiq". Removed variables "pyu1", "pyv1", "ftsol1",
"ftsol2", "ftsol3", "ftsol4", "psrf1", "psrf2", "psrf3", "psrf4"
"mfu", "mfd", "en_u", "en_d", "de_d", "de_u", "coefh" from
"histrac.nc".

Variable "raz_date" of module "conf_gcm_m" has logical type instead of
integer type.

-- Should not change any result at run time:

Modified calls to "IOIPSL_Lionel" procedures because the interfaces of
these procedures have been simplified.

Changed name of variable in module "start_init_orog_m": "masque" to
"mask".

Created a module containing procedure "phyredem".

Removed arguments "punjours", "pdayref" and "ptimestep" of procedure
"iniphysiq".

Renamed procedure "gr_phy_write" to "gr_phy_write_2d". Created
procedure "gr_phy_write_3d".

Removed procedures "ini_undefstd", "moy_undefSTD", "calcul_STDlev",
"calcul_divers".

1 guez 3 !
2     ! $Header: /home/cvsroot/LMDZ4/libf/bibio/initfluxsto.F,v 1.1.1.1 2004/05/19 12:53:05 lmdzadmin Exp $
3     !
4     subroutine initfluxsto
5     . (infile,tstep,t_ops,t_wrt,nq,
6     . fileid,filevid,filedid)
7    
8     USE IOIPSL
9    
10     C
11     C Routine d'initialisation des ecritures des fichiers histoires LMDZ
12     C au format IOIPSL
13     C
14     C Appels succesifs des routines: histbeg
15     C histhori
16     C histver
17     C histdef
18     C histend
19     C
20     C Entree:
21     C
22     C infile: nom du fichier histoire a creer
23     C day0,anne0: date de reference
24     C tstep: duree du pas de temps en seconde
25     C t_ops: frequence de l'operation pour IOIPSL
26     C t_wrt: frequence d'ecriture sur le fichier
27     C nq: nombre de traceurs
28     C
29     C Sortie:
30     C fileid: ID du fichier netcdf cree
31     C filevid:ID du fichier netcdf pour la grille v
32     C
33     C L. Fairhead, LMD, 03/99
34     C
35     C =====================================================================
36     C
37     C Declarations
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, only: annee_ref, day_ref, itau_dyn
46     use ener
47     implicit none
48    
49     C Arguments
50     C
51     character*(*) infile
52     integer*4 itau
53     real tstep, t_ops, t_wrt
54     integer fileid, filevid,filedid
55     integer nq,ndex(1)
56     real nivd(1)
57    
58     C Variables locales
59     C
60     integer tau0
61     real zjulian
62     character*3 str
63     character*10 ctrac
64     integer iq
65 guez 15 real rlong(iip1,jjp1), rlat(iip1,jjp1)
66 guez 3 integer uhoriid, vhoriid, thoriid, zvertiid,dhoriid,dvertiid
67     integer ii,jj
68     integer zan, idayref
69     logical ok_sync
70     C
71     C Initialisations
72     C
73     pi = 4. * atan (1.)
74     str='q '
75     ctrac = 'traceur '
76     ok_sync = .true.
77     C
78     C Appel a histbeg: creation du fichier netcdf et initialisations diverses
79     C
80    
81     zan = annee_ref
82     idayref = day_ref
83     CALL ymds2ju(zan, 1, idayref, 0.0, zjulian)
84     tau0 = itau_dyn
85    
86     do jj = 1, jjp1
87     do ii = 1, iip1
88     rlong(ii,jj) = rlonu(ii) * 180. / pi
89     rlat(ii,jj) = rlatu(jj) * 180. / pi
90     enddo
91     enddo
92    
93 guez 15 call histbeg_totreg(infile, rlong(:,1), rlat(1,:),
94 guez 3 . 1, iip1, 1, jjp1,
95     . tau0, zjulian, tstep, uhoriid, fileid)
96     C
97     C Creation du fichier histoire pour la grille en V (oblige pour l'instant,
98     C IOIPSL ne permet pas de grilles avec des nombres de point differents dans
99     C un meme fichier)
100    
101    
102     do jj = 1, jjm
103     do ii = 1, iip1
104     rlong(ii,jj) = rlonv(ii) * 180. / pi
105     rlat(ii,jj) = rlatv(jj) * 180. / pi
106     enddo
107     enddo
108    
109 guez 15 call histbeg_totreg('fluxstokev.nc', rlong(:,1),
110     . rlat(1,:jjm),1, iip1, 1, jjm,
111 guez 3 . tau0, zjulian, tstep, vhoriid, filevid)
112    
113 guez 15 call histbeg_totreg('defstoke.nc', (/1./), (/1./),
114 guez 3 . 1, 1, 1, 1,
115     . tau0, zjulian, tstep, dhoriid, filedid)
116    
117     C
118     C Appel a histhori pour rajouter les autres grilles horizontales
119     C
120     do jj = 1, jjp1
121     do ii = 1, iip1
122     rlong(ii,jj) = rlonv(ii) * 180. / pi
123     rlat(ii,jj) = rlatu(jj) * 180. / pi
124     enddo
125     enddo
126    
127 guez 15 call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar',
128 guez 3 . 'Grille points scalaires', thoriid)
129    
130     C
131     C Appel a histvert pour la grille verticale
132     C
133     call histvert(fileid, 'sig_s', 'Niveaux sigma',
134     . 'sigma_level',
135     . llm, nivsigs, zvertiid)
136     C Pour le fichier V
137     call histvert(filevid, 'sig_s', 'Niveaux sigma',
138     . 'sigma_level',
139     . llm, nivsigs, zvertiid)
140     c pour le fichier def
141     nivd(1) = 1
142     call histvert(filedid, 'sig_s', 'Niveaux sigma',
143     . 'sigma_level',
144     . 1, nivd, dvertiid)
145    
146     C
147     C Appels a histdef pour la definition des variables a sauvegarder
148    
149     CALL histdef(fileid, "phis", "Surface geop. height", "-",
150 guez 15 . iip1,jjp1,thoriid, 1,1,1, -99,
151 guez 3 . "once", t_ops, t_wrt)
152    
153     CALL histdef(fileid, "aire", "Grid area", "-",
154 guez 15 . iip1,jjp1,thoriid, 1,1,1, -99,
155 guez 3 . "once", t_ops, t_wrt)
156    
157     CALL histdef(filedid, "dtvr", "tps dyn", "s",
158 guez 15 . 1,1,dhoriid, 1,1,1, -99,
159 guez 3 . "once", t_ops, t_wrt)
160    
161     CALL histdef(filedid, "istdyn", "tps stock", "s",
162 guez 15 . 1,1,dhoriid, 1,1,1, -99,
163 guez 3 . "once", t_ops, t_wrt)
164    
165     CALL histdef(filedid, "istphy", "tps stock phy", "s",
166 guez 15 . 1,1,dhoriid, 1,1,1, -99,
167 guez 3 . "once", t_ops, t_wrt)
168    
169    
170     C
171     C Masse
172     C
173     call histdef(fileid, 'masse', 'Masse', 'kg',
174     . iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
175 guez 15 . 'inst(X)', t_ops, t_wrt)
176 guez 3 C
177     C Pbaru
178     C
179     call histdef(fileid, 'pbaru', 'flx de masse zonal', 'kg m/s',
180     . iip1, jjp1, uhoriid, llm, 1, llm, zvertiid,
181 guez 15 . 'inst(X)', t_ops, t_wrt)
182 guez 3
183     C
184     C Pbarv
185     C
186     call histdef(filevid, 'pbarv', 'flx de masse mer', 'kg m/s',
187     . iip1, jjm, vhoriid, llm, 1, llm, zvertiid,
188 guez 15 . 'inst(X)', t_ops, t_wrt)
189 guez 3 C
190     C w
191     C
192     call histdef(fileid, 'w', 'flx de masse vert', 'kg m/s',
193     . iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
194 guez 15 . 'inst(X)', t_ops, t_wrt)
195 guez 3
196     C
197     C Temperature potentielle
198     C
199     call histdef(fileid, 'teta', 'temperature potentielle', '-',
200     . iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
201 guez 15 . 'inst(X)', t_ops, t_wrt)
202 guez 3 C
203    
204     C
205     C Geopotentiel
206     C
207     call histdef(fileid, 'phi', 'geopotentiel instantane', '-',
208     . iip1, jjp1, thoriid, llm, 1, llm, zvertiid,
209 guez 15 . 'inst(X)', t_ops, t_wrt)
210 guez 3 C
211     C Fin
212     C
213     call histend(fileid)
214     call histend(filevid)
215     call histend(filedid)
216     if (ok_sync) then
217     call histsync(fileid)
218     call histsync(filevid)
219     call histsync(filedid)
220     endif
221    
222     return
223     end

  ViewVC Help
Powered by ViewVC 1.1.21