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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 1 month ago) by guez
Original Path: trunk/libf/bibio/initdynav.f90
File size: 3333 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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     ! Routine d'initialisation des écritures des fichiers histoires LMDZ
14     ! au format IOIPSL. Initialisation du fichier histoire moyenne.
15 guez 3
16 guez 30 use calendar, ONLY: ymds2ju
17 guez 39 USE comvert, ONLY : nivsigs
18     USE comgeom, ONLY : rlatu, rlonv
19 guez 27 USE dimens_m, ONLY : llm
20 guez 61 USE histbeg_totreg_m, ONLY : histbeg_totreg
21     USE histdef_m, ONLY : histdef
22     USE histend_m, ONLY : histend
23     USE histvert_m, ONLY : histvert
24 guez 39 USE iniadvtrac_m, ONLY : ttext
25     USE nr_util, ONLY : pi
26 guez 27 USE paramet_m, ONLY : iip1, jjp1
27     USE temps, ONLY : itau_dyn
28 guez 3
29 guez 40 integer, intent(in):: day0, anne0 ! date de reference
30     real, intent(in):: tstep ! frequence d'ecriture
31 guez 56 integer, intent(in):: nq ! nombre de traceurs
32 guez 40 real, intent(in):: t_ops ! frequence de l'operation pour IOIPSL
33     real, intent(in):: t_wrt ! frequence d'ecriture sur le fichier
34 guez 3
35 guez 40 ! Variables locales
36 guez 3 integer thoriid, zvertiid
37     real zjulian
38     integer iq
39 guez 40 real rlong(iip1, jjp1), rlat(iip1, jjp1)
40     integer ii, jj
41 guez 3 integer zan, dayref
42    
43     !----------------------------------------------------
44    
45 guez 40 print *, "Call sequence information: initdynav"
46 guez 3
47     zan = anne0
48     dayref = day0
49     CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
50    
51     do jj = 1, jjp1
52     do ii = 1, iip1
53 guez 40 rlong(ii, jj) = rlonv(ii) * 180. / pi
54     rlat(ii, jj) = rlatu(jj) * 180. / pi
55 guez 3 enddo
56     enddo
57    
58 guez 40 call histbeg_totreg('dyn_hist_ave.nc', rlong(:, 1), rlat(1, :), 1, iip1, &
59 guez 56 1, jjp1, itau_dyn, zjulian, tstep, thoriid, histaveid)
60     call histvert(histaveid, 'sigss', 'Niveaux sigma', 'Pa', llm, nivsigs, &
61 guez 40 zvertiid)
62 guez 3
63 guez 56 call histdef(histaveid, 'u', 'vents u scalaires moyennes', &
64 guez 3 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
65 guez 15 'ave(X)', t_ops, t_wrt)
66 guez 56 call histdef(histaveid, 'v', 'vents v scalaires moyennes', &
67 guez 3 'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
68 guez 15 'ave(X)', t_ops, t_wrt)
69 guez 56 call histdef(histaveid, 'temp', 'temperature moyennee', 'K', &
70 guez 3 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
71 guez 15 'ave(X)', t_ops, t_wrt)
72 guez 56 call histdef(histaveid, 'theta', 'temperature potentielle', 'K', &
73 guez 3 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
74 guez 15 'ave(X)', t_ops, t_wrt)
75 guez 56 call histdef(histaveid, 'phi', 'geopotentiel moyenne', '-', &
76 guez 3 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
77 guez 15 'ave(X)', t_ops, t_wrt)
78 guez 27
79 guez 40 ! Traceurs
80     DO iq = 1, nq
81 guez 56 call histdef(histaveid, ttext(iq), ttext(iq), '-', &
82 guez 3 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
83 guez 15 'ave(X)', t_ops, t_wrt)
84 guez 3 enddo
85 guez 27
86 guez 56 call histdef(histaveid, 'masse', 'masse', 'kg', &
87 guez 3 iip1, jjp1, thoriid, 1, 1, 1, -99, &
88 guez 15 'ave(X)', t_ops, t_wrt)
89 guez 56 call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', &
90 guez 3 iip1, jjp1, thoriid, 1, 1, 1, -99, &
91 guez 15 'ave(X)', t_ops, t_wrt)
92 guez 56 call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', &
93 guez 3 iip1, jjp1, thoriid, 1, 1, 1, -99, &
94 guez 15 'ave(X)', t_ops, t_wrt)
95 guez 27
96 guez 56 call histend(histaveid)
97 guez 3
98     end subroutine initdynav
99    
100     end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21