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

Contents of /trunk/Sources/bibio/initdynav.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: 2889 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 module initdynav_m
2
3 implicit none
4
5 integer histaveid
6
7 contains
8
9 subroutine initdynav(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
14 ! Routine d'initialisation des écritures des fichiers histoires au
15 ! format IOIPSL. Initialisation du fichier histoire moyenne.
16
17 USE dimens_m, ONLY: llm
18 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlonv
19 USE histbeg_totreg_m, ONLY: histbeg_totreg
20 USE histdef_m, ONLY: histdef
21 USE histend_m, ONLY: histend
22 USE histvert_m, ONLY: histvert
23 USE iniadvtrac_m, ONLY: ttext
24 USE nr_util, ONLY: pi
25 USE paramet_m, ONLY: iip1, jjp1
26 USE temps, ONLY: itau_dyn
27 use ymds2ju_m, ONLY: ymds2ju
28
29 real, intent(in):: tstep ! fréquence d'écriture
30 integer, intent(in):: nq ! nombre de traceurs
31 real, intent(in):: t_ops ! fréquence de l'opération pour IOIPSL
32 real, intent(in):: t_wrt ! fréquence d'écriture sur le fichier
33
34 ! Variables locales
35 integer horiid, zvertiid
36 real julian
37 integer iq, l
38
39 !----------------------------------------------------
40
41 print *, "Call sequence information: initdynav"
42
43 CALL ymds2ju(annee_ref, 1, day_ref, 0., julian)
44 call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, &
45 rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, &
46 horiid, histaveid)
47 call histvert(histaveid, 'sigss', 'Niveaux sigma', '', &
48 (/(real(l), l = 1, llm)/), zvertiid)
49
50 call histdef(histaveid, 'u', 'vents u scalaires moyennes', 'm/s', iip1, &
51 jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
52 call histdef(histaveid, 'v', 'vents v scalaires moyennes', 'm/s', iip1, &
53 jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
54 call histdef(histaveid, 'temp', 'temperature moyennee', 'K', iip1, jjp1, &
55 horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
56 call histdef(histaveid, 'theta', 'temperature potentielle', 'K', iip1, &
57 jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
58 call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', iip1, jjp1, &
59 horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
60
61 ! Traceurs
62 DO iq = 1, nq
63 call histdef(histaveid, ttext(iq), ttext(iq), '-', iip1, jjp1, &
64 horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
65 enddo
66
67 call histdef(histaveid, 'masse', 'masse', 'kg', iip1, jjp1, horiid, 1, &
68 1, 1, -99, 'ave(X)', t_ops, t_wrt)
69 call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', iip1, &
70 jjp1, horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
71 call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
72 horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
73
74 call histend(histaveid)
75
76 end subroutine initdynav
77
78 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21