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

Contents of /trunk/Sources/bibio/initfluxsto.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 139 - (show annotations)
Tue May 26 17:46:03 2015 UTC (8 years, 11 months ago) by guez
File size: 6211 byte(s)
dynetat0 read rlonu, rlatu, rlonv, rlatv, cu_2d, cv_2d, aire_2d from
"start.nc" and then these variables were overwritten by
inigeom. Corrected this. Now, inigeom does not compute rlonu, rlatu,
rlonv and rlatv. Moreover, cu_2d, cv_2d, aire_2d are not written to
"restart.nc". Since xprimu, xprimv, xprimm025, xprimp025, rlatu1,
rlatu2, yprimu1, yprimu2 are computed at the same time as rlonu,
rlatu, rlonv, rlatv, and since it would not be convenient to separate
those computations, we decide to write xprimu, xprimv, xprimm025,
xprimp025, rlatu1, rlatu2, yprimu1, yprimu2 into "restart.nc", read
them from "start.nc" and not compute them in inigeom. So, in summary,
"start.nc" contains all the coordinates and their derivatives, and
inigeom only computes the 2D-variables.

Technical details:

Moved variables rlatu, rlonv, rlonu, rlatv, xprimu, xprimv from module
comgeom to module dynetat0_m. Upgraded local variables rlatu1,
yprimu1, rlatu2, yprimu2, xprimm025, xprimp025 of procedure inigeom to
variables of module dynetat0_m.

Removed unused local variable yprimu of procedure inigeom and
corresponding argument yyprimu of fyhyp.

Moved variables clat, clon, grossismx, grossismy, dzoomx, dzoomy,
taux, tauy from module serre to module dynetat0_m (since they are read
from "start.nc"). The default values are now defined in read_serre
instead of in the declarations. Changed name of module serre to
read_serre_m, no more module variable here.

The calls to fxhyp and fyhyp are moved from inigeom to etat0.

Side effects in programs other than gcm: etat0 and read_serre write
variables of module dynetat0; the programs test_fxyp and
test_inter_barxy need more source files.

Removed unused arguments len and nd of cv3_tracer. Removed unused
argument PPSOL of LWU.

Bug fix in test_inter_barxy: forgotten call to read_serre.

1 SUBROUTINE initfluxsto(tstep, t_ops, t_wrt, nq, fileid, filevid, filedid)
2
3 ! From bibio/initfluxsto.F, v 1.1.1.1 2004/05/19 12:53:05
4
5 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
6 ! au format IOIPSL
7 ! Appels succesifs des routines: histbeg
8 ! histhori
9 ! histver
10 ! histdef
11 ! histend
12
13 ! Entree:
14 ! day0, anne0: date de reference
15 ! tstep: duree du pas de temps en seconde
16 ! t_ops: frequence de l'operation pour IOIPSL
17 ! t_wrt: frequence d'ecriture sur le fichier
18 ! nq: nombre de traceurs
19
20 ! Sortie:
21 ! fileid: ID du fichier netcdf cree
22 ! filevid:ID du fichier netcdf pour la grille v
23
24 ! L. Fairhead, 03/99
25
26 USE comconst
27 use conf_gcm_m
28 USE dimens_m
29 USE disvert_m
30 use dynetat0_m, only: day_ref, annee_ref, rlonu, rlatu, rlonv, rlatv
31 USE histbeg_totreg_m, ONLY : histbeg_totreg
32 USE histdef_m, ONLY : histdef
33 USE histend_m, ONLY : histend
34 use histhori_regular_m, only: histhori_regular
35 use histsync_m, only: histsync
36 USE histvert_m, ONLY : histvert
37 USE nr_util, ONLY : pi
38 USE paramet_m
39 USE temps, ONLY : itau_dyn
40 use ymds2ju_m, only: ymds2ju
41
42 IMPLICIT NONE
43
44 ! Arguments
45 INTEGER itau
46 REAL, INTENT (IN) :: tstep
47 REAL t_ops, t_wrt
48 INTEGER fileid, filevid, filedid
49 INTEGER nq, ndex(1)
50 REAL nivd(1)
51
52 ! Variables locales
53 REAL zjulian
54 CHARACTER(len=3) str
55 CHARACTER(len=10) ctrac
56 INTEGER iq
57 REAL rlong(iip1, jjp1), rlat(iip1, jjp1)
58 INTEGER uhoriid, vhoriid, thoriid, zvertiid, dhoriid, dvertiid
59 INTEGER ii, jj, l
60 LOGICAL ok_sync
61
62 !---------------------------------------------------------
63
64 ! Initialisations
65 str = 'q '
66 ctrac = 'traceur '
67 ok_sync = .TRUE.
68
69 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
70
71 CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian)
72
73 DO jj = 1, jjp1
74 DO ii = 1, iip1
75 rlong(ii, jj) = rlonu(ii)*180./pi
76 rlat(ii, jj) = rlatu(jj)*180./pi
77 END DO
78 END DO
79
80 CALL histbeg_totreg('fluxstoke', rlong(:, 1), rlat(1, :), 1, iip1, 1, jjp1, &
81 itau_dyn, zjulian, tstep, uhoriid, fileid)
82
83 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
84 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
85 ! un meme fichier)
86
87 DO jj = 1, jjm
88 DO ii = 1, iip1
89 rlong(ii, jj) = rlonv(ii)*180./pi
90 rlat(ii, jj) = rlatv(jj)*180./pi
91 END DO
92 END DO
93
94 CALL histbeg_totreg('fluxstokev.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, &
95 1, jjm, itau_dyn, zjulian, tstep, vhoriid, filevid)
96
97 CALL histbeg_totreg('defstoke.nc', (/1./), (/1./), 1, 1, 1, 1, itau_dyn, &
98 zjulian, tstep, dhoriid, filedid)
99
100 ! Appel a histhori pour rajouter les autres grilles horizontales
101
102 DO jj = 1, jjp1
103 DO ii = 1, iip1
104 rlong(ii, jj) = rlonv(ii)*180./pi
105 rlat(ii, jj) = rlatu(jj)*180./pi
106 END DO
107 END DO
108
109 CALL histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &
110 'Grille points scalaires', thoriid)
111
112 ! Appel a histvert pour la grille verticale
113
114 CALL histvert(fileid, 'sig_s', 'Niveaux sigma', 'sigma_level', &
115 (/(real(l), l = 1, llm)/), zvertiid)
116 ! Pour le fichier V
117 CALL histvert(filevid, 'sig_s', 'Niveaux sigma', 'sigma_level', &
118 (/(real(l), l = 1, llm)/), zvertiid)
119 ! pour le fichier def
120 nivd(1) = 1
121 CALL histvert(filedid, 'sig_s', 'Niveaux sigma', 'sigma_level', nivd, &
122 dvertiid)
123
124 ! Appels a histdef pour la definition des variables a sauvegarder
125 CALL histdef(fileid, 'phis', 'Surface geop. height', '-', iip1, jjp1, &
126 thoriid, 1, 1, 1, -99, 'once', t_ops, t_wrt)
127 CALL histdef(fileid, 'aire', 'Grid area', '-', iip1, jjp1, thoriid, 1, 1, &
128 1, -99, 'once', t_ops, t_wrt)
129 CALL histdef(filedid, 'dtvr', 'tps dyn', 's', 1, 1, dhoriid, 1, 1, 1, -99, &
130 'once', t_ops, t_wrt)
131 CALL histdef(filedid, 'istdyn', 'tps stock', 's', 1, 1, dhoriid, 1, 1, 1, &
132 -99, 'once', t_ops, t_wrt)
133 CALL histdef(filedid, 'istphy', 'tps stock phy', 's', 1, 1, dhoriid, 1, 1, &
134 1, -99, 'once', t_ops, t_wrt)
135 CALL histdef(fileid, 'masse', 'Masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
136 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
137 CALL histdef(fileid, 'pbaru', 'flx de masse zonal', 'kg m/s', iip1, jjp1, &
138 uhoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
139 CALL histdef(filevid, 'pbarv', 'flx de masse mer', 'kg m/s', iip1, jjm, &
140 vhoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
141 CALL histdef(fileid, 'w', 'flx de masse vert', 'kg m/s', iip1, jjp1, &
142 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
143 CALL histdef(fileid, 'teta', 'temperature potentielle', '-', iip1, jjp1, &
144 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
145 CALL histdef(fileid, 'phi', 'geopotentiel instantane', '-', iip1, jjp1, &
146 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
147
148 CALL histend(fileid)
149 CALL histend(filevid)
150 CALL histend(filedid)
151 IF (ok_sync) THEN
152 CALL histsync(fileid)
153 CALL histsync(filevid)
154 CALL histsync(filedid)
155 END IF
156
157 END SUBROUTINE initfluxsto

  ViewVC Help
Powered by ViewVC 1.1.21