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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 253 - (hide annotations)
Tue Jan 23 15:49:10 2018 UTC (6 years, 4 months ago) by guez
File size: 2857 byte(s)
No need for intermediary variables rlong and rlat in inithist.
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 253 ! Routine d'initialisation des écritures du fichier histoire
15     ! moyenne au format IOIPSL.
16 guez 67
17 guez 62 USE dimens_m, ONLY: llm
18 guez 139 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlonv
19 guez 253 use histbeg_totreg_m, only: histbeg_totreg
20 guez 62 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 guez 92 use ymds2ju_m, ONLY: ymds2ju
28 guez 3
29 guez 67 real, intent(in):: tstep ! fréquence d'écriture
30 guez 56 integer, intent(in):: nq ! nombre de traceurs
31 guez 67 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 guez 62 real julian
36 guez 105 integer iq, l
37 guez 253 integer horiid, zvertiid
38 guez 3
39 guez 253 !-----------------------------------------------------------------------
40 guez 3
41 guez 40 print *, "Call sequence information: initdynav"
42 guez 129 CALL ymds2ju(annee_ref, 1, day_ref, 0., julian)
43 guez 62 call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, &
44     rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, &
45     horiid, histaveid)
46 guez 83 call histvert(histaveid, 'sigss', 'Niveaux sigma', '', &
47     (/(real(l), l = 1, llm)/), zvertiid)
48 guez 3
49 guez 62 call histdef(histaveid, 'u', 'vents u scalaires moyennes', 'm/s', iip1, &
50     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
51     call histdef(histaveid, 'v', 'vents v scalaires moyennes', 'm/s', iip1, &
52     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
53     call histdef(histaveid, 'temp', 'temperature moyennee', 'K', iip1, jjp1, &
54     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
55     call histdef(histaveid, 'theta', 'temperature potentielle', 'K', iip1, &
56     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
57     call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', iip1, jjp1, &
58     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
59 guez 27
60 guez 40 ! Traceurs
61     DO iq = 1, nq
62 guez 62 call histdef(histaveid, ttext(iq), ttext(iq), '-', iip1, jjp1, &
63     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
64 guez 3 enddo
65 guez 27
66 guez 62 call histdef(histaveid, 'masse', 'masse', 'kg', iip1, jjp1, horiid, 1, &
67     1, 1, -99, 'ave(X)', t_ops, t_wrt)
68     call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', iip1, &
69     jjp1, horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
70     call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
71     horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
72 guez 27
73 guez 56 call histend(histaveid)
74 guez 3
75     end subroutine initdynav
76    
77     end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21