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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 144 - (hide annotations)
Wed Jun 10 16:46:46 2015 UTC (8 years, 11 months ago) by guez
File size: 2889 byte(s)
In procedure fxhyp, the convoluted computation of tanh(fa/fb) occurred
three times. Extracted it into a function. Also, the computation of
xmoy and fxm was repeated. So stored the values in arrays instead.

In procedure fxhyp, in the computation of fhyp, there were tests
xtild(i) == 0. and xtild(i) == pi_d. No use to do these tests at each
iteration. We now they are true for i == nmax and i == 2 * nmax,
respectively, and we know they are false for other values of
"i". Similarly, in the computations of ffdx and xxpr, there were the
tests xmoy == 0. and xmoy == pi_d, these could not be true.

Moved files from bibio to dyn3d, following LMDZ.

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 67 ! Routine d'initialisation des écritures des fichiers histoires au
15     ! format IOIPSL. Initialisation du fichier histoire moyenne.
16    
17 guez 62 USE dimens_m, ONLY: llm
18 guez 139 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlonv
19 guez 62 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 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 40 ! Variables locales
35 guez 62 integer horiid, zvertiid
36     real julian
37 guez 105 integer iq, l
38 guez 3
39     !----------------------------------------------------
40    
41 guez 40 print *, "Call sequence information: initdynav"
42 guez 3
43 guez 129 CALL ymds2ju(annee_ref, 1, day_ref, 0., julian)
44 guez 62 call histbeg_totreg('dyn_hist_ave.nc', rlonv * 180. / pi, &
45     rlatu * 180. / pi, 1, iip1, 1, jjp1, itau_dyn, julian, tstep, &
46     horiid, histaveid)
47 guez 83 call histvert(histaveid, 'sigss', 'Niveaux sigma', '', &
48     (/(real(l), l = 1, llm)/), zvertiid)
49 guez 3
50 guez 62 call histdef(histaveid, 'u', 'vents u scalaires moyennes', 'm/s', iip1, &
51     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
52     call histdef(histaveid, 'v', 'vents v scalaires moyennes', 'm/s', iip1, &
53     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
54     call histdef(histaveid, 'temp', 'temperature moyennee', 'K', iip1, jjp1, &
55     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
56     call histdef(histaveid, 'theta', 'temperature potentielle', 'K', iip1, &
57     jjp1, horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
58     call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', iip1, jjp1, &
59     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
60 guez 27
61 guez 40 ! Traceurs
62     DO iq = 1, nq
63 guez 62 call histdef(histaveid, ttext(iq), ttext(iq), '-', iip1, jjp1, &
64     horiid, llm, 1, llm, zvertiid, 'ave(X)', t_ops, t_wrt)
65 guez 3 enddo
66 guez 27
67 guez 62 call histdef(histaveid, 'masse', 'masse', 'kg', iip1, jjp1, horiid, 1, &
68     1, 1, -99, 'ave(X)', t_ops, t_wrt)
69     call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', iip1, &
70     jjp1, horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
71     call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
72     horiid, 1, 1, 1, -99, 'ave(X)', t_ops, t_wrt)
73 guez 27
74 guez 56 call histend(histaveid)
75 guez 3
76     end subroutine initdynav
77    
78     end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21