1 | c======================================================================= |
---|
2 | PROGRAM start2archive |
---|
3 | c======================================================================= |
---|
4 | c |
---|
5 | c |
---|
6 | c Date: 01/1997 |
---|
7 | c ---- |
---|
8 | c |
---|
9 | c Version Venus: 09/2007 |
---|
10 | c Titan: 02/2009 |
---|
11 | c |
---|
12 | c Objet: Passage des fichiers netcdf d etat initial "start" et |
---|
13 | c ----- "startphy" a un fichier netcdf unique "start_archive" |
---|
14 | c |
---|
15 | c======================================================================= |
---|
16 | |
---|
17 | USE filtreg_mod |
---|
18 | USE infotrac |
---|
19 | USE control_mod |
---|
20 | use cpdet_mod, only: tpot2t,ini_cpdet |
---|
21 | |
---|
22 | implicit none |
---|
23 | |
---|
24 | #include "dimensions.h" |
---|
25 | #include "paramet.h" |
---|
26 | #include "comconst.h" |
---|
27 | #include "comdissnew.h" |
---|
28 | #include "comvert.h" |
---|
29 | #include "comgeom.h" |
---|
30 | #include "logic.h" |
---|
31 | #include "temps.h" |
---|
32 | #include "ener.h" |
---|
33 | #include "description.h" |
---|
34 | #include "dimsoil.h" |
---|
35 | #include "netcdf.inc" |
---|
36 | |
---|
37 | c----------------------------------------------------------------------- |
---|
38 | c Declarations |
---|
39 | c----------------------------------------------------------------------- |
---|
40 | |
---|
41 | c variables dynamiques du GCM |
---|
42 | c ----------------------------- |
---|
43 | REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants |
---|
44 | REAL teta(ip1jmp1,llm) ! temperature potentielle |
---|
45 | REAL, ALLOCATABLE, DIMENSION(:,:,:):: q! champs advectes |
---|
46 | REAL pks(ip1jmp1) ! exner (f pour filtre) |
---|
47 | REAL pk(ip1jmp1,llm) |
---|
48 | REAL pkf(ip1jmp1,llm) |
---|
49 | REAL alpha(iip1,jjp1,llm),beta(iip1,jjp1,llm) |
---|
50 | REAL phis(ip1jmp1) ! geopotentiel au sol |
---|
51 | REAL masse(ip1jmp1,llm) ! masse de l'atmosphere |
---|
52 | REAL ps(ip1jmp1) ! pression au sol |
---|
53 | REAL p3d(iip1,jjp1,llm+1) ! pression aux interfaces |
---|
54 | |
---|
55 | c Variable Physiques (grille physique) |
---|
56 | c ------------------------------------ |
---|
57 | integer ngridmx,nlayermx |
---|
58 | parameter (ngridmx=(2+(jjm-1)*iim - 1/jjm)) |
---|
59 | parameter (nlayermx=llm) |
---|
60 | |
---|
61 | real rlat(ngridmx),rlon(ngridmx) |
---|
62 | REAL tsurf(ngridmx),tsoil(ngridmx,nsoilmx) |
---|
63 | REAL albe(ngridmx),radsol(ngridmx),sollw(ngridmx) |
---|
64 | real solsw(ngridmx),fder(ngridmx) |
---|
65 | real sollwdown(ngridmx),dlw(ngridmx) |
---|
66 | REAL zmea(ngridmx), zstd(ngridmx) |
---|
67 | REAL zsig(ngridmx), zgam(ngridmx), zthe(ngridmx) |
---|
68 | REAL zpic(ngridmx), zval(ngridmx) |
---|
69 | |
---|
70 | INTEGER start,length |
---|
71 | PARAMETER (length = 100) |
---|
72 | REAL tab_cntrl_fi(length) ! tableau des parametres de startfi |
---|
73 | REAL tab_cntrl_dyn(length) ! tableau des parametres de start |
---|
74 | INTEGER*4 day_ini_fi |
---|
75 | |
---|
76 | c Variable naturelle / grille scalaire |
---|
77 | c ------------------------------------ |
---|
78 | REAL T(ip1jmp1,llm),us(ip1jmp1,llm),vs(ip1jmp1,llm) |
---|
79 | REAL tsurfS(ip1jmp1),tsoilS(ip1jmp1,nsoilmx) |
---|
80 | real rlatS(ip1jmp1),rlonS(ip1jmp1) |
---|
81 | real albeS(ip1jmp1),radsolS(ip1jmp1),sollwS(ip1jmp1) |
---|
82 | real solswS(ip1jmp1),fderS(ip1jmp1) |
---|
83 | real dlwS(ip1jmp1),sollwdownS(ip1jmp1) |
---|
84 | real zmeaS(ip1jmp1),zstdS(ip1jmp1),zsigS(ip1jmp1) |
---|
85 | real zgamS(ip1jmp1),ztheS(ip1jmp1),zpicS(ip1jmp1) |
---|
86 | real zvalS(ip1jmp1) |
---|
87 | |
---|
88 | c Variables intermediaires : vent naturel, mais pas coord scalaire |
---|
89 | c---------------------------------------------------------------- |
---|
90 | REAL vn(ip1jm,llm),un(ip1jmp1,llm) |
---|
91 | |
---|
92 | c Autres variables |
---|
93 | c ----------------- |
---|
94 | REAL ptotal |
---|
95 | |
---|
96 | CHARACTER*2 str2 |
---|
97 | |
---|
98 | INTEGER ij, l,i,j,isoil,iq |
---|
99 | character*80 fichnom |
---|
100 | integer :: ierr |
---|
101 | |
---|
102 | c Netcdf |
---|
103 | c------- |
---|
104 | integer varid,dimid |
---|
105 | INTEGER nid |
---|
106 | |
---|
107 | c----------------------------------------------------------------------- |
---|
108 | c Initialisations |
---|
109 | c----------------------------------------------------------------------- |
---|
110 | |
---|
111 | c VENUS/TITAN |
---|
112 | |
---|
113 | iflag_trac = 1 |
---|
114 | c----------------------------------------------------------------------- |
---|
115 | c Initialisation des traceurs |
---|
116 | c --------------------------- |
---|
117 | c Choix du nombre de traceurs et du schema pour l advection |
---|
118 | c dans fichier traceur.def, par default ou via INCA |
---|
119 | call infotrac_init |
---|
120 | |
---|
121 | c Allocation de la tableau q : champs advectes |
---|
122 | allocate(q(ip1jmp1,llm,nqtot)) |
---|
123 | |
---|
124 | c======================================================================= |
---|
125 | c Lecture des donnees |
---|
126 | c======================================================================= |
---|
127 | |
---|
128 | fichnom = 'start.nc' |
---|
129 | CALL readstart(fichnom,nqtot,vcov,ucov,teta,q,masse, |
---|
130 | . ps,phis,tab_cntrl_dyn) |
---|
131 | |
---|
132 | fichnom = 'startphy.nc' |
---|
133 | CALL readstartphy(fichnom, |
---|
134 | . rlat,rlon,tsurf,tsoil, |
---|
135 | . albe, solsw, sollw, |
---|
136 | . fder,dlw,sollwdown,radsol, |
---|
137 | . zmea,zstd,zsig,zgam,zthe,zpic,zval, |
---|
138 | . tab_cntrl_fi) |
---|
139 | |
---|
140 | c----------------------------------------------------------------------- |
---|
141 | c Initialisations |
---|
142 | c----------------------------------------------------------------------- |
---|
143 | |
---|
144 | CALL conf_gcm( 99, .TRUE. ) |
---|
145 | call iniconst |
---|
146 | call inigeom |
---|
147 | call inifilr |
---|
148 | call ini_cpdet |
---|
149 | |
---|
150 | CALL pression(ip1jmp1, ap, bp, ps, p3d) |
---|
151 | if (disvert_type==1) then |
---|
152 | CALL exner_hyb( ip1jmp1, ps, p3d,alpha,beta,pks, pk, pkf ) |
---|
153 | else ! we assume that we are in the disvert_type==2 case |
---|
154 | CALL exner_milieu( ip1jmp1, ps, p3d, beta, pks, pk, pkf ) |
---|
155 | endif |
---|
156 | |
---|
157 | c======================================================================= |
---|
158 | c Transformation EN VARIABLE NATURELLE / GRILLE SCALAIRE si necessaire |
---|
159 | c======================================================================= |
---|
160 | c Les variables modeles dependent de la resolution. Il faut donc |
---|
161 | c eliminer les facteurs responsables de cette dependance |
---|
162 | c (pour utiliser newstart) |
---|
163 | c======================================================================= |
---|
164 | |
---|
165 | c----------------------------------------------------------------------- |
---|
166 | c Vent (depend de la resolution horizontale) |
---|
167 | c----------------------------------------------------------------------- |
---|
168 | c |
---|
169 | c ucov --> un et vcov --> vn |
---|
170 | c un --> us et vn --> vs |
---|
171 | c |
---|
172 | c----------------------------------------------------------------------- |
---|
173 | |
---|
174 | call covnat(llm,ucov, vcov, un, vn) |
---|
175 | call wind_scal(un,vn,us,vs) |
---|
176 | |
---|
177 | c----------------------------------------------------------------------- |
---|
178 | c Temperature (depend de la resolution verticale => de "sigma.def") |
---|
179 | c----------------------------------------------------------------------- |
---|
180 | c |
---|
181 | c h --> T |
---|
182 | c |
---|
183 | c----------------------------------------------------------------------- |
---|
184 | ! ADAPTATION GCM POUR CP(T) |
---|
185 | |
---|
186 | call tpot2t(ip1jmp1*llm,teta,T,pk) |
---|
187 | |
---|
188 | c----------------------------------------------------------------------- |
---|
189 | c Variable physique |
---|
190 | c----------------------------------------------------------------------- |
---|
191 | c |
---|
192 | c tsurf --> tsurfS |
---|
193 | c et autres... |
---|
194 | c |
---|
195 | c----------------------------------------------------------------------- |
---|
196 | |
---|
197 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,tsurf,tsurfS) |
---|
198 | call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,tsoil,tsoilS) |
---|
199 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,rlat,rlatS) |
---|
200 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,rlon,rlonS) |
---|
201 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,albe,albeS) |
---|
202 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,radsol,radsolS) |
---|
203 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,sollw,sollwS) |
---|
204 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,solsw,solswS) |
---|
205 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,fder,fderS) |
---|
206 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,dlw,dlwS) |
---|
207 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,sollwdown,sollwdownS) |
---|
208 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zmea,zmeaS) |
---|
209 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zstd,zstdS) |
---|
210 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zsig,zsigS) |
---|
211 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zgam,zgamS) |
---|
212 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zthe,ztheS) |
---|
213 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zpic,zpicS) |
---|
214 | call gr_fi_dyn(1,ngridmx,iip1,jjp1,zval,zvalS) |
---|
215 | |
---|
216 | c======================================================================= |
---|
217 | c Info pour controler |
---|
218 | c======================================================================= |
---|
219 | |
---|
220 | ptotal = 0. |
---|
221 | DO j=1,jjp1 |
---|
222 | DO i=1,iim |
---|
223 | ptotal=ptotal+aire(i+(iim+1)*(j-1))*ps(i+(iim+1)*(j-1))/g |
---|
224 | ENDDO |
---|
225 | ENDDO |
---|
226 | write(*,*)'Ancienne grille : masse de l''atm :',ptotal |
---|
227 | |
---|
228 | c----------------------------------------------------------------------- |
---|
229 | c Passage de "ptotal" par tab_cntrl_fi |
---|
230 | c----------------------------------------------------------------------- |
---|
231 | |
---|
232 | tab_cntrl_fi(length) = ptotal |
---|
233 | |
---|
234 | c======================================================================= |
---|
235 | c Ecriture dans le fichier "start_archive" |
---|
236 | c======================================================================= |
---|
237 | |
---|
238 | c----------------------------------------------------------------------- |
---|
239 | c Ouverture de "start_archive" |
---|
240 | c----------------------------------------------------------------------- |
---|
241 | |
---|
242 | ierr = NF_OPEN ('start_archive.nc', NF_WRITE,nid) |
---|
243 | |
---|
244 | c----------------------------------------------------------------------- |
---|
245 | c si "start_archive" n'existe pas: |
---|
246 | c 1_ ouverture |
---|
247 | c 2_ creation de l'entete dynamique ("ini_archive") |
---|
248 | c----------------------------------------------------------------------- |
---|
249 | c ini_archive: |
---|
250 | c On met dans l'entete le tab_cntrl_dyn (1 a length) |
---|
251 | c On y ajoute les valeurs du tab_cntrl_fi (length+1 a 2*length) |
---|
252 | c----------------------------------------------------------------------- |
---|
253 | |
---|
254 | if (ierr.ne.NF_NOERR) then |
---|
255 | write(*,*)'OK, Could not open file "start_archive.nc"' |
---|
256 | write(*,*)'So let s create a new "start_archive"' |
---|
257 | ierr = NF_CREATE('start_archive.nc', NF_CLOBBER, nid) |
---|
258 | call ini_archive(nid,phis,tab_cntrl_dyn,tab_cntrl_fi) |
---|
259 | else |
---|
260 | write(*,*)'Attention, start_archive.nc existe déjà...' |
---|
261 | call abort |
---|
262 | endif |
---|
263 | |
---|
264 | c----------------------------------------------------------------------- |
---|
265 | c Ecriture des champs |
---|
266 | c----------------------------------------------------------------------- |
---|
267 | |
---|
268 | call write_archive(nid,'u','Vent zonal','m.s-1',3,us) |
---|
269 | call write_archive(nid,'v','Vent merid','m.s-1',3,vs) |
---|
270 | call write_archive(nid,'temp','temperature','K',3,T) |
---|
271 | c----------------------------------------------------------------------- |
---|
272 | c Ecriture du champs q ( q[1,nqtot] ) |
---|
273 | c----------------------------------------------------------------------- |
---|
274 | do iq=1,nqtot |
---|
275 | write(str2,'(i2.2)') iq |
---|
276 | call write_archive(nid,tname(iq),'tracer','kg/kg', |
---|
277 | . 3,q(1,1,iq)) |
---|
278 | end do |
---|
279 | c----------------------------------------------------------------------- |
---|
280 | call write_archive(nid,'masse','Masse','kg',3,masse) |
---|
281 | call write_archive(nid,'ps','Psurf','Pa',2,ps) |
---|
282 | call write_archive(nid,'tsurf','surf T','K',2,tsurfS) |
---|
283 | c----------------------------------------------------------------------- |
---|
284 | c Ecriture du champs tsoil ( Tsoil[1,nsoilmx] ) |
---|
285 | c----------------------------------------------------------------------- |
---|
286 | c "tsoil" Temperature au sol definie dans nsoilmx couches dans le sol |
---|
287 | c Les nsoilmx couches sont lues comme nsoilmx champs |
---|
288 | c nommees Tsoil[1,nsoilmx] |
---|
289 | do isoil=1,nsoilmx |
---|
290 | write(str2,'(i2.2)') isoil |
---|
291 | call write_archive(nid,'Tsoil'//str2,'Ground Temperature ', |
---|
292 | . 'K',2,tsoilS(1,isoil)) |
---|
293 | enddo |
---|
294 | c----------------------------------------------------------------------- |
---|
295 | call write_archive(nid,'rlat','Latitude','rad',2,rlatS) |
---|
296 | call write_archive(nid,'rlon','Longitude','rad',2,rlonS) |
---|
297 | call write_archive(nid,'albe','Albedo','',2,albeS) |
---|
298 | call write_archive(nid,'radsol', |
---|
299 | . 'Net flux at surface','W m-2',2,radsolS) |
---|
300 | call write_archive(nid,'sollw', |
---|
301 | . 'LW flux at surface','W m-2',2,sollwS) |
---|
302 | call write_archive(nid,'solsw', |
---|
303 | . 'SW flux at surface','W m-2',2,solswS) |
---|
304 | call write_archive(nid,'fder','derive','?',2,fderS) |
---|
305 | call write_archive(nid,'dlw','LW derive','?',2,dlwS) |
---|
306 | call write_archive(nid,'sollwdown', |
---|
307 | . 'LW dwn flux at surface','?',2,sollwdownS) |
---|
308 | call write_archive(nid,'zmea','param oro sous-maille','m',2,zmeaS) |
---|
309 | call write_archive(nid,'zstd','param oro sous-maille','m',2,zstdS) |
---|
310 | call write_archive(nid,'zsig','param oro sous-maille','m',2,zsigS) |
---|
311 | call write_archive(nid,'zgam','param oro sous-maille','m',2,zgamS) |
---|
312 | call write_archive(nid,'zthe','param oro sous-maille','m',2,ztheS) |
---|
313 | call write_archive(nid,'zpic','param oro sous-maille','m',2,zpicS) |
---|
314 | call write_archive(nid,'zval','param oro sous-maille','m',2,zvalS) |
---|
315 | |
---|
316 | ierr=NF_CLOSE(nid) |
---|
317 | |
---|
318 | c----------------------------------------------------------------------- |
---|
319 | c Fin |
---|
320 | c----------------------------------------------------------------------- |
---|
321 | |
---|
322 | end |
---|