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

Annotation of /trunk/dyn3d/inithist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 302 - (hide annotations)
Thu Sep 6 13:19:51 2018 UTC (5 years, 8 months ago) by guez
File size: 3412 byte(s)
In procedure physiq, rename dsens and devap to dflux_t and dflux_q so
they correspond to arguments with the same name in pbl_surface. (In
LMDZ, dsens and devap are not used in physiq, the computation of fder
is done inside pbl_surface.)

1 guez 3 module inithist_m
2    
3     implicit none
4    
5 guez 212 integer histid, histvid, histuid
6    
7 guez 3 contains
8    
9 guez 263 subroutine inithist(t_ops, t_wrt)
10 guez 3
11 guez 253 ! From inithist.F, version 1.1.1.1, 2004/05/19 12:53:05
12 guez 61 ! L. Fairhead, LMD, 03/99
13 guez 3
14 guez 253 ! Routine d'initialisation des écritures des fichiers histoires au
15     ! format IOIPSL.
16 guez 3
17 guez 263 use comconst, only: dtvr
18 guez 265 USE dimensions, ONLY: jjm, llm, nqmx
19 guez 92 USE disvert_m, ONLY: presnivs
20 guez 139 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlatv, rlonu, rlonv
21 guez 253 use histbeg_totreg_m, only: histbeg_totreg
22     USE histdef_m, ONLY: histdef
23     USE histend_m, ONLY: histend
24     USE histvert_m, ONLY: histvert
25 guez 92 USE iniadvtrac_m, ONLY: ttext
26     USE nr_util, ONLY: pi
27 guez 61 USE paramet_m, ONLY: iip1, jjp1
28     USE temps, ONLY: itau_dyn
29 guez 253 use ymds2ju_m, ONLY: ymds2ju
30 guez 3
31 guez 61 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 guez 3
34 guez 253 ! Local:
35     real julian
36 guez 3 integer iq
37     integer uhoriid, vhoriid, thoriid, zvertiid
38    
39     !-----------------------------------------------------------------------
40    
41 guez 253 print *, "Call sequence information: inithist"
42     CALL ymds2ju(annee_ref, 1, day_ref, 0., julian)
43 guez 3
44 guez 253 call histbeg_totreg("dyn_histu.nc", rlonu * 180. / pi, rlatu * 180. / pi, &
45 guez 263 1, iip1, 1, jjp1, itau_dyn, julian, dtvr, uhoriid, histuid)
46 guez 3
47 guez 56 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
48     ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
49     ! un meme fichier)
50    
51 guez 253 call histbeg_totreg('dyn_histv.nc', rlonv * 180. / pi, rlatv * 180. / pi, &
52 guez 263 1, iip1, 1, jjm, itau_dyn, julian, dtvr, vhoriid, histvid)
53 guez 61
54 guez 253 call histbeg_totreg("dyn_hist.nc", rlonv * 180. / pi, rlatu * 180. / pi, &
55 guez 263 1, iip1, 1, jjp1, itau_dyn, julian, dtvr, thoriid, histid)
56 guez 3
57 guez 253 call histvert(histid, 'presnivs', 'Niveaux pression', 'mb', presnivs/100., &
58     zvertiid, 'down')
59     call histvert(histvid, 'presnivs', 'Niveaux pression', 'mb', &
60     presnivs/100., zvertiid, 'down')
61     call histvert(histuid, 'presnivs', 'Niveaux pression', 'mb', &
62     presnivs/100., zvertiid, 'down')
63 guez 61
64     call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, &
65     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
66     call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, &
67     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
68 guez 261 call histdef(histid, 'temp', 'temperature', 'K', iip1, jjp1, &
69 guez 61 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
70 guez 261 call histdef(histid, 'theta', 'temperature potentielle', 'K', iip1, jjp1, &
71     thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
72 guez 302 call histdef(histid, 'phi', 'geopotential', 'm2 s-2', iip1, jjp1, thoriid, &
73 guez 61 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
74    
75 guez 39 ! Traceurs
76 guez 263 DO iq = 1, nqmx
77 guez 61 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
78     llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
79 guez 3 enddo
80 guez 61
81     call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
82     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
83     call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, &
84     thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
85    
86 guez 56 call histend(histid)
87     call histend(histuid)
88     call histend(histvid)
89 guez 3
90     end subroutine inithist
91    
92     end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21