1 | module physiq_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | subroutine physiq(ngrid,nlayer,nq, & |
---|
8 | nametrac, & |
---|
9 | firstcall,lastcall, & |
---|
10 | pday,ptime,ptimestep, & |
---|
11 | pplev,pplay,pphi, & |
---|
12 | pu,pv,pt,pq, & |
---|
13 | flxw, & |
---|
14 | pdu,pdv,pdt,pdq,pdpsrf) |
---|
15 | |
---|
16 | use radinc_h, only : L_NSPECTI,L_NSPECTV |
---|
17 | use radcommon_h, only: sigma, gzlat, gzlat_ig, Cmk, grav, BWNV |
---|
18 | use surfdat_h, only: phisfi, zmea, zstd, zsig, zgam, zthe |
---|
19 | use comchem_h, only: nkim, cnames, nlaykim_up, ykim_up, ykim_tot, botCH4 |
---|
20 | use comdiurn_h, only: coslat, sinlat, coslon, sinlon |
---|
21 | use comsaison_h, only: mu0, fract, dist_star, declin, right_ascen |
---|
22 | use comsoil_h, only: nsoilmx, layer, mlayer, inertiedat |
---|
23 | use datafile_mod, only: datadir, corrkdir, banddir |
---|
24 | use geometry_mod, only: latitude, longitude, cell_area |
---|
25 | USE comgeomfi_h, only: totarea, totarea_planet |
---|
26 | USE tracer_h |
---|
27 | use time_phylmdz_mod, only: ecritphy, iphysiq, nday |
---|
28 | use phyetat0_mod, only: phyetat0 |
---|
29 | use phyredem, only: physdem0, physdem1 |
---|
30 | use planetwide_mod, only: planetwide_minval,planetwide_maxval,planetwide_sumval |
---|
31 | use mod_phys_lmdz_para, only : is_master |
---|
32 | use planete_mod, only: apoastr, periastr, year_day, peri_day, & |
---|
33 | obliquit, nres, z0 |
---|
34 | use comcstfi_mod, only: pi, g, rcp, r, rad, mugaz, cpp |
---|
35 | use time_phylmdz_mod, only: daysec |
---|
36 | use logic_mod, only: moyzon_ch |
---|
37 | use moyzon_mod, only: zphibar, zphisbar, zplevbar, zplaybar, & |
---|
38 | zzlevbar, zzlaybar, ztfibar, zqfibar |
---|
39 | use callkeys_mod |
---|
40 | use vertical_layers_mod, only: presnivs, pseudoalt |
---|
41 | use ioipsl_getin_p_mod, only: getin_p |
---|
42 | use mod_phys_lmdz_omp_data, ONLY: is_omp_master |
---|
43 | #ifdef CPP_XIOS |
---|
44 | use xios_output_mod, only: initialize_xios_output, & |
---|
45 | update_xios_timestep, & |
---|
46 | send_xios_field |
---|
47 | use wxios, only: wxios_context_init, xios_context_finalize |
---|
48 | #endif |
---|
49 | use MMP_OPTICS |
---|
50 | use muphy_diag |
---|
51 | implicit none |
---|
52 | |
---|
53 | |
---|
54 | !================================================================== |
---|
55 | ! |
---|
56 | ! Purpose |
---|
57 | ! ------- |
---|
58 | ! Central subroutine for all the physics parameterisations in the |
---|
59 | ! universal model. Originally adapted from the Mars LMDZ model. |
---|
60 | ! |
---|
61 | ! The model can be run without or with tracer transport |
---|
62 | ! depending on the value of "tracer" in file "callphys.def". |
---|
63 | ! |
---|
64 | ! |
---|
65 | ! It includes: |
---|
66 | ! |
---|
67 | ! I. Initialization : |
---|
68 | ! I.1 Firstcall initializations. |
---|
69 | ! I.2 Initialization for every call to physiq. |
---|
70 | ! |
---|
71 | ! II. Compute radiative transfer tendencies (longwave and shortwave) : |
---|
72 | ! II.a Option 1 : Call correlated-k radiative transfer scheme. |
---|
73 | ! II.b Option 2 : Call Newtonian cooling scheme. |
---|
74 | ! II.c Option 3 : Atmosphere has no radiative effect. |
---|
75 | ! |
---|
76 | ! III. Vertical diffusion (turbulent mixing) : |
---|
77 | ! |
---|
78 | ! IV. Dry Convective adjustment : |
---|
79 | ! |
---|
80 | ! V. Tracers |
---|
81 | ! V.1. Microphysics |
---|
82 | ! V.2. Chemistry |
---|
83 | ! V.3. Updates (pressure variations, surface budget). |
---|
84 | ! V.4. Surface Tracer Update. |
---|
85 | ! |
---|
86 | ! VI. Surface and sub-surface soil temperature. |
---|
87 | ! |
---|
88 | ! VII. Perform diagnostics and write output files. |
---|
89 | ! |
---|
90 | ! |
---|
91 | ! arguments |
---|
92 | ! --------- |
---|
93 | ! |
---|
94 | ! INPUT |
---|
95 | ! ----- |
---|
96 | ! |
---|
97 | ! ngrid Size of the horizontal grid. |
---|
98 | ! nlayer Number of vertical layers. |
---|
99 | ! nq Number of advected fields. |
---|
100 | ! nametrac Name of corresponding advected fields. |
---|
101 | ! |
---|
102 | ! firstcall True at the first call. |
---|
103 | ! lastcall True at the last call. |
---|
104 | ! |
---|
105 | ! pday Number of days counted from the North. Spring equinoxe. |
---|
106 | ! ptime Universal time (0<ptime<1): ptime=0.5 at 12:00 UT. |
---|
107 | ! ptimestep timestep (s). |
---|
108 | ! |
---|
109 | ! pplay(ngrid,nlayer) Pressure at the middle of the layers (Pa). |
---|
110 | ! pplev(ngrid,nlayer+1) Intermediate pressure levels (Pa). |
---|
111 | ! pphi(ngrid,nlayer) Geopotential at the middle of the layers (m2.s-2). |
---|
112 | ! |
---|
113 | ! pu(ngrid,nlayer) u, zonal component of the wind (ms-1). |
---|
114 | ! pv(ngrid,nlayer) v, meridional component of the wind (ms-1). |
---|
115 | ! |
---|
116 | ! pt(ngrid,nlayer) Temperature (K). |
---|
117 | ! |
---|
118 | ! pq(ngrid,nlayer,nq) Advected fields. |
---|
119 | ! |
---|
120 | ! pudyn(ngrid,nlayer) \ |
---|
121 | ! pvdyn(ngrid,nlayer) \ Dynamical temporal derivative for the |
---|
122 | ! ptdyn(ngrid,nlayer) / corresponding variables. |
---|
123 | ! pqdyn(ngrid,nlayer,nq) / |
---|
124 | ! flxw(ngrid,nlayer) vertical mass flux (kg/s) at layer lower boundary |
---|
125 | ! |
---|
126 | ! OUTPUT |
---|
127 | ! ------ |
---|
128 | ! |
---|
129 | ! pdu(ngrid,nlayer) \ |
---|
130 | ! pdv(ngrid,nlayer) \ Temporal derivative of the corresponding |
---|
131 | ! pdt(ngrid,nlayer) / variables due to physical processes. |
---|
132 | ! pdq(ngrid,nlayer) / |
---|
133 | ! pdpsrf(ngrid) / |
---|
134 | ! |
---|
135 | ! |
---|
136 | ! Authors |
---|
137 | ! ------- |
---|
138 | ! Frederic Hourdin 15/10/93 |
---|
139 | ! Francois Forget 1994 |
---|
140 | ! Christophe Hourdin 02/1997 |
---|
141 | ! Subroutine completely rewritten by F. Forget (01/2000) |
---|
142 | ! Water ice clouds: Franck Montmessin (update 06/2003) |
---|
143 | ! Radiatively active tracers: J.-B. Madeleine (10/2008-06/2009) |
---|
144 | ! New correlated-k radiative scheme: R. Wordsworth (2009) |
---|
145 | ! Many specifically Martian subroutines removed: R. Wordsworth (2009) |
---|
146 | ! Improved water cycle: R. Wordsworth / B. Charnay (2010) |
---|
147 | ! To F90: R. Wordsworth (2010) |
---|
148 | ! New turbulent diffusion scheme: J. Leconte (2012) |
---|
149 | ! Loops converted to F90 matrix format: J. Leconte (2012) |
---|
150 | ! No more ngridmx/nqmx, F90 commons and adaptation to parallel: A. Spiga (2012) |
---|
151 | ! Purge of the code : M. Turbet (2015) |
---|
152 | ! Fork for Titan : J. Vatant d'Ollone (2017) |
---|
153 | ! + clean of all too-generic (ocean, water, co2 ...) routines |
---|
154 | ! + Titan's chemistry |
---|
155 | ! Microphysical moment model - J.Burgalat / J.Vatant d'Ollone (2017-2018) |
---|
156 | !============================================================================================ |
---|
157 | |
---|
158 | |
---|
159 | ! 0. Declarations : |
---|
160 | ! ------------------ |
---|
161 | |
---|
162 | include "netcdf.inc" |
---|
163 | |
---|
164 | ! Arguments : |
---|
165 | ! ----------- |
---|
166 | |
---|
167 | ! INPUTS: |
---|
168 | ! ------- |
---|
169 | |
---|
170 | |
---|
171 | integer,intent(in) :: ngrid ! Number of atmospheric columns. |
---|
172 | integer,intent(in) :: nlayer ! Number of atmospheric layers. |
---|
173 | integer,intent(in) :: nq ! Number of tracers. |
---|
174 | character*30,intent(in) :: nametrac(nq) ! Names of the tracers taken from dynamics. |
---|
175 | |
---|
176 | logical,intent(in) :: firstcall ! Signals first call to physics. |
---|
177 | logical,intent(in) :: lastcall ! Signals last call to physics. |
---|
178 | |
---|
179 | real,intent(in) :: pday ! Number of elapsed sols since reference Ls=0. |
---|
180 | real,intent(in) :: ptime ! "Universal time", given as fraction of sol (e.g.: 0.5 for noon). |
---|
181 | real,intent(in) :: ptimestep ! Physics timestep (s). |
---|
182 | real,intent(in) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa). |
---|
183 | real,intent(in) :: pplay(ngrid,nlayer) ! Mid-layer pressure (Pa). |
---|
184 | real,intent(in) :: pphi(ngrid,nlayer) ! Geopotential at mid-layer (m2s-2). |
---|
185 | real,intent(in) :: pu(ngrid,nlayer) ! Zonal wind component (m/s). |
---|
186 | real,intent(in) :: pv(ngrid,nlayer) ! Meridional wind component (m/s). |
---|
187 | real,intent(in) :: pt(ngrid,nlayer) ! Temperature (K). |
---|
188 | real,intent(in) :: pq(ngrid,nlayer,nq) ! Tracers (kg/kg_of_air). |
---|
189 | real,intent(in) :: flxw(ngrid,nlayer) ! Vertical mass flux (ks/s) at lower boundary of layer |
---|
190 | |
---|
191 | ! OUTPUTS: |
---|
192 | ! -------- |
---|
193 | |
---|
194 | ! Physical tendencies : |
---|
195 | |
---|
196 | real,intent(out) :: pdu(ngrid,nlayer) ! Zonal wind tendencies (m/s/s). |
---|
197 | real,intent(out) :: pdv(ngrid,nlayer) ! Meridional wind tendencies (m/s/s). |
---|
198 | real,intent(out) :: pdt(ngrid,nlayer) ! Temperature tendencies (K/s). |
---|
199 | real,intent(out) :: pdq(ngrid,nlayer,nq) ! Tracer tendencies (kg/kg_of_air/s). |
---|
200 | real,intent(out) :: pdpsrf(ngrid) ! Surface pressure tendency (Pa/s). |
---|
201 | |
---|
202 | ! Local saved variables: |
---|
203 | ! ---------------------- |
---|
204 | |
---|
205 | integer,save :: day_ini ! Initial date of the run (sol since Ls=0). |
---|
206 | integer,save :: icount ! Counter of calls to physiq during the run. |
---|
207 | !$OMP THREADPRIVATE(day_ini,icount) |
---|
208 | |
---|
209 | real, dimension(:),allocatable,save :: tsurf ! Surface temperature (K). |
---|
210 | real, dimension(:,:),allocatable,save :: tsoil ! Sub-surface temperatures (K). |
---|
211 | real, dimension(:,:),allocatable,save :: albedo ! Surface Spectral albedo. By MT2015. |
---|
212 | real, dimension(:),allocatable,save :: albedo_equivalent ! Spectral Mean albedo. |
---|
213 | |
---|
214 | !$OMP THREADPRIVATE(tsurf,tsoil,albedo,albedo_equivalent) |
---|
215 | |
---|
216 | real,dimension(:),allocatable,save :: albedo_bareground ! Bare Ground Albedo. By MT 2015. |
---|
217 | |
---|
218 | !$OMP THREADPRIVATE(albedo_bareground) |
---|
219 | |
---|
220 | real,dimension(:),allocatable,save :: emis ! Thermal IR surface emissivity. |
---|
221 | real,dimension(:,:),allocatable,save :: dtrad ! Net atmospheric radiative heating rate (K.s-1). |
---|
222 | real,dimension(:),allocatable,save :: fluxrad_sky ! Radiative flux from sky absorbed by surface (W.m-2). |
---|
223 | real,dimension(:),allocatable,save :: fluxrad ! Net radiative surface flux (W.m-2). |
---|
224 | real,dimension(:),allocatable,save :: capcal ! Surface heat capacity (J m-2 K-1). |
---|
225 | real,dimension(:),allocatable,save :: fluxgrd ! Surface conduction flux (W.m-2). |
---|
226 | real,dimension(:,:),allocatable,save :: qsurf ! Tracer on surface (e.g. kg.m-2). |
---|
227 | real,dimension(:,:),allocatable,save :: q2 ! Turbulent Kinetic Energy. |
---|
228 | |
---|
229 | !$OMP THREADPRIVATE(emis,dtrad,fluxrad_sky,fluxrad,capcal,fluxgrd,qsurf,q2) |
---|
230 | |
---|
231 | |
---|
232 | ! Local variables : |
---|
233 | ! ----------------- |
---|
234 | |
---|
235 | real zh(ngrid,nlayer) ! Potential temperature (K). |
---|
236 | real pw(ngrid,nlayer) ! Vertical velocity (m/s). (NOTE : >0 WHEN DOWNWARDS !!) |
---|
237 | |
---|
238 | integer l,ig,ierr,iq,nw,isoil |
---|
239 | |
---|
240 | ! FOR DIAGNOSTIC : |
---|
241 | |
---|
242 | real,dimension(:),allocatable,save :: fluxsurf_lw ! Incident Long Wave (IR) surface flux (W.m-2). |
---|
243 | real,dimension(:),allocatable,save :: fluxsurf_sw ! Incident Short Wave (stellar) surface flux (W.m-2). |
---|
244 | real,dimension(:),allocatable,save :: fluxsurfabs_sw ! Absorbed Short Wave (stellar) flux by the surface (W.m-2). |
---|
245 | real,dimension(:),allocatable,save :: fluxtop_lw ! Outgoing LW (IR) flux to space (W.m-2). |
---|
246 | real,dimension(:),allocatable,save :: fluxabs_sw ! Absorbed SW (stellar) flux (W.m-2). |
---|
247 | real,dimension(:),allocatable,save :: fluxtop_dn ! Incoming SW (stellar) radiation at the top of the atmosphere (W.m-2). |
---|
248 | real,dimension(:),allocatable,save :: fluxdyn ! Horizontal heat transport by dynamics (W.m-2). |
---|
249 | real,dimension(:,:),allocatable,save :: OLR_nu ! Outgoing LW radiation in each band (Normalized to the band width (W/m2/cm-1)). |
---|
250 | real,dimension(:,:),allocatable,save :: OSR_nu ! Outgoing SW radiation in each band (Normalized to the band width (W/m2/cm-1)). |
---|
251 | real,dimension(:,:),allocatable,save :: zdtlw ! LW heating tendencies (K/s). |
---|
252 | real,dimension(:,:),allocatable,save :: zdtsw ! SW heating tendencies (K/s). |
---|
253 | real,dimension(:),allocatable,save :: sensibFlux ! Turbulent flux given by the atmosphere to the surface (W.m-2). |
---|
254 | |
---|
255 | !$OMP THREADPRIVATE(fluxsurf_lw,fluxsurf_sw,fluxsurfabs_sw,fluxtop_lw,fluxabs_sw,fluxtop_dn,fluxdyn,OLR_nu,OSR_nu,& |
---|
256 | |
---|
257 | !$OMP zdtlw,zdtsw,sensibFlux) |
---|
258 | |
---|
259 | real zls ! Solar longitude (radians). |
---|
260 | real zlss ! Sub solar point longitude (radians). |
---|
261 | real zday ! Date (time since Ls=0, calculated in sols). |
---|
262 | real zzlay(ngrid,nlayer) ! Altitude at the middle of the atmospheric layers (ref : local surf). |
---|
263 | real zzlev(ngrid,nlayer+1) ! Altitude at the atmospheric layer boundaries (ref : local surf). |
---|
264 | real zzlay_eff(ngrid,nlayer) ! Effective altitude at the middle of the atmospheric layers (ref : geoid ). |
---|
265 | real zzlev_eff(ngrid,nlayer+1) ! Effective altitude at the atmospheric layer boundaries ( ref : geoid ). |
---|
266 | |
---|
267 | ! TENDENCIES due to various processes : |
---|
268 | |
---|
269 | ! For Surface Temperature : (K/s) |
---|
270 | real zdtsurf(ngrid) ! Cumulated tendencies. |
---|
271 | real zdtsurfmr(ngrid) ! Mass_redistribution routine. |
---|
272 | real zdtsdif(ngrid) ! Turbdiff/vdifc routines. |
---|
273 | |
---|
274 | ! For Atmospheric Temperatures : (K/s) |
---|
275 | real zdtdif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
---|
276 | real zdtmr(ngrid,nlayer) ! Mass_redistribution routine. |
---|
277 | real zdtsw1(ngrid,nlayer), zdtlw1(ngrid,nlayer) ! Callcorrk routine. |
---|
278 | |
---|
279 | ! For Surface Tracers : (kg/m2/s) |
---|
280 | real dqsurf(ngrid,nq) ! Cumulated tendencies. |
---|
281 | real zdqsdif(ngrid,nq) ! Turbdiff/vdifc routines. |
---|
282 | real zdqsurfmr(ngrid,nq) ! Mass_redistribution routine. |
---|
283 | |
---|
284 | ! For Tracers : (kg/kg_of_air/s) |
---|
285 | real zdqadj(ngrid,nlayer,nq) ! Convadj routine. |
---|
286 | real zdqdif(ngrid,nlayer,nq) ! Turbdiff/vdifc routines. |
---|
287 | real zdqevap(ngrid,nlayer) ! Turbdiff routine. |
---|
288 | real zdqmr(ngrid,nlayer,nq) ! Mass_redistribution routine. |
---|
289 | |
---|
290 | real zdqchi(ngrid,nlayer,nq) ! Chemical tendency ( chemistry routine ). |
---|
291 | |
---|
292 | real zdqmufi(ngrid,nlayer,nq) ! Microphysical tendency. |
---|
293 | |
---|
294 | real zdqfibar(ngrid,nlayer,nq) ! For 2D chemistry |
---|
295 | real zdqmufibar(ngrid,nlayer,nq) ! For 2D chemistry |
---|
296 | |
---|
297 | ! For Winds : (m/s/s) |
---|
298 | real zdvadj(ngrid,nlayer),zduadj(ngrid,nlayer) ! Convadj routine. |
---|
299 | real zdumr(ngrid,nlayer),zdvmr(ngrid,nlayer) ! Mass_redistribution routine. |
---|
300 | real zdvdif(ngrid,nlayer),zdudif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
---|
301 | real zdhdif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
---|
302 | real zdhadj(ngrid,nlayer) ! Convadj routine. |
---|
303 | |
---|
304 | ! For Pressure and Mass : |
---|
305 | real zdmassmr(ngrid,nlayer) ! Atmospheric Mass tendency for mass_redistribution (kg_of_air/m2/s). |
---|
306 | real zdmassmr_col(ngrid) ! Atmospheric Column Mass tendency for mass_redistribution (kg_of_air/m2/s). |
---|
307 | real zdpsrfmr(ngrid) ! Pressure tendency for mass_redistribution routine (Pa/s). |
---|
308 | |
---|
309 | |
---|
310 | |
---|
311 | ! Local variables for LOCAL CALCULATIONS: |
---|
312 | ! --------------------------------------- |
---|
313 | real zflubid(ngrid) |
---|
314 | real zplanck(ngrid),zpopsk(ngrid,nlayer) |
---|
315 | real ztim1,ztim2,ztim3, z1,z2 |
---|
316 | real ztime_fin |
---|
317 | real zdh(ngrid,nlayer) |
---|
318 | real gmplanet |
---|
319 | real taux(ngrid),tauy(ngrid) |
---|
320 | |
---|
321 | |
---|
322 | |
---|
323 | ! local variables for DIAGNOSTICS : (diagfi & stat) |
---|
324 | ! ------------------------------------------------- |
---|
325 | real ps(ngrid) ! Surface Pressure. |
---|
326 | real zt(ngrid,nlayer) ! Atmospheric Temperature. |
---|
327 | real zu(ngrid,nlayer),zv(ngrid,nlayer) ! Zonal and Meridional Winds. |
---|
328 | real zq(ngrid,nlayer,nq) ! Atmospheric Tracers. |
---|
329 | real zdtadj(ngrid,nlayer) ! Convadj Diagnostic. |
---|
330 | real zdtdyn(ngrid,nlayer) ! Dynamical Heating (K/s). |
---|
331 | real zdudyn(ngrid,nlayer) ! Dynamical Zonal Wind tendency (m.s-2). |
---|
332 | real,allocatable,dimension(:,:),save :: ztprevious ! Previous loop Atmospheric Temperature (K) ! Useful for Dynamical Heating calculation. |
---|
333 | real,allocatable,dimension(:,:),save :: zuprevious ! Previous loop Zonal Wind (m.s-1) ! Useful for Zonal Wind tendency calculation. |
---|
334 | !$OMP THREADPRIVATE(ztprevious,zuprevious) |
---|
335 | |
---|
336 | real zhorizwind(ngrid,nlayer) ! Horizontal Wind ( sqrt(u**+v*v)) |
---|
337 | |
---|
338 | real vmr(ngrid,nlayer) ! volume mixing ratio |
---|
339 | real time_phys |
---|
340 | |
---|
341 | real ISR,ASR,OLR,GND,DYN,GSR,Ts1,Ts2,Ts3,TsS ! for Diagnostic. |
---|
342 | |
---|
343 | ! to test energy conservation (RW+JL) |
---|
344 | real mass(ngrid,nlayer),massarea(ngrid,nlayer) |
---|
345 | real dEtot, dEtots, AtmToSurf_TurbFlux |
---|
346 | real,save :: dEtotSW, dEtotsSW, dEtotLW, dEtotsLW |
---|
347 | !$OMP THREADPRIVATE(dEtotSW, dEtotsSW, dEtotLW, dEtotsLW) |
---|
348 | real dEzRadsw(ngrid,nlayer),dEzRadlw(ngrid,nlayer),dEzdiff(ngrid,nlayer) |
---|
349 | real dEdiffs(ngrid),dEdiff(ngrid) |
---|
350 | |
---|
351 | !JL12 conservation test for mean flow kinetic energy has been disabled temporarily |
---|
352 | |
---|
353 | real dItot, dItot_tmp, dVtot, dVtot_tmp |
---|
354 | |
---|
355 | real dWtot, dWtot_tmp, dWtots, dWtots_tmp |
---|
356 | |
---|
357 | |
---|
358 | ! For Clear Sky Case. |
---|
359 | real fluxsurf_lw1(ngrid), fluxsurf_sw1(ngrid), fluxsurfabs_sw1(ngrid) ! For SW/LW flux. |
---|
360 | real fluxtop_lw1(ngrid), fluxabs_sw1(ngrid) ! For SW/LW flux. |
---|
361 | real albedo_equivalent1(ngrid) ! For Equivalent albedo calculation. |
---|
362 | real tf, ntf |
---|
363 | |
---|
364 | real,allocatable,dimension(:,:),save :: qsurf_hist |
---|
365 | !$OMP THREADPRIVATE(qsurf_hist) |
---|
366 | |
---|
367 | ! Miscellaneous : |
---|
368 | character(len=10) :: tmp1 |
---|
369 | character(len=10) :: tmp2 |
---|
370 | |
---|
371 | ! Local variables for Titan chemistry and microphysics |
---|
372 | ! ---------------------------------------------------- |
---|
373 | |
---|
374 | real,save :: ctimestep ! Chemistry timestep (s) |
---|
375 | !$OMP THREADPRIVATE(ctimestep) |
---|
376 | |
---|
377 | ! Chemical tracers in molar fraction |
---|
378 | real, dimension(ngrid,nlayer,nkim) :: ychim ! (mol/mol) |
---|
379 | real, dimension(ngrid,nlayer,nkim) :: ychimbar ! For 2D chemistry |
---|
380 | |
---|
381 | ! Molar fraction tendencies ( chemistry, condensation and evaporation ) for tracers (mol/mol/s) |
---|
382 | real, dimension(:,:,:), allocatable, save :: dycchi ! NB : Only for chem tracers. Saved since chemistry is not called every step. |
---|
383 | !$OMP THREADPRIVATE(dycchi) |
---|
384 | real, dimension(ngrid,nlayer,nq) :: dyccond ! Condensation rate. NB : for all tracers, as we want to use indx on it. |
---|
385 | real, dimension(ngrid,nlayer,nq) :: dyccondbar ! For 2D chemistry |
---|
386 | real, dimension(ngrid) :: dycevapCH4 ! Surface "pseudo-evaporation" rate (forcing constant surface humidity). |
---|
387 | |
---|
388 | ! Saturation profiles |
---|
389 | real, dimension(ngrid,nlayer,nkim) :: ysat ! (mol/mol) |
---|
390 | |
---|
391 | ! Surface methane |
---|
392 | real, dimension(:), allocatable, save :: tankCH4 ! Depth of surface methane tank (m) |
---|
393 | !$OMP THREADPRIVATE(tankCH4) |
---|
394 | |
---|
395 | real :: i2e(ngrid,nlayer) ! int 2 ext factor ( X.kg-1 -> X.m-3 for diags ) |
---|
396 | |
---|
397 | real,save,dimension(:,:,:), allocatable :: tpq ! Tracers for decoupled microphysical tests ( temporary in 01/18 ) |
---|
398 | !$OMP THREADPRIVATE(tpq) |
---|
399 | real,dimension(ngrid,nlayer,nq) :: dtpq ! (temporary in 01/18) |
---|
400 | |
---|
401 | |
---|
402 | !----------------------------------------------------------------------------- |
---|
403 | ! Interface to calmufi |
---|
404 | ! --> needed in order to pass assumed-shape arrays. Otherwise we must put calmufi in a module |
---|
405 | ! (to have an explicit interface generated by the compiler). |
---|
406 | ! Or one can put calmufi in MMP_GCM module (in muphytitan). |
---|
407 | INTERFACE |
---|
408 | SUBROUTINE calmufi(dt, plev, zlev, play, zlay, g3d, temp, pq, zdqfi, zdq) |
---|
409 | REAL(kind=8), INTENT(IN) :: dt !! Physics timestep (s). |
---|
410 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: plev !! Pressure levels (Pa). |
---|
411 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlev !! Altitude levels (m). |
---|
412 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: play !! Pressure layers (Pa). |
---|
413 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlay !! Altitude at the center of each layer (m). |
---|
414 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: g3d !! Latitude-Altitude depending gravitational acceleration (m.s-2). |
---|
415 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: temp !! Temperature at the center of each layer (K). |
---|
416 | REAL(kind=8), DIMENSION(:,:,:), INTENT(IN) :: pq !! Tracers (\(kg.kg^{-1}}\)). |
---|
417 | REAL(kind=8), DIMENSION(:,:,:), INTENT(IN) :: zdqfi !! Tendency from former processes for tracers (\(kg.kg^{-1}}\)). |
---|
418 | REAL(kind=8), DIMENSION(:,:,:), INTENT(OUT) :: zdq !! Microphysical tendency for tracers (\(kg.kg^{-1}}\)). |
---|
419 | END SUBROUTINE calmufi |
---|
420 | END INTERFACE |
---|
421 | |
---|
422 | !================================================================================================== |
---|
423 | |
---|
424 | ! ----------------- |
---|
425 | ! I. INITIALISATION |
---|
426 | ! ----------------- |
---|
427 | |
---|
428 | ! -------------------------------- |
---|
429 | ! I.1 First Call Initialisation. |
---|
430 | ! -------------------------------- |
---|
431 | if (firstcall) then |
---|
432 | allocate(tpq(ngrid,nlayer,nq)) |
---|
433 | tpq(:,:,:) = pq(:,:,:) |
---|
434 | |
---|
435 | ! Initialisation of nmicro as well as tracers names, indexes ... |
---|
436 | if (ngrid.ne.1) then ! Already done in rcm1d |
---|
437 | call initracer2(nq,nametrac) ! WARNING JB (27/03/2018): should be wrapped in an OMP SINGLE statement (see module notes) |
---|
438 | endif |
---|
439 | |
---|
440 | ! Allocate saved arrays. |
---|
441 | ALLOCATE(tsurf(ngrid)) |
---|
442 | ALLOCATE(tsoil(ngrid,nsoilmx)) |
---|
443 | ALLOCATE(albedo(ngrid,L_NSPECTV)) |
---|
444 | ALLOCATE(albedo_equivalent(ngrid)) |
---|
445 | ALLOCATE(albedo_bareground(ngrid)) |
---|
446 | ALLOCATE(emis(ngrid)) |
---|
447 | ALLOCATE(dtrad(ngrid,nlayer)) |
---|
448 | ALLOCATE(fluxrad_sky(ngrid)) |
---|
449 | ALLOCATE(fluxrad(ngrid)) |
---|
450 | ALLOCATE(capcal(ngrid)) |
---|
451 | ALLOCATE(fluxgrd(ngrid)) |
---|
452 | ALLOCATE(qsurf(ngrid,nq)) |
---|
453 | ALLOCATE(q2(ngrid,nlayer+1)) |
---|
454 | ALLOCATE(tankCH4(ngrid)) |
---|
455 | ALLOCATE(ztprevious(ngrid,nlayer)) |
---|
456 | ALLOCATE(zuprevious(ngrid,nlayer)) |
---|
457 | ALLOCATE(qsurf_hist(ngrid,nq)) |
---|
458 | ALLOCATE(fluxsurf_lw(ngrid)) |
---|
459 | ALLOCATE(fluxsurf_sw(ngrid)) |
---|
460 | ALLOCATE(fluxsurfabs_sw(ngrid)) |
---|
461 | ALLOCATE(fluxtop_lw(ngrid)) |
---|
462 | ALLOCATE(fluxabs_sw(ngrid)) |
---|
463 | ALLOCATE(fluxtop_dn(ngrid)) |
---|
464 | ALLOCATE(fluxdyn(ngrid)) |
---|
465 | ALLOCATE(OLR_nu(ngrid,L_NSPECTI)) |
---|
466 | ALLOCATE(OSR_nu(ngrid,L_NSPECTV)) |
---|
467 | ALLOCATE(sensibFlux(ngrid)) |
---|
468 | ALLOCATE(zdtlw(ngrid,nlayer)) |
---|
469 | ALLOCATE(zdtsw(ngrid,nlayer)) |
---|
470 | |
---|
471 | ! This is defined in comsaison_h |
---|
472 | ALLOCATE(mu0(ngrid)) |
---|
473 | ALLOCATE(fract(ngrid)) |
---|
474 | ! This is defined in radcommon_h |
---|
475 | ALLOCATE(gzlat(ngrid,nlayer)) |
---|
476 | ALLOCATE(gzlat_ig(nlayer)) |
---|
477 | ALLOCATE(Cmk(nlayer)) |
---|
478 | |
---|
479 | ! Variables set to 0 |
---|
480 | ! ~~~~~~~~~~~~~~~~~~ |
---|
481 | dtrad(:,:) = 0.0 |
---|
482 | fluxrad(:) = 0.0 |
---|
483 | zdtsw(:,:) = 0.0 |
---|
484 | zdtlw(:,:) = 0.0 |
---|
485 | |
---|
486 | ! Initialize setup for correlated-k radiative transfer |
---|
487 | ! JVO 17 : Was in callcorrk firstcall, but we need spectral intervals for microphysics |
---|
488 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
489 | |
---|
490 | if (corrk) then |
---|
491 | |
---|
492 | call system('rm -f surf_vals_long.out') |
---|
493 | |
---|
494 | write( tmp1, '(i3)' ) L_NSPECTI |
---|
495 | write( tmp2, '(i3)' ) L_NSPECTV |
---|
496 | banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)) |
---|
497 | banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)) |
---|
498 | |
---|
499 | call setspi ! Basic infrared properties. |
---|
500 | call setspv ! Basic visible properties. |
---|
501 | call sugas_corrk ! Set up gaseous absorption properties. |
---|
502 | |
---|
503 | OLR_nu(:,:) = 0. |
---|
504 | OSR_nu(:,:) = 0. |
---|
505 | |
---|
506 | endif |
---|
507 | |
---|
508 | ! Initialize names and timestep for chemistry |
---|
509 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
510 | |
---|
511 | if ( callchim ) then |
---|
512 | |
---|
513 | if ( moyzon_ch .and. ngrid.eq.1 ) then |
---|
514 | print *, "moyzon_ch=",moyzon_ch," and ngrid=1" |
---|
515 | print *, "Please desactivate zonal mean for 1D !" |
---|
516 | stop |
---|
517 | endif |
---|
518 | |
---|
519 | allocate(dycchi(ngrid,nlayer,nkim)) ! only for chemical tracers |
---|
520 | |
---|
521 | ! Chemistry timestep |
---|
522 | ctimestep = ptimestep*REAL(ichim) |
---|
523 | |
---|
524 | zdqchi(:,:,:) = 0.0 |
---|
525 | |
---|
526 | endif |
---|
527 | |
---|
528 | ! Initialize microphysics. |
---|
529 | ! ~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
530 | |
---|
531 | IF ( callmufi ) THEN |
---|
532 | ! WARNING JB (27/03/2018): inimufi AND mmp_initialize_optics should be wrapped in an OMP SINGLE statement. |
---|
533 | call inimufi(ptimestep) |
---|
534 | ! Optical coupling of YAMMS is plugged but inactivated for now |
---|
535 | ! as long as the microphysics only isn't fully debugged -- JVO 01/18 |
---|
536 | ! NOTE JB (03/18): mmp_optic_file is initialized in inimufi => mmp_initialize_optics must be called AFTER inimufi |
---|
537 | IF (.NOT.uncoupl_optic_haze) call mmp_initialize_optics(mmp_optic_file) |
---|
538 | |
---|
539 | ! initialize microphysics diagnostics arrays. |
---|
540 | call ini_diag_arrays(ngrid,nlayer,nice) |
---|
541 | |
---|
542 | ENDIF |
---|
543 | |
---|
544 | |
---|
545 | #ifdef CPP_XIOS |
---|
546 | ! Initialize XIOS context |
---|
547 | write(*,*) "physiq: call wxios_context_init" |
---|
548 | CALL wxios_context_init |
---|
549 | #endif |
---|
550 | |
---|
551 | ! Read 'startfi.nc' file. |
---|
552 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
553 | call phyetat0(startphy_file,ngrid,nlayer,"startfi.nc",0,0,nsoilmx,nq, & |
---|
554 | day_ini,time_phys,tsurf,tsoil,emis,q2,qsurf,tankCH4) |
---|
555 | if (.not.startphy_file) then |
---|
556 | ! additionnal "academic" initialization of physics |
---|
557 | if (is_master) write(*,*) "Physiq: initializing tsurf(:) to pt(:,1) !!" |
---|
558 | tsurf(:)=pt(:,1) |
---|
559 | if (is_master) write(*,*) "Physiq: initializing tsoil(:) to pt(:,1) !!" |
---|
560 | do isoil=1,nsoilmx |
---|
561 | tsoil(:,isoil)=tsurf(:) |
---|
562 | enddo |
---|
563 | if (is_master) write(*,*) "Physiq: initializing day_ini to pdat !" |
---|
564 | day_ini=pday |
---|
565 | endif |
---|
566 | |
---|
567 | if (pday.ne.day_ini) then |
---|
568 | write(*,*) "ERROR in physiq.F90:" |
---|
569 | write(*,*) "bad synchronization between physics and dynamics" |
---|
570 | write(*,*) "dynamics day: ",pday |
---|
571 | write(*,*) "physics day: ",day_ini |
---|
572 | stop |
---|
573 | endif |
---|
574 | write (*,*) 'In physiq day_ini =', day_ini |
---|
575 | |
---|
576 | ! Initialize albedo calculation. |
---|
577 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
578 | albedo(:,:)=0.0 |
---|
579 | albedo_bareground(:)=0.0 |
---|
580 | call surfini(ngrid,nq,qsurf,albedo,albedo_bareground) |
---|
581 | |
---|
582 | ! Initialize orbital calculation. |
---|
583 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
584 | call iniorbit(apoastr,periastr,year_day,peri_day,obliquit) |
---|
585 | |
---|
586 | if(tlocked)then |
---|
587 | print*,'Planet is tidally locked at resonance n=',nres |
---|
588 | print*,'Make sure you have the right rotation rate!!!' |
---|
589 | endif |
---|
590 | |
---|
591 | |
---|
592 | ! Initialize soil. |
---|
593 | ! ~~~~~~~~~~~~~~~~ |
---|
594 | if (callsoil) then |
---|
595 | |
---|
596 | call soil(ngrid,nsoilmx,firstcall,lastcall,inertiedat, & |
---|
597 | ptimestep,tsurf,tsoil,capcal,fluxgrd) |
---|
598 | |
---|
599 | else ! else of 'callsoil'. |
---|
600 | |
---|
601 | print*,'WARNING! Thermal conduction in the soil turned off' |
---|
602 | capcal(:)=1.e6 |
---|
603 | fluxgrd(:)=intheat |
---|
604 | print*,'Flux from ground = ',intheat,' W m^-2' |
---|
605 | |
---|
606 | endif ! end of 'callsoil'. |
---|
607 | |
---|
608 | icount=1 |
---|
609 | |
---|
610 | |
---|
611 | ! Initialize surface history variable. |
---|
612 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
613 | qsurf_hist(:,:)=qsurf(:,:) |
---|
614 | |
---|
615 | ! Initialize variable for dynamical heating and zonal wind tendency diagnostic |
---|
616 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
617 | ztprevious(:,:)=pt(:,:) |
---|
618 | zuprevious(:,:)=pu(:,:) |
---|
619 | |
---|
620 | |
---|
621 | if(meanOLR)then |
---|
622 | call system('rm -f rad_bal.out') ! to record global radiative balance. |
---|
623 | call system('rm -f tem_bal.out') ! to record global mean/max/min temperatures. |
---|
624 | call system('rm -f h2o_bal.out') ! to record global hydrological balance. |
---|
625 | endif |
---|
626 | |
---|
627 | if (ngrid.ne.1) then |
---|
628 | ! Note : no need to create a restart file in 1d. |
---|
629 | call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid,nlayer,nq, & |
---|
630 | ptimestep,pday+nday,time_phys,cell_area, & |
---|
631 | albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe) |
---|
632 | endif |
---|
633 | |
---|
634 | ! XIOS outputs |
---|
635 | #ifdef CPP_XIOS |
---|
636 | |
---|
637 | write(*,*) "physiq: call initialize_xios_output" |
---|
638 | call initialize_xios_output(pday,ptime,ptimestep,daysec, & |
---|
639 | presnivs,pseudoalt) |
---|
640 | #endif |
---|
641 | write(*,*) "physiq: end of firstcall" |
---|
642 | endif ! end of 'firstcall' |
---|
643 | |
---|
644 | ! ------------------------------------------------------ |
---|
645 | ! I.2 Initializations done at every physical timestep: |
---|
646 | ! ------------------------------------------------------ |
---|
647 | |
---|
648 | #ifdef CPP_XIOS |
---|
649 | ! update XIOS time/calendar |
---|
650 | call update_xios_timestep |
---|
651 | #endif |
---|
652 | |
---|
653 | ! Initialize various variables |
---|
654 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
655 | |
---|
656 | pdt(:,:) = 0.0 |
---|
657 | zdtsurf(:) = 0.0 |
---|
658 | pdq(:,:,:) = 0.0 |
---|
659 | dqsurf(:,:) = 0.0 |
---|
660 | pdu(:,:) = 0.0 |
---|
661 | pdv(:,:) = 0.0 |
---|
662 | pdpsrf(:) = 0.0 |
---|
663 | zflubid(:) = 0.0 |
---|
664 | taux(:) = 0.0 |
---|
665 | tauy(:) = 0.0 |
---|
666 | |
---|
667 | zday=pday+ptime ! Compute time, in sols (and fraction thereof). |
---|
668 | |
---|
669 | ! Compute Stellar Longitude (Ls), and orbital parameters. |
---|
670 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
671 | if (season) then |
---|
672 | call stellarlong(zday,zls) |
---|
673 | else |
---|
674 | call stellarlong(float(day_ini),zls) |
---|
675 | end if |
---|
676 | |
---|
677 | call orbite(zls,dist_star,declin,right_ascen) |
---|
678 | |
---|
679 | if (tlocked) then |
---|
680 | zlss=Mod(-(2.*pi*(zday/year_day)*nres - right_ascen),2.*pi) |
---|
681 | elseif (diurnal) then |
---|
682 | zlss=-2.*pi*(zday-.5) |
---|
683 | else if(diurnal .eqv. .false.) then |
---|
684 | zlss=9999. |
---|
685 | endif |
---|
686 | |
---|
687 | |
---|
688 | ! Compute variations of g with latitude (oblate case). |
---|
689 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
690 | if (oblate .eqv. .false.) then |
---|
691 | gzlat(:,:) = g |
---|
692 | else if (flatten .eq. 0.0 .or. J2 .eq. 0.0 .or. Rmean .eq. 0.0 .or. MassPlanet .eq. 0.0) then |
---|
693 | print*,'I need values for flatten, J2, Rmean and MassPlanet to compute gzlat (else set oblate=.false.)' |
---|
694 | call abort |
---|
695 | else |
---|
696 | gmplanet = MassPlanet*grav*1e24 |
---|
697 | do ig=1,ngrid |
---|
698 | gzlat(ig,:)= gmplanet/(Rmean**2) * (1.D0 + 0.75 *J2 - 2.0*flatten/3. + (2.*flatten - 15./4.* J2) * cos(2. * (pi/2. - latitude(ig)))) |
---|
699 | end do |
---|
700 | endif |
---|
701 | |
---|
702 | ! Compute altitudes with the geopotential coming from the dynamics. |
---|
703 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
704 | |
---|
705 | if (eff_gz .eqv. .false.) then |
---|
706 | |
---|
707 | do l=1,nlayer |
---|
708 | zzlay(:,l) = pphi(:,l) / gzlat(:,l) ! Reference = local surface |
---|
709 | enddo |
---|
710 | |
---|
711 | else ! In this case we consider variations of g with altitude |
---|
712 | |
---|
713 | do l=1,nlayer |
---|
714 | zzlay(:,l) = g*rad*rad / ( g*rad - ( pphi(:,l) + phisfi(:) )) - rad |
---|
715 | gzlat(:,l) = g*rad*rad / ( rad + zzlay(:,l) )**2 |
---|
716 | end do |
---|
717 | |
---|
718 | endif ! if eff_gz |
---|
719 | |
---|
720 | zzlev(:,1)=0. |
---|
721 | zzlev(:,nlayer+1)=1.e7 ! Dummy top of last layer above 10000 km... |
---|
722 | ! JVO 19 : This altitude is indeed dummy for the GCM and fits ptop=0 |
---|
723 | ! but for upper chemistry that's a pb -> we anyway redefine it just after .. |
---|
724 | |
---|
725 | do l=2,nlayer |
---|
726 | do ig=1,ngrid |
---|
727 | z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l)) |
---|
728 | z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l)) |
---|
729 | zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2) |
---|
730 | enddo |
---|
731 | enddo |
---|
732 | |
---|
733 | ! Effective altitudes ( eg needed for chemistry ) with correct g, and with reference to the geoid |
---|
734 | ! JVO 19 : We shall always have correct altitudes in chemistry no matter what's in physics |
---|
735 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
736 | if (moyzon_ch) then ! Zonal averages |
---|
737 | |
---|
738 | zzlaybar(1,:)=g*rad*rad/(g*rad-(zphibar(1,:)+zphisbar(1)))-rad ! reference = geoid |
---|
739 | zzlevbar(1,1)=zphisbar(1)/g |
---|
740 | DO l=2,nlayer |
---|
741 | z1=(zplaybar(1,l-1)+zplevbar(1,l))/(zplaybar(1,l-1)-zplevbar(1,l)) |
---|
742 | z2=(zplevbar(1,l) +zplaybar(1,l))/(zplevbar(1,l) -zplaybar(1,l)) |
---|
743 | zzlevbar(1,l)=(z1*zzlaybar(1,l-1)+z2*zzlaybar(1,l))/(z1+z2) |
---|
744 | ENDDO |
---|
745 | zzlevbar(1,nlayer+1)=zzlaybar(1,nlayer)+(zzlaybar(1,nlayer)-zzlevbar(1,nlayer)) |
---|
746 | |
---|
747 | DO ig=2,ngrid |
---|
748 | if (latitude(ig).ne.latitude(ig-1)) then |
---|
749 | DO l=1,nlayer |
---|
750 | zzlaybar(ig,l)=g*rad*rad/(g*rad-(zphibar(ig,l)+zphisbar(ig)))-rad |
---|
751 | ENDDO |
---|
752 | zzlevbar(ig,1)=zphisbar(ig)/g |
---|
753 | DO l=2,nlayer |
---|
754 | z1=(zplaybar(ig,l-1)+zplevbar(ig,l))/ (zplaybar(ig,l-1)-zplevbar(ig,l)) |
---|
755 | z2=(zplevbar(ig,l) +zplaybar(ig,l))/(zplevbar(ig,l) -zplaybar(ig,l)) |
---|
756 | zzlevbar(ig,l)=(z1*zzlaybar(ig,l-1)+z2*zzlaybar(ig,l))/(z1+z2) |
---|
757 | ENDDO |
---|
758 | zzlevbar(ig,nlayer+1)=zzlaybar(ig,nlayer)+(zzlaybar(ig,nlayer)-zzlevbar(ig,nlayer)) |
---|
759 | else |
---|
760 | zzlaybar(ig,:)=zzlaybar(ig-1,:) |
---|
761 | zzlevbar(ig,:)=zzlevbar(ig-1,:) |
---|
762 | endif |
---|
763 | ENDDO |
---|
764 | |
---|
765 | else ! if not moyzon |
---|
766 | |
---|
767 | DO l=1,nlayer |
---|
768 | zzlay_eff(ig,l)=g*rad*rad/(g*rad-(pphi(ig,l)+phisfi(ig)))-rad ! reference = geoid |
---|
769 | ENDDO |
---|
770 | zzlev_eff(ig,1)=phisfi(ig)/g |
---|
771 | DO l=2,nlayer |
---|
772 | z1=(pplay(ig,l-1)+pplev(ig,l))/ (pplay(ig,l-1)-pplev(ig,l)) |
---|
773 | z2=(pplev(ig,l) +pplay(ig,l))/(pplev(ig,l) -pplay(ig,l)) |
---|
774 | zzlev_eff(ig,l)=(z1*zzlay_eff(ig,l-1)+z2*zzlay_eff(ig,l))/(z1+z2) |
---|
775 | ENDDO |
---|
776 | zzlev_eff(ig,nlayer+1)=zzlay_eff(ig,nlayer)+(zzlay_eff(ig,nlayer)-zzlev_eff(ig,nlayer)) |
---|
777 | |
---|
778 | endif ! moyzon |
---|
779 | |
---|
780 | ! ------------------------------------------------------------------------------------- |
---|
781 | ! Compute potential temperature |
---|
782 | ! Note : Potential temperature calculation may not be the same in physiq and dynamic... |
---|
783 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
784 | do l=1,nlayer |
---|
785 | do ig=1,ngrid |
---|
786 | zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp |
---|
787 | zh(ig,l)=pt(ig,l)/zpopsk(ig,l) |
---|
788 | mass(ig,l) = (pplev(ig,l) - pplev(ig,l+1))/gzlat(ig,l) |
---|
789 | massarea(ig,l)=mass(ig,l)*cell_area(ig) |
---|
790 | enddo |
---|
791 | enddo |
---|
792 | |
---|
793 | ! Compute vertical velocity (m/s) from vertical mass flux |
---|
794 | ! w = F / (rho*area) and rho = P/(r*T) |
---|
795 | ! But first linearly interpolate mass flux to mid-layers |
---|
796 | do l=1,nlayer-1 |
---|
797 | pw(:,l)=0.5*(flxw(:,l)+flxw(:,l+1)) |
---|
798 | enddo |
---|
799 | pw(:,nlayer)=0.5*flxw(:,nlayer) ! since flxw(nlayer+1)=0 |
---|
800 | do l=1,nlayer |
---|
801 | pw(:,l)=(pw(:,l)*r*pt(:,l)) / (pplay(:,l)*cell_area(:)) |
---|
802 | enddo |
---|
803 | |
---|
804 | !--------------------------------- |
---|
805 | ! II. Compute radiative tendencies |
---|
806 | !--------------------------------- |
---|
807 | |
---|
808 | if (callrad) then |
---|
809 | if( mod(icount-1,iradia).eq.0.or.lastcall) then |
---|
810 | |
---|
811 | ! Compute local stellar zenith angles |
---|
812 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
813 | if (tlocked) then |
---|
814 | ! JL14 corrects tidally resonant (and inclined) cases. nres=omega_rot/omega_orb |
---|
815 | ztim1=SIN(declin) |
---|
816 | ztim2=COS(declin)*COS(zlss) |
---|
817 | ztim3=COS(declin)*SIN(zlss) |
---|
818 | |
---|
819 | call stelang(ngrid,sinlon,coslon,sinlat,coslat, & |
---|
820 | ztim1,ztim2,ztim3,mu0,fract, flatten) |
---|
821 | |
---|
822 | elseif (diurnal) then |
---|
823 | ztim1=SIN(declin) |
---|
824 | ztim2=COS(declin)*COS(2.*pi*(zday-.5)) |
---|
825 | ztim3=-COS(declin)*SIN(2.*pi*(zday-.5)) |
---|
826 | |
---|
827 | call stelang(ngrid,sinlon,coslon,sinlat,coslat, & |
---|
828 | ztim1,ztim2,ztim3,mu0,fract, flatten) |
---|
829 | else if(diurnal .eqv. .false.) then |
---|
830 | |
---|
831 | call mucorr(ngrid,declin,latitude,mu0,fract,10000.,rad,flatten) |
---|
832 | ! WARNING: this function appears not to work in 1D |
---|
833 | |
---|
834 | endif |
---|
835 | |
---|
836 | ! Eclipse incoming sunlight (e.g. Saturn ring shadowing). |
---|
837 | if(rings_shadow) then |
---|
838 | call call_rings(ngrid, ptime, pday, diurnal) |
---|
839 | endif |
---|
840 | |
---|
841 | |
---|
842 | if (corrk) then |
---|
843 | |
---|
844 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
845 | ! II.a Call correlated-k radiative transfer scheme |
---|
846 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
847 | |
---|
848 | call call_profilgases(nlayer) |
---|
849 | |
---|
850 | ! standard callcorrk |
---|
851 | call callcorrk(ngrid,nlayer,pq,nq,qsurf,zday, & |
---|
852 | albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt, & |
---|
853 | tsurf,fract,dist_star, & |
---|
854 | zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw, & |
---|
855 | fluxsurfabs_sw,fluxtop_lw, & |
---|
856 | fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu, & |
---|
857 | lastcall) |
---|
858 | |
---|
859 | ! Radiative flux from the sky absorbed by the surface (W.m-2). |
---|
860 | GSR=0.0 |
---|
861 | fluxrad_sky(:)=emis(:)*fluxsurf_lw(:)+fluxsurfabs_sw(:) |
---|
862 | |
---|
863 | !if(noradsurf)then ! no lower surface; SW flux just disappears |
---|
864 | ! GSR = SUM(fluxsurf_sw(:)*cell_area(:))/totarea |
---|
865 | ! fluxrad_sky(:)=emis(:)*fluxsurf_lw(:) |
---|
866 | ! print*,'SW lost in deep atmosphere = ',GSR,' W m^-2' |
---|
867 | !endif |
---|
868 | |
---|
869 | ! Net atmospheric radiative heating rate (K.s-1) |
---|
870 | dtrad(:,:)=zdtsw(:,:)+zdtlw(:,:) |
---|
871 | |
---|
872 | elseif(newtonian)then |
---|
873 | |
---|
874 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
875 | ! II.b Call Newtonian cooling scheme |
---|
876 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
877 | call newtrelax(ngrid,nlayer,mu0,sinlat,zpopsk,pt,pplay,pplev,dtrad,firstcall) |
---|
878 | |
---|
879 | zdtsurf(:) = +(pt(:,1)-tsurf(:))/ptimestep |
---|
880 | ! e.g. surface becomes proxy for 1st atmospheric layer ? |
---|
881 | |
---|
882 | else |
---|
883 | |
---|
884 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
885 | ! II.c Atmosphere has no radiative effect |
---|
886 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
887 | fluxtop_dn(:) = fract(:)*mu0(:)*Fat1AU/dist_star**2 |
---|
888 | if(ngrid.eq.1)then ! / by 4 globally in 1D case... |
---|
889 | fluxtop_dn(1) = fract(1)*Fat1AU/dist_star**2/2.0 |
---|
890 | endif |
---|
891 | fluxsurf_sw(:) = fluxtop_dn(:) |
---|
892 | print*,'------------WARNING---WARNING------------' ! by MT2015. |
---|
893 | print*,'You are in corrk=false mode, ' |
---|
894 | print*,'and the surface albedo is taken equal to the first visible spectral value' |
---|
895 | |
---|
896 | fluxsurfabs_sw(:) = fluxtop_dn(:)*(1.-albedo(:,1)) |
---|
897 | fluxrad_sky(:) = fluxsurfabs_sw(:) |
---|
898 | fluxtop_lw(:) = emis(:)*sigma*tsurf(:)**4 |
---|
899 | |
---|
900 | dtrad(:,:)=0.0 ! no atmospheric radiative heating |
---|
901 | |
---|
902 | endif ! end of corrk |
---|
903 | |
---|
904 | endif ! of if(mod(icount-1,iradia).eq.0) |
---|
905 | |
---|
906 | |
---|
907 | ! Transformation of the radiative tendencies |
---|
908 | ! ------------------------------------------ |
---|
909 | zplanck(:)=tsurf(:)*tsurf(:) |
---|
910 | zplanck(:)=emis(:)*sigma*zplanck(:)*zplanck(:) |
---|
911 | fluxrad(:)=fluxrad_sky(:)-zplanck(:) |
---|
912 | pdt(:,:)=pdt(:,:)+dtrad(:,:) |
---|
913 | |
---|
914 | ! Test of energy conservation |
---|
915 | !---------------------------- |
---|
916 | if(enertest)then |
---|
917 | call planetwide_sumval(cpp*massarea(:,:)*zdtsw(:,:)/totarea_planet,dEtotSW) |
---|
918 | call planetwide_sumval(cpp*massarea(:,:)*zdtlw(:,:)/totarea_planet,dEtotLW) |
---|
919 | !call planetwide_sumval(fluxsurf_sw(:)*(1.-albedo_equivalent(:))*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk |
---|
920 | call planetwide_sumval(fluxsurfabs_sw(:)*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk |
---|
921 | call planetwide_sumval((fluxsurf_lw(:)*emis(:)-zplanck(:))*cell_area(:)/totarea_planet,dEtotsLW) |
---|
922 | dEzRadsw(:,:)=cpp*mass(:,:)*zdtsw(:,:) |
---|
923 | dEzRadlw(:,:)=cpp*mass(:,:)*zdtlw(:,:) |
---|
924 | if (is_master) then |
---|
925 | print*,'---------------------------------------------------------------' |
---|
926 | print*,'In corrk SW atmospheric heating =',dEtotSW,' W m-2' |
---|
927 | print*,'In corrk LW atmospheric heating =',dEtotLW,' W m-2' |
---|
928 | print*,'atmospheric net rad heating (SW+LW) =',dEtotLW+dEtotSW,' W m-2' |
---|
929 | print*,'In corrk SW surface heating =',dEtotsSW,' W m-2' |
---|
930 | print*,'In corrk LW surface heating =',dEtotsLW,' W m-2' |
---|
931 | print*,'surface net rad heating (SW+LW) =',dEtotsLW+dEtotsSW,' W m-2' |
---|
932 | endif |
---|
933 | endif ! end of 'enertest' |
---|
934 | |
---|
935 | endif ! of if (callrad) |
---|
936 | |
---|
937 | |
---|
938 | |
---|
939 | ! -------------------------------------------- |
---|
940 | ! III. Vertical diffusion (turbulent mixing) : |
---|
941 | ! -------------------------------------------- |
---|
942 | |
---|
943 | if (calldifv) then |
---|
944 | |
---|
945 | zflubid(:)=fluxrad(:)+fluxgrd(:) |
---|
946 | |
---|
947 | ! JL12 the following if test is temporarily there to allow us to compare the old vdifc with turbdiff. |
---|
948 | if (UseTurbDiff) then |
---|
949 | |
---|
950 | call turbdiff(ngrid,nlayer,nq, & |
---|
951 | ptimestep,capcal,lwrite, & |
---|
952 | pplay,pplev,zzlay,zzlev,z0, & |
---|
953 | pu,pv,pt,zpopsk,pq,tsurf,emis,qsurf, & |
---|
954 | pdt,pdq,zflubid, & |
---|
955 | zdudif,zdvdif,zdtdif,zdtsdif, & |
---|
956 | sensibFlux,q2,zdqdif,zdqsdif, & |
---|
957 | taux,tauy,lastcall) |
---|
958 | |
---|
959 | else |
---|
960 | |
---|
961 | zdh(:,:)=pdt(:,:)/zpopsk(:,:) |
---|
962 | |
---|
963 | call vdifc(ngrid,nlayer,nq,zpopsk, & |
---|
964 | ptimestep,capcal,lwrite, & |
---|
965 | pplay,pplev,zzlay,zzlev,z0, & |
---|
966 | pu,pv,zh,pq,tsurf,emis,qsurf, & |
---|
967 | zdh,pdq,zflubid, & |
---|
968 | zdudif,zdvdif,zdhdif,zdtsdif, & |
---|
969 | sensibFlux,q2,zdqdif,zdqsdif, & |
---|
970 | taux,tauy,lastcall) |
---|
971 | |
---|
972 | zdtdif(:,:)=zdhdif(:,:)*zpopsk(:,:) ! for diagnostic only |
---|
973 | zdqevap(:,:)=0. |
---|
974 | |
---|
975 | end if !end of 'UseTurbDiff' |
---|
976 | |
---|
977 | |
---|
978 | pdv(:,:)=pdv(:,:)+zdvdif(:,:) |
---|
979 | pdu(:,:)=pdu(:,:)+zdudif(:,:) |
---|
980 | pdt(:,:)=pdt(:,:)+zdtdif(:,:) |
---|
981 | zdtsurf(:)=zdtsurf(:)+zdtsdif(:) |
---|
982 | |
---|
983 | if (tracer) then |
---|
984 | pdq(:,:,:)=pdq(:,:,:)+ zdqdif(:,:,:) |
---|
985 | dqsurf(:,:)=dqsurf(:,:) + zdqsdif(:,:) |
---|
986 | end if ! of if (tracer) |
---|
987 | |
---|
988 | |
---|
989 | ! test energy conservation |
---|
990 | !------------------------- |
---|
991 | if(enertest)then |
---|
992 | |
---|
993 | dEzdiff(:,:)=cpp*mass(:,:)*zdtdif(:,:) |
---|
994 | do ig = 1, ngrid |
---|
995 | dEdiff(ig)=SUM(dEzdiff (ig,:))+ sensibFlux(ig)! subtract flux to the ground |
---|
996 | dEzdiff(ig,1)= dEzdiff(ig,1)+ sensibFlux(ig)! subtract flux to the ground |
---|
997 | enddo |
---|
998 | |
---|
999 | call planetwide_sumval(dEdiff(:)*cell_area(:)/totarea_planet,dEtot) |
---|
1000 | dEdiffs(:)=capcal(:)*zdtsdif(:)-zflubid(:)-sensibFlux(:) |
---|
1001 | call planetwide_sumval(dEdiffs(:)*cell_area(:)/totarea_planet,dEtots) |
---|
1002 | call planetwide_sumval(sensibFlux(:)*cell_area(:)/totarea_planet,AtmToSurf_TurbFlux) |
---|
1003 | |
---|
1004 | if (is_master) then |
---|
1005 | |
---|
1006 | if (UseTurbDiff) then |
---|
1007 | print*,'In TurbDiff sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' |
---|
1008 | print*,'In TurbDiff non-cons atm nrj change =',dEtot,' W m-2' |
---|
1009 | print*,'In TurbDiff (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' |
---|
1010 | else |
---|
1011 | print*,'In vdifc sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' |
---|
1012 | print*,'In vdifc non-cons atm nrj change =',dEtot,' W m-2' |
---|
1013 | print*,'In vdifc (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' |
---|
1014 | end if |
---|
1015 | endif ! end of 'is_master' |
---|
1016 | |
---|
1017 | ! JL12 : note that the black body radiative flux emitted by the surface has been updated by the implicit scheme but not given back elsewhere. |
---|
1018 | endif ! end of 'enertest' |
---|
1019 | |
---|
1020 | else ! calldifv |
---|
1021 | |
---|
1022 | if(.not.newtonian)then |
---|
1023 | |
---|
1024 | zdtsurf(:) = zdtsurf(:) + (fluxrad(:) + fluxgrd(:))/capcal(:) |
---|
1025 | |
---|
1026 | endif |
---|
1027 | |
---|
1028 | endif ! end of 'calldifv' |
---|
1029 | |
---|
1030 | |
---|
1031 | !---------------------------------- |
---|
1032 | ! IV. Dry convective adjustment : |
---|
1033 | !---------------------------------- |
---|
1034 | |
---|
1035 | if(calladj) then |
---|
1036 | |
---|
1037 | zdh(:,:) = pdt(:,:)/zpopsk(:,:) |
---|
1038 | zduadj(:,:)=0.0 |
---|
1039 | zdvadj(:,:)=0.0 |
---|
1040 | zdhadj(:,:)=0.0 |
---|
1041 | |
---|
1042 | |
---|
1043 | call convadj(ngrid,nlayer,nq,ptimestep, & |
---|
1044 | pplay,pplev,zpopsk, & |
---|
1045 | pu,pv,zh,pq, & |
---|
1046 | pdu,pdv,zdh,pdq, & |
---|
1047 | zduadj,zdvadj,zdhadj, & |
---|
1048 | zdqadj) |
---|
1049 | |
---|
1050 | pdu(:,:) = pdu(:,:) + zduadj(:,:) |
---|
1051 | pdv(:,:) = pdv(:,:) + zdvadj(:,:) |
---|
1052 | pdt(:,:) = pdt(:,:) + zdhadj(:,:)*zpopsk(:,:) |
---|
1053 | zdtadj(:,:) = zdhadj(:,:)*zpopsk(:,:) ! for diagnostic only |
---|
1054 | |
---|
1055 | if(tracer) then |
---|
1056 | pdq(:,:,:) = pdq(:,:,:) + zdqadj(:,:,:) |
---|
1057 | end if |
---|
1058 | |
---|
1059 | ! Test energy conservation |
---|
1060 | if(enertest)then |
---|
1061 | call planetwide_sumval(cpp*massarea(:,:)*zdtadj(:,:)/totarea_planet,dEtot) |
---|
1062 | if (is_master) print*,'In convadj atmospheric energy change =',dEtot,' W m-2' |
---|
1063 | endif |
---|
1064 | |
---|
1065 | |
---|
1066 | endif ! end of 'calladj' |
---|
1067 | |
---|
1068 | |
---|
1069 | !--------------------------------------------- |
---|
1070 | ! V. Specific parameterizations for tracers |
---|
1071 | !--------------------------------------------- |
---|
1072 | |
---|
1073 | if (tracer) then |
---|
1074 | |
---|
1075 | ! ------------------- |
---|
1076 | ! V.1 Microphysics |
---|
1077 | ! ------------------- |
---|
1078 | |
---|
1079 | ! JVO 05/18 : We must call microphysics before chemistry, for condensation ! |
---|
1080 | |
---|
1081 | if (callmufi) then |
---|
1082 | #ifdef USE_QTEST |
---|
1083 | dtpq(:,:,:) = 0.0 ! we want tpq to go only through mufi |
---|
1084 | call calmufi(ptimestep,pplev,zzlev,pplay,zzlay,gzlat,pt,tpq,dtpq,zdqmufi) |
---|
1085 | tpq(:,:,:) = tpq(:,:,:) + zdqmufi(:,:,:)*ptimestep ! only manipulation of tpq->*ptimestep here |
---|
1086 | #else |
---|
1087 | call calmufi(ptimestep,pplev,zzlev,pplay,zzlay,gzlat,pt,pq,pdq,zdqmufi) ! JVO 19 : To be fixed, what altitude do we need ? |
---|
1088 | pdq(:,:,:) = pdq(:,:,:) + zdqmufi(:,:,:) |
---|
1089 | #endif |
---|
1090 | |
---|
1091 | ! Microphysics condensation for 2D fields to sent non-saturated fields to photochem |
---|
1092 | if ( callclouds .and. moyzon_ch .and. mod(icount-1,ichim).eq.0 ) then |
---|
1093 | zdqfibar(:,:,:) = 0.0 ! We work in zonal average -> forget processes other than condensation |
---|
1094 | call calmufi(ptimestep,zplevbar,zzlevbar,zplaybar,zzlaybar, & |
---|
1095 | gzlat,ztfibar,zqfibar,zdqfibar,zdqmufibar) |
---|
1096 | endif |
---|
1097 | |
---|
1098 | endif |
---|
1099 | |
---|
1100 | ! ----------------- |
---|
1101 | ! V.2. Chemistry |
---|
1102 | ! ----------------- |
---|
1103 | ! NB : Must be call last ( brings fields back to an equilibrium ) |
---|
1104 | |
---|
1105 | ! Known bug ? ( JVO 18 ) : If you try to use chimi_indx instead of iq+nmicro |
---|
1106 | ! it leads to weird results / crash on dev mod ( ok in debug ) ... Why ? Idk ... |
---|
1107 | |
---|
1108 | if (callchim) then |
---|
1109 | |
---|
1110 | ! o. Convert updated tracers to molar fraction |
---|
1111 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1112 | do iq = 1,nkim |
---|
1113 | ychim(:,:,iq) = ( pq(:,:,iq+nmicro) + pdq(:,:,iq+nmicro) ) / rat_mmol(iq+nmicro) |
---|
1114 | enddo |
---|
1115 | |
---|
1116 | ! JVO 05/18 : We update zonally averaged fields with condensation |
---|
1117 | ! as it is compulsory to have correct photochem production. But for other |
---|
1118 | ! processes ( convadj ... ) we miss them in any case as we work in zonally/diurnal |
---|
1119 | ! mean -> no fine diurnal/short time evolution, only seasonal evolution only. |
---|
1120 | if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) then |
---|
1121 | do iq = 1,nkim |
---|
1122 | ychimbar(:,:,iq) = zqfibar(:,:,iq+nmicro) / rat_mmol(iq+nmicro) |
---|
1123 | if ( callclouds ) then |
---|
1124 | ychimbar(:,:,iq) = ychimbar(:,:,iq) + ( zdqmufibar(:,:,iq+nmicro) / rat_mmol(iq+nmicro) ) |
---|
1125 | endif |
---|
1126 | enddo |
---|
1127 | endif |
---|
1128 | |
---|
1129 | ! i. Condensation of the 3D tracers after the transport |
---|
1130 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1131 | |
---|
1132 | call calc_ysat(ngrid,nlayer,pplay/100.0,pt,ysat) ! Compute saturation profiles for every grid point |
---|
1133 | |
---|
1134 | dyccond(:,:,:) = 0.0 ! Default value -> no condensation |
---|
1135 | |
---|
1136 | do iq=1,nkim |
---|
1137 | where ( ychim(:,:,iq).gt.ysat(:,:,iq) ) & |
---|
1138 | dyccond(:,:,iq+nmicro) = ( -ychim(:,:,iq)+ysat(:,:,iq) ) / ptimestep |
---|
1139 | enddo |
---|
1140 | |
---|
1141 | if ( callclouds ) dyccond(:,:,ices_indx) = 0.0 ! Condensation have been calculated in the cloud microphysics |
---|
1142 | |
---|
1143 | do iq=1,nkim |
---|
1144 | ychim(:,:,iq) = ychim(:,:,iq) + dyccond(:,:,iq+nmicro) ! update molar ychim for following calchim |
---|
1145 | |
---|
1146 | pdq(:,:,iq+nmicro) = pdq(:,:,iq+nmicro) + dyccond(:,:,iq+nmicro)*rat_mmol(iq+nmicro) ! convert tendencies to mass mixing ratio |
---|
1147 | enddo |
---|
1148 | |
---|
1149 | |
---|
1150 | ! ii. 2D zonally averaged fields need to condense and evap before photochem |
---|
1151 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1152 | if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) then |
---|
1153 | |
---|
1154 | call calc_ysat(ngrid,nlayer,zplaybar/100.0,ztfibar,ysat) ! Compute saturation profiles for every grid point for the zon-ave fields |
---|
1155 | |
---|
1156 | dyccondbar(:,:,:) = 0.0 ! Default value -> no condensation |
---|
1157 | |
---|
1158 | do iq = 1,nkim |
---|
1159 | where ( ychimbar(:,:,iq).gt.ysat(:,:,iq) ) & |
---|
1160 | dyccondbar(:,:,iq+nmicro) = ( -ychimbar(:,:,iq)+ysat(:,:,iq) ) / ptimestep |
---|
1161 | enddo |
---|
1162 | |
---|
1163 | if ( callclouds ) dyccondbar(:,:,ices_indx) = 0.0 ! Condensation have been calculated in the cloud microphysics |
---|
1164 | |
---|
1165 | do iq=1,nkim |
---|
1166 | ychimbar(:,:,iq) = ychimbar(:,:,iq) + dyccondbar(:,:,iq+nmicro) |
---|
1167 | enddo |
---|
1168 | |
---|
1169 | ! Pseudo-evap ( forcing constant surface humidity ) |
---|
1170 | do ig=1,ngrid |
---|
1171 | if ( ychimbar(ig,1,7+nmicro) .lt. botCH4 ) ychimbar(ig,1,7+nmicro) = botCH4 |
---|
1172 | enddo |
---|
1173 | |
---|
1174 | endif ! if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) |
---|
1175 | |
---|
1176 | ! iii. Photochemistry ( must be call after condensation (and evap of 2D) ) |
---|
1177 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1178 | if( mod(icount-1,ichim).eq.0. ) then |
---|
1179 | |
---|
1180 | print *, "We enter in the chemistry ..." |
---|
1181 | |
---|
1182 | if (moyzon_ch) then ! 2D zonally averaged chemistry |
---|
1183 | |
---|
1184 | ! Here we send zonal average fields ( corrected with cond ) from dynamics to chem. module |
---|
1185 | call calchim(ngrid,ychimbar,declin,ctimestep,ztfibar,zphibar,zphisbar, & |
---|
1186 | zplaybar,zplevbar,zzlaybar,zzlevbar,dycchi) |
---|
1187 | else ! 3D chemistry (or 1D run) |
---|
1188 | call calchim(ngrid,ychim,declin,ctimestep,pt,pphi,phisfi, & |
---|
1189 | pplay,pplev,zzlay_eff,zzlev_eff,dycchi) |
---|
1190 | endif ! if moyzon |
---|
1191 | |
---|
1192 | endif |
---|
1193 | |
---|
1194 | ! iv. Surface pseudo-evaporation |
---|
1195 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1196 | do ig=1,ngrid |
---|
1197 | if ( (ychim(ig,1,7+nmicro)+dycchi(ig,1,7+nmicro)*ptimestep) .lt. botCH4 ) then ! +dycchi because ychim not yet updated |
---|
1198 | dycevapCH4(ig) = ( -ychim(ig,1,7+nmicro)+botCH4 ) / ptimestep - dycchi(ig,1,7+nmicro) |
---|
1199 | else |
---|
1200 | dycevapCH4(ig) = 0.0 |
---|
1201 | endif |
---|
1202 | enddo |
---|
1203 | |
---|
1204 | pdq(:,1,7+nmicro) = pdq(:,1,7+nmicro) + dycevapCH4(:)*rat_mmol(7+nmicro) |
---|
1205 | |
---|
1206 | ! v. Updates and positivity check |
---|
1207 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1208 | do iq=1,nkim |
---|
1209 | zdqchi(:,:,iq+nmicro) = dycchi(:,:,iq)*rat_mmol(iq+nmicro) ! convert tendencies to mass mixing ratio |
---|
1210 | |
---|
1211 | where( (pq(:,:,iq+nmicro) + ( pdq(:,:,iq+nmicro)+zdqchi(:,:,iq+nmicro) )*ptimestep ) .LT. 0.) & ! When using zonal means we set the same tendency |
---|
1212 | zdqchi(:,:,iq+nmicro) = 1.e-30 - pdq(:,:,iq+nmicro) - pq(:,:,iq+nmicro)/ptimestep ! everywhere in longitude -> could lead to negs ! |
---|
1213 | enddo |
---|
1214 | |
---|
1215 | pdq(:,:,:) = pdq(:,:,:) + zdqchi(:,:,:) |
---|
1216 | |
---|
1217 | endif ! end of 'callchim' |
---|
1218 | |
---|
1219 | ! --------------- |
---|
1220 | ! V.3 Updates |
---|
1221 | ! --------------- |
---|
1222 | |
---|
1223 | ! Updating Atmospheric Mass and Tracers budgets. |
---|
1224 | if(mass_redistrib) then |
---|
1225 | |
---|
1226 | zdmassmr(:,:) = mass(:,:) * zdqevap(:,:) |
---|
1227 | |
---|
1228 | do ig = 1, ngrid |
---|
1229 | zdmassmr_col(ig)=SUM(zdmassmr(ig,:)) |
---|
1230 | enddo |
---|
1231 | |
---|
1232 | call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr) |
---|
1233 | call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col) |
---|
1234 | call writediagfi(ngrid,"mass","mass","kg/m2",3,mass) |
---|
1235 | |
---|
1236 | call mass_redistribution(ngrid,nlayer,nq,ptimestep, & |
---|
1237 | capcal,pplay,pplev,pt,tsurf,pq,qsurf, & |
---|
1238 | pu,pv,pdt,zdtsurf,pdq,pdu,pdv,zdmassmr, & |
---|
1239 | zdtmr,zdtsurfmr,zdpsrfmr,zdumr,zdvmr,zdqmr,zdqsurfmr) |
---|
1240 | |
---|
1241 | pdq(:,:,:) = pdq(:,:,:) + zdqmr(:,:,:) |
---|
1242 | dqsurf(:,:) = dqsurf(:,:) + zdqsurfmr(:,:) |
---|
1243 | pdt(:,:) = pdt(:,:) + zdtmr(:,:) |
---|
1244 | pdu(:,:) = pdu(:,:) + zdumr(:,:) |
---|
1245 | pdv(:,:) = pdv(:,:) + zdvmr(:,:) |
---|
1246 | pdpsrf(:) = pdpsrf(:) + zdpsrfmr(:) |
---|
1247 | zdtsurf(:) = zdtsurf(:) + zdtsurfmr(:) |
---|
1248 | |
---|
1249 | endif |
---|
1250 | |
---|
1251 | ! ----------------------------- |
---|
1252 | ! V.4. Surface Tracer Update |
---|
1253 | ! ----------------------------- |
---|
1254 | |
---|
1255 | qsurf(:,:) = qsurf(:,:) + ptimestep*dqsurf(:,:) |
---|
1256 | |
---|
1257 | ! Add qsurf to qsurf_hist, which is what we save in diagfi.nc. At the same time, we set the water |
---|
1258 | ! content of ocean gridpoints back to zero, in order to avoid rounding errors in vdifc, rain. |
---|
1259 | qsurf_hist(:,:) = qsurf(:,:) |
---|
1260 | |
---|
1261 | endif! end of if 'tracer' |
---|
1262 | |
---|
1263 | |
---|
1264 | !------------------------------------------------ |
---|
1265 | ! VI. Surface and sub-surface soil temperature |
---|
1266 | !------------------------------------------------ |
---|
1267 | |
---|
1268 | |
---|
1269 | ! Increment surface temperature |
---|
1270 | |
---|
1271 | tsurf(:)=tsurf(:)+ptimestep*zdtsurf(:) |
---|
1272 | |
---|
1273 | ! Compute soil temperatures and subsurface heat flux. |
---|
1274 | if (callsoil) then |
---|
1275 | call soil(ngrid,nsoilmx,.false.,lastcall,inertiedat, & |
---|
1276 | ptimestep,tsurf,tsoil,capcal,fluxgrd) |
---|
1277 | endif |
---|
1278 | |
---|
1279 | |
---|
1280 | ! Test energy conservation |
---|
1281 | if(enertest)then |
---|
1282 | call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf(:)/totarea_planet,dEtots) |
---|
1283 | if (is_master) print*,'Surface energy change =',dEtots,' W m-2' |
---|
1284 | endif |
---|
1285 | |
---|
1286 | |
---|
1287 | !--------------------------------------------------- |
---|
1288 | ! VII. Perform diagnostics and write output files |
---|
1289 | !--------------------------------------------------- |
---|
1290 | |
---|
1291 | ! Note : For output only: the actual model integration is performed in the dynamics. |
---|
1292 | |
---|
1293 | ! Temperature, zonal and meridional winds. |
---|
1294 | zt(:,:) = pt(:,:) + pdt(:,:)*ptimestep |
---|
1295 | zu(:,:) = pu(:,:) + pdu(:,:)*ptimestep |
---|
1296 | zv(:,:) = pv(:,:) + pdv(:,:)*ptimestep |
---|
1297 | |
---|
1298 | ! Diagnostic. |
---|
1299 | zdtdyn(:,:) = (pt(:,:)-ztprevious(:,:)) / ptimestep |
---|
1300 | ztprevious(:,:) = zt(:,:) |
---|
1301 | |
---|
1302 | zdudyn(:,:) = (pu(:,:)-zuprevious(:,:)) / ptimestep |
---|
1303 | zuprevious(:,:) = zu(:,:) |
---|
1304 | |
---|
1305 | if(firstcall)then |
---|
1306 | zdtdyn(:,:)=0.0 |
---|
1307 | zdudyn(:,:)=0.0 |
---|
1308 | endif |
---|
1309 | |
---|
1310 | ! Horizotal wind |
---|
1311 | zhorizwind(:,:) = sqrt( zu(:,:)*zu(:,:) + zv(:,:)*zv(:,:) ) |
---|
1312 | |
---|
1313 | ! Dynamical heating diagnostic. |
---|
1314 | do ig=1,ngrid |
---|
1315 | fluxdyn(ig)= SUM(zdtdyn(ig,:) *mass(ig,:))*cpp |
---|
1316 | enddo |
---|
1317 | |
---|
1318 | ! Tracers. |
---|
1319 | zq(:,:,:) = pq(:,:,:) + pdq(:,:,:)*ptimestep |
---|
1320 | |
---|
1321 | ! Surface pressure. |
---|
1322 | ps(:) = pplev(:,1) + pdpsrf(:)*ptimestep |
---|
1323 | |
---|
1324 | |
---|
1325 | |
---|
1326 | ! Surface and soil temperature information |
---|
1327 | call planetwide_sumval(cell_area(:)*tsurf(:)/totarea_planet,Ts1) |
---|
1328 | call planetwide_minval(tsurf(:),Ts2) |
---|
1329 | call planetwide_maxval(tsurf(:),Ts3) |
---|
1330 | if(callsoil)then |
---|
1331 | TsS = SUM(cell_area(:)*tsoil(:,nsoilmx))/totarea ! mean temperature at bottom soil layer |
---|
1332 | if (is_master) then |
---|
1333 | print*,' ave[Tsurf] min[Tsurf] max[Tsurf] ave[Tdeep]' |
---|
1334 | print*,Ts1,Ts2,Ts3,TsS |
---|
1335 | end if |
---|
1336 | else |
---|
1337 | if (is_master) then |
---|
1338 | print*,' ave[Tsurf] min[Tsurf] max[Tsurf]' |
---|
1339 | print*,Ts1,Ts2,Ts3 |
---|
1340 | endif |
---|
1341 | end if |
---|
1342 | |
---|
1343 | |
---|
1344 | ! Check the energy balance of the simulation during the run |
---|
1345 | if(corrk)then |
---|
1346 | |
---|
1347 | call planetwide_sumval(cell_area(:)*fluxtop_dn(:)/totarea_planet,ISR) |
---|
1348 | call planetwide_sumval(cell_area(:)*fluxabs_sw(:)/totarea_planet,ASR) |
---|
1349 | call planetwide_sumval(cell_area(:)*fluxtop_lw(:)/totarea_planet,OLR) |
---|
1350 | call planetwide_sumval(cell_area(:)*fluxgrd(:)/totarea_planet,GND) |
---|
1351 | call planetwide_sumval(cell_area(:)*fluxdyn(:)/totarea_planet,DYN) |
---|
1352 | do ig=1,ngrid |
---|
1353 | if(fluxtop_dn(ig).lt.0.0)then |
---|
1354 | print*,'fluxtop_dn has gone crazy' |
---|
1355 | print*,'fluxtop_dn=',fluxtop_dn(ig) |
---|
1356 | print*,'temp= ',pt(ig,:) |
---|
1357 | print*,'pplay= ',pplay(ig,:) |
---|
1358 | call abort |
---|
1359 | endif |
---|
1360 | end do |
---|
1361 | |
---|
1362 | if(ngrid.eq.1)then |
---|
1363 | DYN=0.0 |
---|
1364 | endif |
---|
1365 | |
---|
1366 | if (is_master) then |
---|
1367 | print*,' ISR ASR OLR GND DYN [W m^-2]' |
---|
1368 | print*, ISR,ASR,OLR,GND,DYN |
---|
1369 | endif |
---|
1370 | |
---|
1371 | if(enertest .and. is_master)then |
---|
1372 | print*,'SW flux/heating difference SW++ - ASR = ',dEtotSW+dEtotsSW-ASR,' W m-2' |
---|
1373 | print*,'LW flux/heating difference LW++ - OLR = ',dEtotLW+dEtotsLW+OLR,' W m-2' |
---|
1374 | print*,'LW energy balance LW++ + ASR = ',dEtotLW+dEtotsLW+ASR,' W m-2' |
---|
1375 | endif |
---|
1376 | |
---|
1377 | if(meanOLR .and. is_master)then |
---|
1378 | if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then |
---|
1379 | ! to record global radiative balance |
---|
1380 | open(92,file="rad_bal.out",form='formatted',position='append') |
---|
1381 | write(92,*) zday,ISR,ASR,OLR |
---|
1382 | close(92) |
---|
1383 | open(93,file="tem_bal.out",form='formatted',position='append') |
---|
1384 | if(callsoil)then |
---|
1385 | write(93,*) zday,Ts1,Ts2,Ts3,TsS |
---|
1386 | else |
---|
1387 | write(93,*) zday,Ts1,Ts2,Ts3 |
---|
1388 | endif |
---|
1389 | close(93) |
---|
1390 | endif |
---|
1391 | endif |
---|
1392 | |
---|
1393 | endif ! end of 'corrk' |
---|
1394 | |
---|
1395 | |
---|
1396 | ! Diagnostic to test radiative-convective timescales in code. |
---|
1397 | if(testradtimes)then |
---|
1398 | open(38,file="tau_phys.out",form='formatted',position='append') |
---|
1399 | ig=1 |
---|
1400 | do l=1,nlayer |
---|
1401 | write(38,*) -1./pdt(ig,l),pt(ig,l),pplay(ig,l) |
---|
1402 | enddo |
---|
1403 | close(38) |
---|
1404 | print*,'As testradtimes enabled,' |
---|
1405 | print*,'exiting physics on first call' |
---|
1406 | call abort |
---|
1407 | endif |
---|
1408 | |
---|
1409 | |
---|
1410 | if (is_master) print*,'--> Ls =',zls*180./pi |
---|
1411 | |
---|
1412 | |
---|
1413 | !---------------------------------------------------------------------- |
---|
1414 | ! Writing NetCDF file "RESTARTFI" at the end of the run |
---|
1415 | !---------------------------------------------------------------------- |
---|
1416 | |
---|
1417 | ! Note: 'restartfi' is stored just before dynamics are stored |
---|
1418 | ! in 'restart'. Between now and the writting of 'restart', |
---|
1419 | ! there will have been the itau=itau+1 instruction and |
---|
1420 | ! a reset of 'time' (lastacll = .true. when itau+1= itaufin) |
---|
1421 | ! thus we store for time=time+dtvr |
---|
1422 | |
---|
1423 | |
---|
1424 | |
---|
1425 | if(lastcall) then |
---|
1426 | ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec) |
---|
1427 | |
---|
1428 | if (ngrid.ne.1) then |
---|
1429 | write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin |
---|
1430 | |
---|
1431 | call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq, & |
---|
1432 | ptimestep,ztime_fin, & |
---|
1433 | tsurf,tsoil,emis,q2,qsurf_hist,tankCH4) |
---|
1434 | endif |
---|
1435 | endif ! end of 'lastcall' |
---|
1436 | |
---|
1437 | |
---|
1438 | !----------------------------------------------------------------------------------------------------- |
---|
1439 | ! OUTPUT in netcdf file "DIAGFI.NC", containing any variable for diagnostic |
---|
1440 | ! |
---|
1441 | ! Note 1 : output with period "ecritphy", set in "run.def" |
---|
1442 | ! |
---|
1443 | ! Note 2 : writediagfi can also be called from any other subroutine for any variable, |
---|
1444 | ! but its preferable to keep all the calls in one place ... |
---|
1445 | !----------------------------------------------------------------------------------------------------- |
---|
1446 | |
---|
1447 | |
---|
1448 | call writediagfi(ngrid,"Ls","solar longitude","deg",0,zls*180./pi) |
---|
1449 | call writediagfi(ngrid,"Lss","sub solar longitude","deg",0,zlss*180./pi) |
---|
1450 | call writediagfi(ngrid,"RA","right ascension","deg",0,right_ascen*180./pi) |
---|
1451 | call writediagfi(ngrid,"Declin","solar declination","deg",0,declin*180./pi) |
---|
1452 | call writediagfi(ngrid,"tsurf","Surface temperature","K",2,tsurf) |
---|
1453 | call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps) |
---|
1454 | call writediagfi(ngrid,"temp","temperature","K",3,zt) |
---|
1455 | call writediagfi(ngrid,"teta","potential temperature","K",3,zh) |
---|
1456 | call writediagfi(ngrid,"u","Zonal wind","m.s-1",3,zu) |
---|
1457 | call writediagfi(ngrid,"v","Meridional wind","m.s-1",3,zv) |
---|
1458 | call writediagfi(ngrid,"w","Vertical wind","m.s-1",3,pw) |
---|
1459 | call writediagfi(ngrid,"p","Pressure","Pa",3,pplay) |
---|
1460 | |
---|
1461 | ! Subsurface temperatures |
---|
1462 | ! call writediagsoil(ngrid,"tsurf","Surface temperature","K",2,tsurf) |
---|
1463 | ! call writediagsoil(ngrid,"temp","temperature","K",3,tsoil) |
---|
1464 | |
---|
1465 | |
---|
1466 | |
---|
1467 | ! Total energy balance diagnostics |
---|
1468 | if(callrad.and.(.not.newtonian))then |
---|
1469 | |
---|
1470 | !call writediagfi(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent) |
---|
1471 | call writediagfi(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn) |
---|
1472 | call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw) |
---|
1473 | call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw) |
---|
1474 | |
---|
1475 | ! call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1) |
---|
1476 | ! call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1) |
---|
1477 | ! call writediagfi(ngrid,"fluxsurfsw","sw surface flux.","W m-2",2,fluxsurf_sw) |
---|
1478 | ! call writediagfi(ngrid,"fluxsurflw","lw back radiation.","W m-2",2,fluxsurf_lw) |
---|
1479 | ! call writediagfi(ngrid,"fluxsurfswcs","sw surface flux (cs).","W m-2",2,fluxsurf_sw1) |
---|
1480 | ! call writediagfi(ngrid,"fluxsurflwcs","lw back radiation (cs).","W m-2",2,fluxsurf_lw1) |
---|
1481 | |
---|
1482 | |
---|
1483 | call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd) |
---|
1484 | |
---|
1485 | call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn) |
---|
1486 | |
---|
1487 | endif ! end of 'callrad' |
---|
1488 | |
---|
1489 | if(enertest) then |
---|
1490 | |
---|
1491 | if (calldifv) then |
---|
1492 | |
---|
1493 | call writediagfi(ngrid,"q2","turbulent kinetic energy","J.kg^-1",3,q2) |
---|
1494 | call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux) |
---|
1495 | |
---|
1496 | ! call writediagfi(ngrid,"dEzdiff","turbulent diffusion heating (-sensible flux)","w.m^-2",3,dEzdiff) |
---|
1497 | ! call writediagfi(ngrid,"dEdiff","integrated turbulent diffusion heating (-sensible flux)","w.m^-2",2,dEdiff) |
---|
1498 | ! call writediagfi(ngrid,"dEdiffs","In TurbDiff (correc rad+latent heat) surf nrj change","w.m^-2",2,dEdiffs) |
---|
1499 | |
---|
1500 | endif |
---|
1501 | |
---|
1502 | if (corrk) then |
---|
1503 | call writediagfi(ngrid,"dEzradsw","radiative heating","w.m^-2",3,dEzradsw) |
---|
1504 | call writediagfi(ngrid,"dEzradlw","radiative heating","w.m^-2",3,dEzradlw) |
---|
1505 | endif |
---|
1506 | |
---|
1507 | endif ! end of 'enertest' |
---|
1508 | |
---|
1509 | ! Temporary inclusions for winds diagnostics. |
---|
1510 | call writediagfi(ngrid,"zdudif","Turbdiff tend. zon. wind","m s-2",3,zdudif) |
---|
1511 | call writediagfi(ngrid,"zdudyn","Dyn. tend. zon. wind","m s-2",3,zdudyn) |
---|
1512 | |
---|
1513 | ! Temporary inclusions for heating diagnostics. |
---|
1514 | call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw) |
---|
1515 | call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw) |
---|
1516 | call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad) |
---|
1517 | call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn) |
---|
1518 | |
---|
1519 | ! For Debugging. |
---|
1520 | call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi) |
---|
1521 | |
---|
1522 | |
---|
1523 | ! Output tracers. |
---|
1524 | if (tracer) then |
---|
1525 | |
---|
1526 | if (callmufi) then |
---|
1527 | ! Microphysical tracers are expressed in unit/m3. |
---|
1528 | ! convert X.kg-1 --> X.m-3 (whereas for optics was -> X.m-2) |
---|
1529 | i2e(:,:) = ( pplev(:,1:nlayer)-pplev(:,2:nlayer+1) ) / gzlat(:,1:nlayer) /(zzlev(:,2:nlayer+1)-zzlev(:,1:nlayer)) |
---|
1530 | |
---|
1531 | #ifdef USE_QTEST |
---|
1532 | ! Microphysical tracers passed through dyn+phys(except mufi) |
---|
1533 | call writediagfi(ngrid,"mu_m0as_dp","Dynphys only spherical mode 0th order moment",'m-3',3,zq(:,:,micro_indx(1))*i2e) |
---|
1534 | call writediagfi(ngrid,"mu_m3as_dp","Dynphys only spherical mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(2))*i2e) |
---|
1535 | call writediagfi(ngrid,"mu_m0af_dp","Dynphys only fractal mode 0th order moment",'m-3',3,zq(:,:,micro_indx(3))*i2e) |
---|
1536 | call writediagfi(ngrid,"mu_m3af_dp","Dynphys only fractal mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(4))*i2e) |
---|
1537 | ! Microphysical tracers passed through mufi only |
---|
1538 | call writediagfi(ngrid,"mu_m0as_mo","Mufi only spherical mode 0th order moment",'m-3',3,tpq(:,:,micro_indx(1))*i2e) |
---|
1539 | call writediagfi(ngrid,"mu_m3as_mo","Mufi only spherical mode 3rd order moment",'m3/m3',3,tpq(:,:,micro_indx(2))*i2e) |
---|
1540 | call writediagfi(ngrid,"mu_m0af_mo","Mufi only fractal mode 0th order moment",'m-3',3,tpq(:,:,micro_indx(3))*i2e) |
---|
1541 | call writediagfi(ngrid,"mu_m3af_mo","Mufi only fractal mode 3rd order moment",'m3/m3',3,tpq(:,:,micro_indx(4))*i2e) |
---|
1542 | #else |
---|
1543 | call writediagfi(ngrid,"mu_m0as","Spherical mode 0th order moment",'m-3',3,zq(:,:,micro_indx(1))*i2e) |
---|
1544 | call writediagfi(ngrid,"mu_m3as","Spherical mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(2))*i2e) |
---|
1545 | call writediagfi(ngrid,"mu_m0af","Fractal mode 0th order moment",'m-3',3,zq(:,:,micro_indx(3))*i2e) |
---|
1546 | call writediagfi(ngrid,"mu_m3af","Fractal mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(4))*i2e) |
---|
1547 | #endif |
---|
1548 | |
---|
1549 | ! Microphysical diagnostics |
---|
1550 | call writediagfi(ngrid,"mmd_aer_prec","Total aerosols precipitations",'m',2,mmd_aer_prec) |
---|
1551 | call writediagfi(ngrid,"mmd_aer_s_flux","Spherical aerosols sedimentation flux",'kg.m-2.s-1',3,mmd_aer_s_flux) |
---|
1552 | call writediagfi(ngrid,"mmd_aer_f_flux","Fractal aerosols sedimentation flux",'kg.m-2.s-1',3,mmd_aer_f_flux) |
---|
1553 | call writediagfi(ngrid,"mmd_rc_sph","Spherical mode caracteristic radius",'m',3,mmd_rc_sph) |
---|
1554 | call writediagfi(ngrid,"mmd_rc_fra","Fractal mode caracteristic radius",'m',3,mmd_rc_fra) |
---|
1555 | |
---|
1556 | endif ! end of 'callmufi' |
---|
1557 | |
---|
1558 | ! Chemical tracers |
---|
1559 | if (callchim) then |
---|
1560 | do iq=1,nkim |
---|
1561 | call writediagfi(ngrid,cnames(iq),cnames(iq),'mol/mol',3,zq(:,:,iq+nmicro)/rat_mmol(iq+nmicro)) |
---|
1562 | enddo |
---|
1563 | call writediagfi(ngrid,"evapCH4","Surface CH4 pseudo-evaporation rate",'mol/mol/s',2,dycevapCH4) |
---|
1564 | endif |
---|
1565 | |
---|
1566 | endif ! end of 'tracer' |
---|
1567 | |
---|
1568 | ! XIOS outputs |
---|
1569 | #ifdef CPP_XIOS |
---|
1570 | ! Send fields to XIOS: (NB these fields must also be defined as |
---|
1571 | ! <field id="..." /> in context_lmdz_physics.xml to be correctly used) |
---|
1572 | CALL send_xios_field("ls",zls*180./pi) |
---|
1573 | CALL send_xios_field("lss",zlss*180./pi) |
---|
1574 | CALL send_xios_field("RA",right_ascen*180./pi) |
---|
1575 | CALL send_xios_field("Declin",declin*180./pi) |
---|
1576 | |
---|
1577 | ! Total energy balance diagnostics |
---|
1578 | if (callrad.and.(.not.newtonian)) then |
---|
1579 | CALL send_xios_field("ISR_TOA",fluxtop_dn) |
---|
1580 | CALL send_xios_field("OLR_TOA",fluxtop_lw) |
---|
1581 | endif |
---|
1582 | |
---|
1583 | CALL send_xios_field("area",cell_area) |
---|
1584 | CALL send_xios_field("pphi",pphi) |
---|
1585 | CALL send_xios_field("pphis",phisfi) |
---|
1586 | |
---|
1587 | CALL send_xios_field("ps",ps) |
---|
1588 | CALL send_xios_field("tsurf",tsurf) |
---|
1589 | |
---|
1590 | if(enertest) then |
---|
1591 | if (calldifv) then |
---|
1592 | CALL send_xios_field("sensibFlux",sensibFlux) |
---|
1593 | endif |
---|
1594 | endif |
---|
1595 | |
---|
1596 | CALL send_xios_field("temp",zt) |
---|
1597 | CALL send_xios_field("teta",zh) |
---|
1598 | CALL send_xios_field("u",zu) |
---|
1599 | CALL send_xios_field("v",zv) |
---|
1600 | CALL send_xios_field("w",pw) |
---|
1601 | CALL send_xios_field("p",pplay) |
---|
1602 | |
---|
1603 | ! Winds diagnostics. |
---|
1604 | CALL send_xios_field("dudif",zdudif) |
---|
1605 | CALL send_xios_field("dudyn",zdudyn) |
---|
1606 | |
---|
1607 | CALL send_xios_field("horizwind",zhorizwind) |
---|
1608 | |
---|
1609 | ! Heating diagnostics. |
---|
1610 | CALL send_xios_field("dtsw",zdtsw) |
---|
1611 | CALL send_xios_field("dtlw",zdtlw) |
---|
1612 | CALL send_xios_field("dtrad",dtrad) |
---|
1613 | CALL send_xios_field("dtdyn",zdtdyn) |
---|
1614 | CALL send_xios_field("dtdif",zdtdif) |
---|
1615 | |
---|
1616 | ! Chemical tracers |
---|
1617 | if (callchim) then |
---|
1618 | |
---|
1619 | ! Advected fields |
---|
1620 | do iq=1,nkim |
---|
1621 | CALL send_xios_field(trim(cnames(iq)),zq(:,:,iq+nmicro)/rat_mmol(iq+nmicro)) ! kg/kg -> mol/mol |
---|
1622 | enddo |
---|
1623 | |
---|
1624 | ! Upper chemistry fields |
---|
1625 | do iq=1,nkim |
---|
1626 | CALL send_xios_field(trim(cnames(iq))//"_up",ykim_up(iq,:,:)) ! mol/mol |
---|
1627 | enddo |
---|
1628 | |
---|
1629 | ! Append fields in ykim_tot for output on the total vertical grid (0->1300km) |
---|
1630 | do iq=1,nkim |
---|
1631 | |
---|
1632 | ! GCM levels |
---|
1633 | do l=1,nlayer |
---|
1634 | ykim_tot(iq,:,l) = zq(:,l,iq+nmicro)/rat_mmol(iq+nmicro) |
---|
1635 | enddo |
---|
1636 | ! Upper levels |
---|
1637 | do l=1,nlaykim_up |
---|
1638 | ykim_tot(iq,:,nlayer+l) = ykim_up(iq,:,l) |
---|
1639 | enddo |
---|
1640 | |
---|
1641 | CALL send_xios_field(trim(cnames(iq))//"_tot",ykim_tot(iq,:,:)) ! mol/mol |
---|
1642 | |
---|
1643 | enddo |
---|
1644 | |
---|
1645 | ! Condensation tendencies ( mol/mol/s ) |
---|
1646 | CALL send_xios_field("dqcond_CH4",dyccond(:,:,7+nmicro)) |
---|
1647 | CALL send_xios_field("dqcond_C2H2",dyccond(:,:,10+nmicro)) |
---|
1648 | CALL send_xios_field("dqcond_C2H4",dyccond(:,:,12+nmicro)) |
---|
1649 | CALL send_xios_field("dqcond_C2H6",dyccond(:,:,14+nmicro)) |
---|
1650 | CALL send_xios_field("dqcond_C3H6",dyccond(:,:,17+nmicro)) |
---|
1651 | CALL send_xios_field("dqcond_C4H4",dyccond(:,:,21+nmicro)) |
---|
1652 | CALL send_xios_field("dqcond_CH3CCH",dyccond(:,:,24+nmicro)) |
---|
1653 | CALL send_xios_field("dqcond_C3H8",dyccond(:,:,25+nmicro)) |
---|
1654 | CALL send_xios_field("dqcond_C4H2",dyccond(:,:,26+nmicro)) |
---|
1655 | CALL send_xios_field("dqcond_C4H6",dyccond(:,:,27+nmicro)) |
---|
1656 | CALL send_xios_field("dqcond_C4H10",dyccond(:,:,28+nmicro)) |
---|
1657 | CALL send_xios_field("dqcond_AC6H6",dyccond(:,:,29+nmicro)) |
---|
1658 | CALL send_xios_field("dqcond_HCN",dyccond(:,:,36+nmicro)) |
---|
1659 | CALL send_xios_field("dqcond_CH3CN",dyccond(:,:,40+nmicro)) |
---|
1660 | CALL send_xios_field("dqcond_HC3N",dyccond(:,:,42+nmicro)) |
---|
1661 | CALL send_xios_field("dqcond_NCCN",dyccond(:,:,43+nmicro)) |
---|
1662 | CALL send_xios_field("dqcond_C4N2",dyccond(:,:,44+nmicro)) |
---|
1663 | |
---|
1664 | ! Pseudo-evaporation flux (mol/mol/s) |
---|
1665 | CALL send_xios_field("evapCH4",dycevapCH4(:)) |
---|
1666 | |
---|
1667 | endif ! of 'if callchim' |
---|
1668 | |
---|
1669 | ! Microphysical tracers |
---|
1670 | if (callmufi) then |
---|
1671 | CALL send_xios_field("mu_m0as",zq(:,:,micro_indx(1))*i2e) |
---|
1672 | CALL send_xios_field("mu_m3as",zq(:,:,micro_indx(2))*i2e) |
---|
1673 | CALL send_xios_field("mu_m0af",zq(:,:,micro_indx(3))*i2e) |
---|
1674 | CALL send_xios_field("mu_m3af",zq(:,:,micro_indx(4))*i2e) |
---|
1675 | endif |
---|
1676 | |
---|
1677 | if (lastcall.and.is_omp_master) then |
---|
1678 | write(*,*) "physiq: call xios_context_finalize" |
---|
1679 | call xios_context_finalize |
---|
1680 | endif |
---|
1681 | #endif |
---|
1682 | |
---|
1683 | icount=icount+1 |
---|
1684 | |
---|
1685 | end subroutine physiq |
---|
1686 | |
---|
1687 | end module physiq_mod |
---|