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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 76 - (hide annotations)
Fri Nov 15 18:45:49 2013 UTC (10 years, 6 months ago) by guez
Original Path: trunk/bibio/initdynav.f90
File size: 2942 byte(s)
Moved everything out of libf.
1 guez 3 module initdynav_m
2    
3     implicit none
4    
5 guez 56 integer histaveid
6    
7 guez 3 contains
8    
9 guez 56 subroutine initdynav(day0, anne0, 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 30 use calendar, ONLY: ymds2ju
18 guez 62 USE comgeom, ONLY: rlatu, rlonv
19     USE dimens_m, ONLY: llm
20 guez 67 USE disvert_m, ONLY: nivsigs
21 guez 62 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     USE iniadvtrac_m, ONLY: ttext
26     USE nr_util, ONLY: pi
27     USE paramet_m, ONLY: iip1, jjp1
28     USE temps, ONLY: itau_dyn
29 guez 3
30 guez 67 integer, intent(in):: day0, anne0 ! date de référence
31     real, intent(in):: tstep ! fréquence d'écriture
32 guez 56 integer, intent(in):: nq ! nombre de traceurs
33 guez 67 real, intent(in):: t_ops ! fréquence de l'opération pour IOIPSL
34     real, intent(in):: t_wrt ! fréquence d'écriture sur le fichier
35 guez 3
36 guez 40 ! Variables locales
37 guez 62 integer horiid, zvertiid
38     real julian
39 guez 67 integer iq, ii, jj
40 guez 3
41     !----------------------------------------------------
42    
43 guez 40 print *, "Call sequence information: initdynav"
44 guez 3
45 guez 62 CALL ymds2ju(anne0, 1, day0, 0., julian)
46     call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, &
47     rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, &
48     horiid, histaveid)
49 guez 67 call histvert(histaveid, 'sigss', 'Niveaux sigma', 'Pa', nivsigs, 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