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

Contents of /trunk/bibio/initdynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 92 - (show annotations)
Wed Mar 26 18:16:05 2014 UTC (10 years, 1 month ago) by guez
File size: 2940 byte(s)
Extracted procedures that were in module calendar into separate files.

1 module initdynav_m
2
3 implicit none
4
5 integer histaveid
6
7 contains
8
9 subroutine initdynav(day0, anne0, tstep, nq, t_ops, t_wrt)
10
11 ! From initdynav.F, version 1.1.1.1, 2004/05/19 12:53:05
12 ! L. Fairhead, LMD
13
14 ! Routine d'initialisation des écritures des fichiers histoires au
15 ! format IOIPSL. Initialisation du fichier histoire moyenne.
16
17 USE comgeom, ONLY: rlatu, rlonv
18 USE dimens_m, ONLY: llm
19 USE histbeg_totreg_m, ONLY: histbeg_totreg
20 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 use ymds2ju_m, ONLY: ymds2ju
28
29 integer, intent(in):: day0, anne0 ! date de référence
30 real, intent(in):: tstep ! fréquence d'écriture
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 integer horiid, zvertiid
37 real julian
38 integer iq, ii, jj, l
39
40 !----------------------------------------------------
41
42 print *, "Call sequence information: initdynav"
43
44 CALL ymds2ju(anne0, 1, day0, 0., julian)
45 call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, &
46 rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, &
47 horiid, histaveid)
48 call histvert(histaveid, 'sigss', 'Niveaux sigma', '', &
49 (/(real(l), l = 1, llm)/), zvertiid)
50
51 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
62 ! Traceurs
63 DO iq = 1, nq
64 call histdef(histaveid, ttext(iq), ttext(iq), '-', iip1, jjp1, &
65 horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
66 enddo
67
68 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
75 call histend(histaveid)
76
77 end subroutine initdynav
78
79 end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21