1 | MODULE timeloop_gcm_mod |
---|
2 | USE profiling_mod |
---|
3 | USE icosa |
---|
4 | USE disvert_mod |
---|
5 | USE trace |
---|
6 | USE omp_para |
---|
7 | USE euler_scheme_mod |
---|
8 | USE explicit_scheme_mod |
---|
9 | USE hevi_scheme_mod |
---|
10 | IMPLICIT NONE |
---|
11 | PRIVATE |
---|
12 | |
---|
13 | INTEGER, PARAMETER :: itau_sync=10 |
---|
14 | TYPE(t_message),SAVE :: req_ps0, req_mass0, req_theta_rhodz0, req_u0, req_q0, req_W0, req_geopot0 |
---|
15 | LOGICAL, SAVE :: positive_theta |
---|
16 | INTEGER :: itau_prof, id_timeloop, id_dyn, id_phys, id_dissip, id_adv, id_diags |
---|
17 | PUBLIC :: init_timeloop, timeloop |
---|
18 | |
---|
19 | CONTAINS |
---|
20 | |
---|
21 | SUBROUTINE init_timeloop |
---|
22 | USE dissip_gcm_mod |
---|
23 | USE observable_mod |
---|
24 | USE caldyn_mod |
---|
25 | USE etat0_mod |
---|
26 | USE guided_mod |
---|
27 | USE advect_tracer_mod |
---|
28 | USE check_conserve_mod |
---|
29 | USE output_field_mod |
---|
30 | USE theta2theta_rhodz_mod |
---|
31 | USE sponge_mod |
---|
32 | |
---|
33 | CHARACTER(len=255) :: def |
---|
34 | |
---|
35 | CALL register_id('timeloop',id_timeloop) |
---|
36 | CALL register_id('dyn',id_dyn) |
---|
37 | CALL register_id('dissip',id_dissip) |
---|
38 | CALL register_id('phys',id_phys) |
---|
39 | CALL register_id('adv',id_adv) |
---|
40 | CALL register_id('diags',id_diags) |
---|
41 | |
---|
42 | CALL init_caldyn |
---|
43 | |
---|
44 | ! IF (xios_output) itau_out=1 |
---|
45 | IF (.NOT. enable_io) itau_out=HUGE(itau_out) |
---|
46 | |
---|
47 | itau_prof=100 |
---|
48 | CALL getin('itau_profiling',itau_prof) |
---|
49 | |
---|
50 | positive_theta=.FALSE. |
---|
51 | CALL getin('positive_theta',positive_theta) |
---|
52 | IF(positive_theta .AND. nqtot<1) THEN |
---|
53 | PRINT *, 'nqtot must be >0 if positive_theta is .TRUE.' |
---|
54 | STOP |
---|
55 | END IF |
---|
56 | |
---|
57 | def='ARK2.3' |
---|
58 | CALL getin('time_scheme',def) |
---|
59 | SELECT CASE (TRIM(def)) |
---|
60 | CASE('euler') |
---|
61 | scheme_family=explicit |
---|
62 | scheme=euler |
---|
63 | nb_stage=1 |
---|
64 | CASE ('runge_kutta') |
---|
65 | scheme_family=explicit |
---|
66 | scheme=rk4 |
---|
67 | nb_stage=4 |
---|
68 | CASE ('RK2.5') |
---|
69 | scheme_family=explicit |
---|
70 | scheme=rk25 |
---|
71 | nb_stage=5 |
---|
72 | CASE ('leapfrog_matsuno') |
---|
73 | scheme_family=explicit |
---|
74 | scheme=mlf |
---|
75 | matsuno_period=5 |
---|
76 | CALL getin('matsuno_period',matsuno_period) |
---|
77 | nb_stage=matsuno_period+1 |
---|
78 | CASE('ARK2.3') |
---|
79 | scheme_family=hevi |
---|
80 | scheme=ark23 |
---|
81 | nb_stage=3 |
---|
82 | CALL set_coefs_ark23(dt) |
---|
83 | CASE('ARK3.3') |
---|
84 | scheme_family=hevi |
---|
85 | scheme=ark33 |
---|
86 | nb_stage=3 |
---|
87 | CALL set_coefs_ark33(dt) |
---|
88 | CASE ('none') |
---|
89 | nb_stage=0 |
---|
90 | CASE default |
---|
91 | PRINT*,'Bad selector for variable scheme : <', TRIM(def), & |
---|
92 | ' > options are <euler>, <runge_kutta>, <leapfrog_matsuno>,<RK2.5>,<ARK2.3>' |
---|
93 | STOP |
---|
94 | END SELECT |
---|
95 | |
---|
96 | ! Time-independant orography |
---|
97 | CALL allocate_field(f_phis,field_t,type_real,name='phis') |
---|
98 | ! Model state at current time step |
---|
99 | CALL allocate_field(f_ps,field_t,type_real, name='ps') |
---|
100 | CALL allocate_field(f_mass,field_t,type_real,llm,name='mass') |
---|
101 | CALL allocate_field(f_rhodz,field_t,type_real,llm,name='rhodz') |
---|
102 | CALL allocate_field(f_theta_rhodz,field_t,type_real,llm,nqdyn,name='theta_rhodz') |
---|
103 | CALL allocate_field(f_u,field_u,type_real,llm,name='u') |
---|
104 | CALL allocate_field(f_geopot,field_t,type_real,llm+1,name='geopot') |
---|
105 | CALL allocate_field(f_W,field_t,type_real,llm+1,name='W') |
---|
106 | CALL allocate_field(f_q,field_t,type_real,llm,nqtot,'q') |
---|
107 | ! Mass fluxes |
---|
108 | CALL allocate_field(f_hflux,field_u,type_real,llm) ! instantaneous mass fluxes |
---|
109 | CALL allocate_field(f_hfluxt,field_u,type_real,llm) ! mass "fluxes" accumulated in time |
---|
110 | CALL allocate_field(f_wflux,field_t,type_real,llm+1) ! vertical mass fluxes |
---|
111 | CALL allocate_field(f_wfluxt,field_t,type_real,llm+1,name='wfluxt') |
---|
112 | |
---|
113 | SELECT CASE(scheme_family) |
---|
114 | CASE(explicit) |
---|
115 | ! Trends |
---|
116 | CALL allocate_field(f_dps,field_t,type_real,name='dps') |
---|
117 | CALL allocate_field(f_dmass,field_t,type_real,llm, name='dmass') |
---|
118 | CALL allocate_field(f_dtheta_rhodz,field_t,type_real,llm,nqdyn,name='dtheta_rhodz') |
---|
119 | CALL allocate_field(f_du,field_u,type_real,llm,name='du') |
---|
120 | ! Model state at previous time step (RK/MLF) |
---|
121 | CALL allocate_field(f_psm1,field_t,type_real,name='psm1') |
---|
122 | CALL allocate_field(f_massm1,field_t,type_real,llm, name='massm1') |
---|
123 | CALL allocate_field(f_theta_rhodzm1,field_t,type_real,llm,nqdyn,name='theta_rhodzm1') |
---|
124 | CALL allocate_field(f_um1,field_u,type_real,llm,name='um1') |
---|
125 | CASE(hevi) |
---|
126 | ! Trends |
---|
127 | CALL allocate_fields(nb_stage,f_dps_slow, field_t,type_real,name='dps_slow') |
---|
128 | CALL allocate_fields(nb_stage,f_dmass_slow, field_t,type_real,llm, name='dmass_slow') |
---|
129 | CALL allocate_fields(nb_stage,f_dtheta_rhodz_slow, field_t,type_real,llm,nqdyn,name='dtheta_rhodz_fast') |
---|
130 | CALL allocate_fields(nb_stage,f_du_slow, field_u,type_real,llm,name='du_slow') |
---|
131 | CALL allocate_fields(nb_stage,f_du_fast, field_u,type_real,llm,name='du_fast') |
---|
132 | CALL allocate_fields(nb_stage,f_dW_slow, field_t,type_real,llm+1,name='dW_slow') |
---|
133 | CALL allocate_fields(nb_stage,f_dW_fast, field_t,type_real,llm+1,name='dW_fast') |
---|
134 | CALL allocate_fields(nb_stage,f_dPhi_slow, field_t,type_real,llm+1,name='dPhi_slow') |
---|
135 | CALL allocate_fields(nb_stage,f_dPhi_fast, field_t,type_real,llm+1,name='dPhi_fast') |
---|
136 | f_dps => f_dps_slow(:,1) |
---|
137 | f_du => f_du_slow(:,1) |
---|
138 | f_dtheta_rhodz => f_dtheta_rhodz_slow(:,1) |
---|
139 | END SELECT |
---|
140 | |
---|
141 | SELECT CASE(scheme) |
---|
142 | CASE(mlf) |
---|
143 | ! Model state 2 time steps ago (MLF) |
---|
144 | CALL allocate_field(f_psm2,field_t,type_real) |
---|
145 | CALL allocate_field(f_massm2,field_t,type_real,llm) |
---|
146 | CALL allocate_field(f_theta_rhodzm2,field_t,type_real,llm,nqdyn) |
---|
147 | CALL allocate_field(f_um2,field_u,type_real,llm) |
---|
148 | END SELECT |
---|
149 | |
---|
150 | CALL init_theta2theta_rhodz |
---|
151 | CALL init_dissip |
---|
152 | CALL init_sponge |
---|
153 | CALL init_observable |
---|
154 | CALL init_guided |
---|
155 | CALL init_advect_tracer |
---|
156 | CALL init_check_conserve |
---|
157 | |
---|
158 | CALL etat0(f_ps,f_mass,f_phis,f_theta_rhodz,f_u, f_geopot,f_W, f_q) |
---|
159 | |
---|
160 | CALL transfert_request(f_phis,req_i0) |
---|
161 | CALL transfert_request(f_phis,req_i1) |
---|
162 | |
---|
163 | CALL init_message(f_ps,req_i0,req_ps0) |
---|
164 | CALL init_message(f_mass,req_i0,req_mass0) |
---|
165 | CALL init_message(f_theta_rhodz,req_i0,req_theta_rhodz0) |
---|
166 | CALL init_message(f_u,req_e0_vect,req_u0) |
---|
167 | CALL init_message(f_q,req_i0,req_q0) |
---|
168 | CALL init_message(f_geopot,req_i0,req_geopot0) |
---|
169 | CALL init_message(f_W,req_i0,req_W0) |
---|
170 | |
---|
171 | END SUBROUTINE init_timeloop |
---|
172 | |
---|
173 | SUBROUTINE timeloop |
---|
174 | USE dissip_gcm_mod |
---|
175 | USE sponge_mod |
---|
176 | USE observable_mod |
---|
177 | USE etat0_mod |
---|
178 | USE guided_mod |
---|
179 | USE caldyn_mod |
---|
180 | USE advect_tracer_mod |
---|
181 | USE diagflux_mod |
---|
182 | USE physics_mod |
---|
183 | USE mpipara |
---|
184 | USE transfert_mod |
---|
185 | USE check_conserve_mod |
---|
186 | USE xios_mod |
---|
187 | USE output_field_mod |
---|
188 | USE write_etat0_mod |
---|
189 | USE restart_mod |
---|
190 | USE checksum_mod |
---|
191 | USE explicit_scheme_mod |
---|
192 | USE hevi_scheme_mod |
---|
193 | REAL(rstd),POINTER :: rhodz(:,:), mass(:,:), ps(:) |
---|
194 | |
---|
195 | REAL(rstd) :: adv_over_out ! ratio itau_adv/itau_out |
---|
196 | INTEGER :: ind, it,i,j,l,n, stage |
---|
197 | LOGICAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating mass fluxes in time |
---|
198 | LOGICAL, PARAMETER :: check_rhodz=.FALSE. |
---|
199 | INTEGER(kind=8) :: start_clock, stop_clock, rate_clock |
---|
200 | LOGICAL,SAVE :: first_physic=.TRUE. |
---|
201 | !$OMP THREADPRIVATE(first_physic) |
---|
202 | |
---|
203 | CALL switch_omp_distrib_level |
---|
204 | CALL caldyn_BC(f_phis, f_geopot, f_wflux) ! set constant values in first/last interfaces |
---|
205 | |
---|
206 | !$OMP BARRIER |
---|
207 | DO ind=1,ndomain |
---|
208 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
209 | CALL swap_dimensions(ind) |
---|
210 | CALL swap_geometry(ind) |
---|
211 | rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind) |
---|
212 | IF(caldyn_eta==eta_mass) THEN |
---|
213 | CALL compute_rhodz(.TRUE., ps, rhodz) ! save rhodz for transport scheme before dynamics update ps |
---|
214 | ELSE |
---|
215 | DO l=ll_begin,ll_end |
---|
216 | rhodz(:,l)=mass(:,l) |
---|
217 | ENDDO |
---|
218 | END IF |
---|
219 | END DO |
---|
220 | !$OMP BARRIER |
---|
221 | fluxt_zero=.TRUE. |
---|
222 | |
---|
223 | IF(positive_theta) CALL copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q) |
---|
224 | IF(diagflux_on) THEN |
---|
225 | adv_over_out = itau_adv*(1./itau_out) |
---|
226 | ELSE |
---|
227 | adv_over_out = 0. |
---|
228 | END IF |
---|
229 | |
---|
230 | CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,itau0) |
---|
231 | |
---|
232 | Call trace_on |
---|
233 | |
---|
234 | IF (xios_output) THEN ! we must call update_calendar before any XIOS output |
---|
235 | IF (is_omp_master) CALL xios_update_calendar(1) |
---|
236 | END IF |
---|
237 | |
---|
238 | ! IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q) |
---|
239 | IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W) |
---|
240 | |
---|
241 | CALL write_output_fields_basic(.TRUE., f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q) |
---|
242 | |
---|
243 | !$OMP MASTER |
---|
244 | CALL SYSTEM_CLOCK(start_clock, rate_clock) |
---|
245 | !$OMP END MASTER |
---|
246 | |
---|
247 | DO it=itau0+1,itau0+itaumax |
---|
248 | IF (is_master) CALL print_iteration(it, itau0, itaumax, start_clock, rate_clock) |
---|
249 | |
---|
250 | CALL enter_profile(id_timeloop) |
---|
251 | |
---|
252 | IF (xios_output) THEN |
---|
253 | |
---|
254 | IF(it>itau0+1 .AND. is_omp_master) CALL xios_update_calendar(it-itau0) |
---|
255 | ELSE |
---|
256 | CALL update_time_counter(dt*it) |
---|
257 | END IF |
---|
258 | |
---|
259 | IF (it==itau0+1 .OR. MOD(it,itau_sync)==0) THEN |
---|
260 | CALL send_message(f_ps,req_ps0) |
---|
261 | CALL wait_message(req_ps0) |
---|
262 | CALL send_message(f_mass,req_mass0) |
---|
263 | CALL wait_message(req_mass0) |
---|
264 | CALL send_message(f_theta_rhodz,req_theta_rhodz0) |
---|
265 | CALL wait_message(req_theta_rhodz0) |
---|
266 | CALL send_message(f_u,req_u0) |
---|
267 | CALL wait_message(req_u0) |
---|
268 | CALL send_message(f_q,req_q0) |
---|
269 | CALL wait_message(req_q0) |
---|
270 | IF(.NOT. hydrostatic) THEN |
---|
271 | CALL send_message(f_geopot,req_geopot0) |
---|
272 | CALL wait_message(req_geopot0) |
---|
273 | CALL send_message(f_W,req_W0) |
---|
274 | CALL wait_message(req_W0) |
---|
275 | END IF |
---|
276 | ENDIF |
---|
277 | |
---|
278 | CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q) |
---|
279 | |
---|
280 | CALL enter_profile(id_dyn) |
---|
281 | SELECT CASE(scheme_family) |
---|
282 | CASE(explicit) |
---|
283 | CALL explicit_scheme(it, fluxt_zero) |
---|
284 | CASE(hevi) |
---|
285 | CALL HEVI_scheme(it, fluxt_zero) |
---|
286 | END SELECT |
---|
287 | CALL exit_profile(id_dyn) |
---|
288 | |
---|
289 | CALL enter_profile(id_dissip) |
---|
290 | IF (MOD(it,itau_dissip)==0) THEN |
---|
291 | |
---|
292 | IF(caldyn_eta==eta_mass) THEN |
---|
293 | !ym flush ps |
---|
294 | !$OMP BARRIER |
---|
295 | DO ind=1,ndomain |
---|
296 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
297 | CALL swap_dimensions(ind) |
---|
298 | CALL swap_geometry(ind) |
---|
299 | mass=f_mass(ind); ps=f_ps(ind); |
---|
300 | CALL compute_rhodz(.TRUE., ps, mass) |
---|
301 | END DO |
---|
302 | ENDIF |
---|
303 | |
---|
304 | CALL enter_profile(id_diags) |
---|
305 | CALL check_conserve_detailed(it, AAM_dyn, & |
---|
306 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
307 | CALL exit_profile(id_diags) |
---|
308 | |
---|
309 | CALL dissip(f_ps,f_mass,f_phis,f_geopot,f_theta_rhodz,f_u, f_dtheta_rhodz,f_du) |
---|
310 | |
---|
311 | CALL euler_scheme(.FALSE.) ! update only u, theta |
---|
312 | IF (iflag_sponge > 0) THEN |
---|
313 | CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz) |
---|
314 | CALL euler_scheme(.FALSE.) ! update only u, theta |
---|
315 | END IF |
---|
316 | |
---|
317 | CALL enter_profile(id_diags) |
---|
318 | CALL check_conserve_detailed(it, AAM_dissip, & |
---|
319 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
320 | CALL exit_profile(id_diags) |
---|
321 | END IF |
---|
322 | CALL exit_profile(id_dissip) |
---|
323 | |
---|
324 | CALL enter_profile(id_adv) |
---|
325 | IF(MOD(it,itau_adv)==0) THEN |
---|
326 | CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz, & ! update q and rhodz after RK step |
---|
327 | adv_over_out, f_masst,f_qmasst,f_massfluxt, f_qfluxt) ! accumulate mass and fluxes if diagflux_on |
---|
328 | fluxt_zero=.TRUE. ! restart accumulation of hfluxt and qfluxt at next call to explicit_scheme / HEVI_scheme |
---|
329 | ! At this point advect_tracer has obtained the halos of u and rhodz, |
---|
330 | ! needed for correct computation of kinetic energy |
---|
331 | IF(diagflux_on) CALL diagflux_energy(adv_over_out, f_phis,f_rhodz,f_theta_rhodz,f_u, f_geopot,f_theta,f_buf_i, f_hfluxt) |
---|
332 | |
---|
333 | IF (check_rhodz) THEN ! check that rhodz is consistent with ps |
---|
334 | DO ind=1,ndomain |
---|
335 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
336 | CALL swap_dimensions(ind) |
---|
337 | CALL swap_geometry(ind) |
---|
338 | rhodz=f_rhodz(ind); ps=f_ps(ind); |
---|
339 | CALL compute_rhodz(.FALSE., ps, rhodz) |
---|
340 | END DO |
---|
341 | ENDIF |
---|
342 | IF(positive_theta) CALL copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q) |
---|
343 | END IF |
---|
344 | CALL exit_profile(id_adv) |
---|
345 | |
---|
346 | CALL enter_profile(id_diags) |
---|
347 | ! IF (MOD(it,itau_physics)==0) THEN |
---|
348 | CALL check_conserve_detailed(it, AAM_dyn, & |
---|
349 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
350 | CALL enter_profile(id_phys) |
---|
351 | CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q) |
---|
352 | CALL exit_profile(id_phys) |
---|
353 | CALL check_conserve_detailed(it, AAM_phys, & |
---|
354 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
355 | !$OMP MASTER |
---|
356 | IF (first_physic) CALL SYSTEM_CLOCK(start_clock) |
---|
357 | !$OMP END MASTER |
---|
358 | first_physic=.FALSE. |
---|
359 | ! END IF |
---|
360 | |
---|
361 | IF (MOD(it,itau_check_conserv)==0) THEN |
---|
362 | CALL check_conserve_detailed(it, AAM_dyn, & |
---|
363 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
364 | CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) |
---|
365 | ENDIF |
---|
366 | |
---|
367 | IF (mod(it,itau_out)==0 ) THEN |
---|
368 | CALL transfert_request(f_u,req_e1_vect) |
---|
369 | CALL write_output_fields_basic(.FALSE.,f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q) |
---|
370 | ENDIF |
---|
371 | CALL exit_profile(id_diags) |
---|
372 | |
---|
373 | CALL exit_profile(id_timeloop) |
---|
374 | END DO |
---|
375 | |
---|
376 | ! CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q) |
---|
377 | CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W) |
---|
378 | |
---|
379 | CALL check_conserve_detailed(it, AAM_dyn, & |
---|
380 | f_ps,f_dps,f_u,f_theta_rhodz,f_phis) |
---|
381 | CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) |
---|
382 | |
---|
383 | !$OMP MASTER |
---|
384 | CALL SYSTEM_CLOCK(stop_clock) |
---|
385 | |
---|
386 | IF (mpi_rank==0) THEN |
---|
387 | PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock |
---|
388 | ENDIF |
---|
389 | !$OMP END MASTER |
---|
390 | |
---|
391 | ! CONTAINS |
---|
392 | END SUBROUTINE timeloop |
---|
393 | |
---|
394 | SUBROUTINE print_iteration(it,itau0,itaumax,start_clock,rate_clock) |
---|
395 | INTEGER :: it, itau0, itaumax, throughput |
---|
396 | INTEGER(kind=8) :: start_clock, stop_clock, rate_clock |
---|
397 | REAL :: per_step,total, elapsed |
---|
398 | WRITE(*,'(A,I7,A,F14.1)') "It No :",it," t :",dt*it |
---|
399 | IF(MOD(it,10)==0) THEN |
---|
400 | CALL SYSTEM_CLOCK(stop_clock) |
---|
401 | elapsed = (stop_clock-start_clock)*1./rate_clock |
---|
402 | per_step = elapsed/(it-itau0) |
---|
403 | throughput = dt/per_step |
---|
404 | total = per_step*itaumax |
---|
405 | WRITE(*,'(A,I5,A,F6.2,A,I6)') 'Time spent (s):',INT(elapsed), & |
---|
406 | ' -- ms/step : ', 1000*per_step, & |
---|
407 | ' -- Throughput :', throughput |
---|
408 | WRITE(*,'(A,I5,A,I5)') 'Whole job (min) :', INT(total/60.), & |
---|
409 | ' -- Completion in (min) : ', INT((total-elapsed)/60.) |
---|
410 | END IF |
---|
411 | IF(MOD(it,itau_prof)==0) CALL print_profile |
---|
412 | |
---|
413 | END SUBROUTINE print_iteration |
---|
414 | |
---|
415 | SUBROUTINE copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q) |
---|
416 | TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:) |
---|
417 | REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:) |
---|
418 | INTEGER :: ind, iq |
---|
419 | DO ind=1, ndomain |
---|
420 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
421 | CALL swap_dimensions(ind) |
---|
422 | CALL swap_geometry(ind) |
---|
423 | theta_rhodz=f_theta_rhodz(ind) |
---|
424 | rhodz=f_rhodz(ind) |
---|
425 | q=f_q(ind) |
---|
426 | DO iq=1, nqdyn |
---|
427 | q(:,:,iq) = theta_rhodz(:,:,iq)/rhodz(:,:) |
---|
428 | END DO |
---|
429 | END DO |
---|
430 | END SUBROUTINE copy_theta_to_q |
---|
431 | |
---|
432 | SUBROUTINE copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q) |
---|
433 | TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:) |
---|
434 | REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:) |
---|
435 | INTEGER :: ind, iq |
---|
436 | DO ind=1, ndomain |
---|
437 | IF (.NOT. assigned_domain(ind)) CYCLE |
---|
438 | CALL swap_dimensions(ind) |
---|
439 | CALL swap_geometry(ind) |
---|
440 | theta_rhodz=f_theta_rhodz(ind) |
---|
441 | rhodz=f_rhodz(ind) |
---|
442 | q=f_q(ind) |
---|
443 | DO iq=1,nqdyn |
---|
444 | theta_rhodz(:,:,iq) = rhodz(:,:)*q(:,:,iq) |
---|
445 | END DO |
---|
446 | END DO |
---|
447 | END SUBROUTINE copy_q_to_theta |
---|
448 | |
---|
449 | END MODULE timeloop_gcm_mod |
---|