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

Annotation of /trunk/Sources/dyn3d/initdynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 129 - (hide annotations)
Fri Feb 13 18:22:38 2015 UTC (9 years, 3 months ago) by guez
Original Path: trunk/bibio/initdynav.f
File size: 2911 byte(s)
Removed arguments day0, anne0 of procedures initdynav and
inithist. Use directly day_ref, annee_ref instead.

Moved variables annee_ref, day_ref of module temps to module
dynetat0_m. Deleted variables dayref and anneeref of module conf_gcm_m
and removed them from namelist conf_gcm_nml. These variables were
troubling intermediary on the way to annee_ref and day_ref. Gave as
default values to annee_ref and day_ref the default values of dayref
and anneeref. Moved the test on raz_date from main unit gcm to
procedure dynetat0. Created namelist dynetat0_nml. Read annee_ref and
day_ref from standard input in dynetat0 and redefine them from
"start.nc" if not raz_date. Rationale: 1 - Choose the best programming
from the point of view of program gcm only, forgetting program ce0l. 2
- The normal case is to define annee_ref and day_ref from "start.nc"
so put them in module dynetat0_m rather than in conf_gcm_m. 3 - Try to
always read the same namelists in the same order regardless of choices
in previous namelists. Downsides: 1 -We now need the file "dynetat0.f"
for the program ce0l, because dynetat0_m is used by dynredem0. 2 - We
need to define annee_ref and day_ref from procedure etat0.

Removed useless variable day_end of module temps.

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

  ViewVC Help
Powered by ViewVC 1.1.21