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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 30 - (hide annotations)
Thu Apr 1 09:07:28 2010 UTC (14 years, 1 month ago) by guez
Original Path: trunk/libf/bibio/initdynav.f90
File size: 4009 byte(s)
Imported Source files of the external library "IOIPSL_Lionel" into
"libf/IOIPSL".

Split "cray.f90" into "scopy.f90" and "ssum.f90".

Rewrote "leapfrog" in order to have a clearer algorithmic structure.

1 guez 3 module initdynav_m
2    
3     ! This module is clean: no C preprocessor directive, no include line
4    
5     implicit none
6    
7     contains
8    
9 guez 27 subroutine initdynav(day0, anne0, tstep, nq, fileid, 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     ! L. Fairhead, LMD
13 guez 3
14 guez 30 USE histcom, ONLY: histbeg_totreg, histdef, histend, histvert
15     use calendar, ONLY: ymds2ju
16 guez 27 USE dimens_m, ONLY : llm
17     USE paramet_m, ONLY : iip1, jjp1
18     USE comconst, ONLY : pi
19     USE comvert, ONLY : nivsigs
20     USE comgeom, ONLY : rlatu, rlonv
21     USE temps, ONLY : itau_dyn
22     USE iniadvtrac_m, ONLY : ttext
23 guez 3
24     ! Routine d'initialisation des ecritures des fichiers histoires LMDZ
25     ! au format IOIPSL. Initialisation du fichier histoire moyenne.
26    
27     ! Appels succesifs des routines: histbeg
28     ! histhori
29     ! histver
30     ! histdef
31     ! histend
32    
33     ! Entree:
34     ! day0,anne0: date de reference
35     ! tstep : frequence d'ecriture
36     ! t_ops: frequence de l'operation pour IOIPSL
37     ! t_wrt: frequence d'ecriture sur le fichier
38     ! nq: nombre de traceurs
39    
40     ! Sortie:
41     ! fileid: ID du fichier netcdf cree
42    
43     ! Arguments
44     integer day0, anne0
45     real, intent(in):: tstep, t_ops, t_wrt
46     integer fileid
47     integer nq
48     integer thoriid, zvertiid
49    
50     ! Variables locales
51    
52     real zjulian
53     integer iq
54     real rlong(iip1,jjp1), rlat(iip1,jjp1)
55     integer ii,jj
56     integer zan, dayref
57    
58     !----------------------------------------------------
59    
60     ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
61    
62     zan = anne0
63     dayref = day0
64     CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
65    
66     do jj = 1, jjp1
67     do ii = 1, iip1
68     rlong(ii,jj) = rlonv(ii) * 180. / pi
69     rlat(ii,jj) = rlatu(jj) * 180. / pi
70     enddo
71     enddo
72    
73 guez 27 call histbeg_totreg('dyn_hist_ave.nc', rlong(:,1), rlat(1,:), &
74 guez 3 1, iip1, 1, jjp1, &
75 guez 27 itau_dyn, zjulian, tstep, thoriid, fileid)
76 guez 3
77 guez 27
78 guez 3 ! Appel a histvert pour la grille verticale
79 guez 27
80 guez 3 call histvert(fileid, 'sigss', 'Niveaux sigma','Pa', &
81     llm, nivsigs, zvertiid)
82 guez 27
83 guez 3 ! Appels a histdef pour la definition des variables a sauvegarder
84 guez 27
85 guez 3 ! Vents U
86 guez 27
87 guez 3 write(6,*)'inithistave',tstep
88     call histdef(fileid, 'u', 'vents u scalaires moyennes', &
89     'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
90 guez 15 'ave(X)', t_ops, t_wrt)
91 guez 3
92 guez 27
93 guez 3 ! Vents V
94 guez 27
95 guez 3 call histdef(fileid, 'v', 'vents v scalaires moyennes', &
96     'm/s', iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
97 guez 15 'ave(X)', t_ops, t_wrt)
98 guez 3
99 guez 27
100 guez 3 ! Temperature
101 guez 27
102 guez 3 call histdef(fileid, 'temp', 'temperature moyennee', 'K', &
103     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
104 guez 15 'ave(X)', t_ops, t_wrt)
105 guez 27
106 guez 3 ! Temperature potentielle
107 guez 27
108 guez 3 call histdef(fileid, 'theta', 'temperature potentielle', 'K', &
109     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
110 guez 15 'ave(X)', t_ops, t_wrt)
111 guez 3
112    
113 guez 27
114 guez 3 ! Geopotentiel
115 guez 27
116 guez 3 call histdef(fileid, 'phi', 'geopotentiel moyenne', '-', &
117     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
118 guez 15 'ave(X)', t_ops, t_wrt)
119 guez 27
120 guez 3 ! Traceurs
121 guez 27
122 guez 3 DO iq=1,nq
123     call histdef(fileid, ttext(iq), ttext(iq), '-', &
124     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
125 guez 15 'ave(X)', t_ops, t_wrt)
126 guez 3 enddo
127 guez 27
128 guez 3 ! Masse
129 guez 27
130 guez 3 call histdef(fileid, 'masse', 'masse', 'kg', &
131     iip1, jjp1, thoriid, 1, 1, 1, -99, &
132 guez 15 'ave(X)', t_ops, t_wrt)
133 guez 27
134 guez 3 ! Pression au sol
135 guez 27
136 guez 3 call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
137     iip1, jjp1, thoriid, 1, 1, 1, -99, &
138 guez 15 'ave(X)', t_ops, t_wrt)
139 guez 27
140 guez 3 ! Pression au sol
141 guez 27
142 guez 3 call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
143     iip1, jjp1, thoriid, 1, 1, 1, -99, &
144 guez 15 'ave(X)', t_ops, t_wrt)
145 guez 27
146 guez 3 call histend(fileid)
147    
148     end subroutine initdynav
149    
150     end module initdynav_m

  ViewVC Help
Powered by ViewVC 1.1.21