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

Contents of /trunk/bibio/inithist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 129 - (show annotations)
Fri Feb 13 18:22:38 2015 UTC (9 years, 3 months ago) by guez
File size: 4131 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 module inithist_m
2
3 implicit none
4
5 contains
6
7 subroutine inithist(tstep, nq, t_ops, t_wrt)
8
9 ! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05
10 ! L. Fairhead, LMD, 03/99
11
12 ! Routine d'initialisation des écritures des fichiers histoires
13 ! LMDZ au format IOIPSL.
14
15 USE comgeom, ONLY: rlatu, rlatv, rlonu, rlonv
16 USE com_io_dyn, ONLY: histid, histuid, histvid
17 USE dimens_m, ONLY: jjm, llm
18 USE disvert_m, ONLY: presnivs
19 use dynetat0_m, only: day_ref, annee_ref
20 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 USE ymds2ju_m, ONLY: ymds2ju
29
30 real, intent(in):: tstep ! durée du pas de temps en secondes
31 integer, intent(in):: nq ! nombre de traceurs
32 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
35 ! Variables locales:
36 real zjulian
37 integer iq
38 real rlong(iip1, jjp1), rlat(iip1, jjp1)
39 integer uhoriid, vhoriid, thoriid, zvertiid
40 integer ii, jj
41
42 !-----------------------------------------------------------------------
43
44 CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian)
45
46 do jj = 1, jjp1
47 do ii = 1, iip1
48 rlong(ii, jj) = rlonu(ii) * 180. / pi
49 rlat(ii, jj) = rlatu(jj) * 180. / pi
50 enddo
51 enddo
52
53 call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, &
54 jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid)
55
56 do jj = 1, jjm
57 do ii = 1, iip1
58 rlong(ii, jj) = rlonv(ii) * 180. / pi
59 rlat(ii, jj) = rlatv(jj) * 180. / pi
60 enddo
61 enddo
62
63 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
64 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
65 ! un meme fichier)
66
67 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, &
68 1, jjm, itau_dyn, zjulian, tstep, vhoriid, histvid)
69
70 do jj = 1, jjp1
71 do ii = 1, iip1
72 rlong(ii, jj) = rlonv(ii) * 180. / pi
73 rlat(ii, jj) = rlatu(jj) * 180. / pi
74 enddo
75 enddo
76
77 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, &
78 jjp1, itau_dyn, zjulian, tstep, thoriid, histid)
79
80 ! Appel a histvert pour la grille verticale
81
82 call histvert(histid, 'presnivs', 'Niveaux pression','mb', presnivs/100., &
83 zvertiid,'down')
84 call histvert(histvid, 'presnivs', 'Niveaux pression','mb', &
85 presnivs/100., zvertiid,'down')
86 call histvert(histuid, 'presnivs', 'Niveaux pression','mb', &
87 presnivs/100., zvertiid,'down')
88
89 ! Appels a histdef pour la definition des variables a sauvegarder
90
91 call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, &
92 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
93 call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, &
94 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
95 call histdef(histid, 'teta', 'temperature potentielle', '-', iip1, jjp1, &
96 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
97 call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, &
98 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
99
100 ! Traceurs
101 DO iq=1, nq
102 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
103 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
104 enddo
105
106 call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
107 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
108 call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, &
109 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
110 call histdef(histid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
111 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
112
113 call histend(histid)
114 call histend(histuid)
115 call histend(histvid)
116
117 end subroutine inithist
118
119 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21