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

Contents of /trunk/dyn3d/inithist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 254 - (show annotations)
Mon Feb 5 10:39:38 2018 UTC (6 years, 3 months ago) by guez
File size: 3496 byte(s)
Move Sources/* to root directory.
1 module inithist_m
2
3 implicit none
4
5 integer histid, histvid, histuid
6
7 contains
8
9 subroutine inithist(tstep, nq, t_ops, t_wrt)
10
11 ! From inithist.F, version 1.1.1.1, 2004/05/19 12:53:05
12 ! L. Fairhead, LMD, 03/99
13
14 ! Routine d'initialisation des écritures des fichiers histoires au
15 ! format IOIPSL.
16
17 USE dimens_m, ONLY: jjm, llm
18 USE disvert_m, ONLY: presnivs
19 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlatv, rlonu, rlonv
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 ! Local:
36 real julian
37 integer iq
38 integer uhoriid, vhoriid, thoriid, zvertiid
39
40 !-----------------------------------------------------------------------
41
42 print *, "Call sequence information: inithist"
43 CALL ymds2ju(annee_ref, 1, day_ref, 0., julian)
44
45 call histbeg_totreg("dyn_histu.nc", rlonu * 180. / pi, rlatu * 180. / pi, &
46 1, iip1, 1, jjp1, itau_dyn, julian, tstep, uhoriid, histuid)
47
48 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
49 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
50 ! un meme fichier)
51
52 call histbeg_totreg('dyn_histv.nc', rlonv * 180. / pi, rlatv * 180. / pi, &
53 1, iip1, 1, jjm, itau_dyn, julian, tstep, vhoriid, histvid)
54
55 call histbeg_totreg("dyn_hist.nc", rlonv * 180. / pi, rlatu * 180. / pi, &
56 1, iip1, 1, jjp1, itau_dyn, julian, tstep, thoriid, histid)
57
58 call histvert(histid, 'presnivs', 'Niveaux pression', 'mb', presnivs/100., &
59 zvertiid, 'down')
60 call histvert(histvid, 'presnivs', 'Niveaux pression', 'mb', &
61 presnivs/100., zvertiid, 'down')
62 call histvert(histuid, 'presnivs', 'Niveaux pression', 'mb', &
63 presnivs/100., zvertiid, 'down')
64
65 call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, &
66 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
67 call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, &
68 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
69 call histdef(histid, 'teta', 'temperature potentielle', '-', iip1, jjp1, &
70 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
71 call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, &
72 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
73
74 ! Traceurs
75 DO iq = 1, nq
76 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
77 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
78 enddo
79
80 call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
81 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
82 call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, &
83 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
84 call histdef(histid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
85 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
86
87 call histend(histid)
88 call histend(histuid)
89 call histend(histvid)
90
91 end subroutine inithist
92
93 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21