source: codes/icosagcm/trunk/src/time/timeloop_gcm.F90 @ 1015

Last change on this file since 1015 was 1015, checked in by ymipsl, 4 years ago

1j+1j=2j

  • the iteration frequency of syncronization procedure at frontier (every 10 time step before) must be a sub multiple of the restart period (daily). So we fix itau_sync to be close of 10 and be a sub multiple of number of iteration in a day.
  • The value of rhodz for advection must be recomputed at a sub-multiple multiple of the restart period, so for now we decide to recomputing it each time at the end of the advection.

YM

File size: 18.5 KB
Line 
1MODULE 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 :: sync_it=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
19CONTAINS
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=1000
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    IF (scheme_family /= hevi) THEN
97       CALL abort_acc("scheme_family /= hevi")
98    END IF
99
100    ! Time-independant orography
101    CALL allocate_field(f_phis,field_t,type_real,name='phis')
102    ! Model state at current time step
103    CALL allocate_field(f_ps,field_t,type_real, name='ps')
104    CALL allocate_field(f_mass,field_t,type_real,llm,name='mass')
105    CALL allocate_field(f_rhodz,field_t,type_real,llm,name='rhodz')
106    CALL allocate_field(f_theta_rhodz,field_t,type_real,llm,nqdyn,name='theta_rhodz')
107    CALL allocate_field(f_u,field_u,type_real,llm,name='u')
108    CALL allocate_field(f_geopot,field_t,type_real,llm+1,name='geopot')
109    CALL allocate_field(f_W,field_t,type_real,llm+1,name='W') ! used only if .not. hydrostatic
110    CALL allocate_field(f_q,field_t,type_real,llm,nqtot,'q')
111    ! Mass fluxes
112    CALL allocate_field(f_hflux,field_u,type_real,llm, ondevice=.TRUE.)    ! instantaneous mass fluxes
113    CALL allocate_field(f_hfluxt,field_u,type_real,llm,ondevice=.TRUE.)   ! mass "fluxes" accumulated in time
114    CALL allocate_field(f_wflux,field_t,type_real,llm+1)  ! vertical mass fluxes
115    CALL allocate_field(f_wfluxt,field_t,type_real,llm+1,name='wfluxt',ondevice=.TRUE.)
116   
117    SELECT CASE(scheme_family)
118    CASE(explicit)
119       ! Trends
120       CALL allocate_field(f_dps,field_t,type_real,name='dps')
121       CALL allocate_field(f_dmass,field_t,type_real,llm, name='dmass')
122       CALL allocate_field(f_dtheta_rhodz,field_t,type_real,llm,nqdyn,name='dtheta_rhodz')
123       CALL allocate_field(f_du,field_u,type_real,llm,name='du')
124       ! Model state at previous time step (RK/MLF)
125       CALL allocate_field(f_psm1,field_t,type_real,name='psm1')
126       CALL allocate_field(f_massm1,field_t,type_real,llm, name='massm1')
127       CALL allocate_field(f_theta_rhodzm1,field_t,type_real,llm,nqdyn,name='theta_rhodzm1')
128       CALL allocate_field(f_um1,field_u,type_real,llm,name='um1')
129    CASE(hevi)
130       ! Trends
131       CALL allocate_fields(nb_stage,f_dps_slow, field_t,type_real,name='dps_slow', ondevice=.TRUE.)
132       CALL allocate_fields(nb_stage,f_dmass_slow, field_t,type_real,llm, name='dmass_slow', ondevice=.TRUE.)
133       CALL allocate_fields(nb_stage,f_dtheta_rhodz_slow, field_t,type_real,llm,nqdyn,name='dtheta_rhodz_fast', ondevice=.TRUE.)
134       CALL allocate_fields(nb_stage,f_du_slow, field_u,type_real,llm,name='du_slow', ondevice=.TRUE.)
135       CALL allocate_fields(nb_stage,f_du_fast, field_u,type_real,llm,name='du_fast', ondevice=.TRUE.)
136       CALL allocate_fields(nb_stage,f_dW_slow, field_t,type_real,llm+1,name='dW_slow')
137       CALL allocate_fields(nb_stage,f_dW_fast, field_t,type_real,llm+1,name='dW_fast')
138       CALL allocate_fields(nb_stage,f_dPhi_slow, field_t,type_real,llm+1,name='dPhi_slow')
139       CALL allocate_fields(nb_stage,f_dPhi_fast, field_t,type_real,llm+1,name='dPhi_fast')
140       f_dps => f_dps_slow(:,1)
141       f_du => f_du_slow(:,1)
142       f_dtheta_rhodz => f_dtheta_rhodz_slow(:,1)
143    END SELECT
144
145    SELECT CASE(scheme)
146    CASE(mlf)
147       ! Model state 2 time steps ago (MLF)
148       CALL allocate_field(f_psm2,field_t,type_real)
149       CALL allocate_field(f_massm2,field_t,type_real,llm)
150       CALL allocate_field(f_theta_rhodzm2,field_t,type_real,llm,nqdyn)
151       CALL allocate_field(f_um2,field_u,type_real,llm)
152    END SELECT
153
154    CALL init_theta2theta_rhodz
155    CALL init_dissip
156    CALL init_sponge
157    CALL init_observable
158    CALL init_guided
159    CALL init_advect_tracer
160    CALL init_check_conserve
161
162    CALL etat0(f_ps,f_mass,f_phis,f_theta_rhodz,f_u, f_geopot,f_W, f_q)
163
164    CALL transfert_request(f_phis,req_i0) 
165    CALL transfert_request(f_phis,req_i1) 
166
167    CALL init_message(f_ps,req_i0,req_ps0)
168    CALL init_message(f_mass,req_i0,req_mass0)
169    CALL init_message(f_theta_rhodz,req_i0,req_theta_rhodz0)
170    CALL init_message(f_u,req_e0_vect,req_u0)
171    CALL init_message(f_q,req_i0,req_q0)
172    CALL init_message(f_geopot,req_i0,req_geopot0)
173    CALL init_message(f_W,req_i0,req_W0)
174
175  END SUBROUTINE init_timeloop
176 
177  SUBROUTINE timeloop
178    USE abort_mod
179    USE dissip_gcm_mod
180    USE sponge_mod
181    USE observable_mod
182    USE etat0_mod
183    USE guided_mod
184    USE caldyn_mod
185    USE advect_tracer_mod
186    USE diagflux_mod
187    USE physics_mod
188    USE mpipara
189    USE transfert_mod
190    USE check_conserve_mod
191    USE xios_mod
192    USE output_field_mod
193    USE write_etat0_mod
194    USE restart_mod
195    USE checksum_mod
196    USE explicit_scheme_mod
197    USE hevi_scheme_mod
198    REAL(rstd),POINTER :: rhodz(:,:), mass(:,:), ps(:)
199
200    REAL(rstd) :: adv_over_out ! ratio itau_adv/itau_out
201    INTEGER :: ind, it,l
202    LOGICAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating mass fluxes in time
203    LOGICAL, PARAMETER :: check_rhodz=.FALSE.
204    INTEGER(kind=8) :: start_clock, stop_clock, rate_clock
205    INTEGER :: itau_sync   ! best iteration for synchronisation and ensure 1+1=2
206    INTEGER :: i
207   
208    LOGICAL,SAVE :: first_physic=.TRUE.
209    !$OMP THREADPRIVATE(first_physic)   
210
211    itau_sync=1
212    DO i=2,3*sync_it
213      IF (MOD(86400,INT(i*dt))==0 .AND. ABS((sync_it-itau_sync)*1./sync_it )/sync_it < (sync_it-itau_sync)*1./sync_it) itau_sync=i
214    ENDDO
215    IF (is_master) PRINT*,"Synchronize frontier every itau_sync =",itau_sync
216     
217    CALL switch_omp_distrib_level
218    CALL caldyn_BC(f_phis, f_geopot, f_wflux) ! set constant values in first/last interfaces
219
220    !$OMP BARRIER
221    DO ind=1,ndomain
222       IF (.NOT. assigned_domain(ind)) CYCLE
223       CALL swap_dimensions(ind)
224       CALL swap_geometry(ind)
225       rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
226       IF(caldyn_eta==eta_mass) THEN
227          CALL compute_rhodz(.TRUE., ps, rhodz, ondevice=.FALSE.) ! save rhodz for transport scheme before dynamics update ps
228       ELSE
229          DO l=ll_begin,ll_end
230             rhodz(:,l)=mass(:,l)
231          ENDDO
232       END IF
233    END DO
234    !$OMP BARRIER
235    fluxt_zero=.TRUE.
236
237    IF(positive_theta) CALL copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
238    IF(diagflux_on) THEN
239       adv_over_out = itau_adv*(1./itau_out)
240    ELSE
241       adv_over_out = 0.
242    END IF
243
244    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,itau0) 
245
246    Call trace_on
247
248    IF (xios_output) THEN ! we must call update_calendar before any XIOS output
249      IF (is_omp_master) CALL xios_update_calendar(1)
250    END IF
251   
252!    IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q)
253    IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W)
254   
255    CALL write_output_fields_basic(.TRUE., f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
256
257    !$OMP MASTER
258    CALL SYSTEM_CLOCK(start_clock, rate_clock)
259    !$OMP END MASTER   
260    call update_device_field(f_ps)
261    call update_device_field(f_mass)
262    CALL update_device_field(f_theta_rhodz)
263    CALL update_device_field(f_u)
264    CALL update_device_field(f_q)
265    CALL update_device_field(f_geopot)
266    CALL update_device_field(f_wflux)
267    CALL update_device_field(f_rhodz)
268    call reset_profiling() 
269
270
271    DO it=itau0+1,itau0+itaumax
272
273       CALL print_iteration(it, itau0, itaumax, start_clock, rate_clock)
274
275       CALL enter_profile(id_timeloop)
276
277       IF (xios_output) THEN
278
279          IF(it>itau0+1 .AND. is_omp_master) CALL xios_update_calendar(it-itau0)
280       ELSE
281          CALL update_time_counter(dt*it)
282       END IF
283
284       IF (it==itau0+1 .OR. MOD(it-1,itau_sync)==0) THEN
285          CALL send_message(f_ps,req_ps0)
286          CALL wait_message(req_ps0)
287          CALL send_message(f_mass,req_mass0)
288          CALL wait_message(req_mass0)
289          CALL send_message(f_theta_rhodz,req_theta_rhodz0) 
290          CALL wait_message(req_theta_rhodz0)
291          CALL send_message(f_u,req_u0)
292          CALL wait_message(req_u0)
293          CALL send_message(f_q,req_q0) 
294          CALL wait_message(req_q0)
295          IF(.NOT. hydrostatic) THEN
296             CALL send_message(f_geopot,req_geopot0)
297             CALL wait_message(req_geopot0)
298             CALL send_message(f_W,req_W0)
299             CALL wait_message(req_W0)
300          END IF
301       ENDIF
302
303!       CALL checksum(f_ps)
304!       CALL checksum(f_theta_rhodz)
305!       CALL checksum(f_u)
306!       CALL checksum(f_q)
307!       CALL checksum(f_mass)
308!       CALL checksum(f_geopot)
309!       CALL checksum(f_rhodz)
310!       CALL checksum(f_wflux)
311
312       CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q)
313
314       CALL enter_profile(id_dyn)
315       SELECT CASE(scheme_family)
316       CASE(explicit)
317          CALL abort_acc("explicit_scheme")
318          CALL explicit_scheme(it, fluxt_zero)
319       CASE(hevi)
320          CALL HEVI_scheme(it, fluxt_zero)
321       END SELECT
322       CALL exit_profile(id_dyn)
323
324!       CALL checksum(f_ps)
325!       CALL checksum(f_theta_rhodz)
326!       CALL checksum(f_u)
327!       CALL checksum(f_q)
328!       CALL checksum(f_mass)
329!       CALL checksum(f_geopot)
330!       CALL checksum(f_rhodz)
331!       CALL checksum(f_wflux)
332
333       CALL enter_profile(id_dissip)
334       IF (MOD(it,itau_dissip)==0) THEN
335         
336          IF(caldyn_eta==eta_mass) THEN
337             !ym flush ps
338             !$OMP BARRIER
339             DO ind=1,ndomain
340                IF (.NOT. assigned_domain(ind)) CYCLE
341                CALL swap_dimensions(ind)
342                CALL swap_geometry(ind)
343                mass=f_mass(ind); ps=f_ps(ind);
344                CALL compute_rhodz(.TRUE., ps, mass, ondevice=.TRUE.)
345             END DO
346          ENDIF
347         
348          CALL enter_profile(id_diags)
349          CALL check_conserve_detailed(it, AAM_dyn, &
350               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
351          CALL exit_profile(id_diags)
352       
353          CALL dissip(f_ps,f_mass,f_phis,f_geopot,f_theta_rhodz,f_u, f_dtheta_rhodz,f_du)
354         
355          CALL euler_scheme(.FALSE.)  ! update only u, theta
356          IF (iflag_sponge > 0) THEN
357             CALL abort_acc("iflag_sponge>0")
358             CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz)
359             CALL euler_scheme(.FALSE.)  ! update only u, theta
360          END IF
361
362          CALL enter_profile(id_diags)
363          CALL check_conserve_detailed(it, AAM_dissip, &
364               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
365          CALL exit_profile(id_diags)
366       END IF
367       CALL exit_profile(id_dissip)
368!       CALL checksum(f_ps)
369!       CALL checksum(f_theta_rhodz)
370!       CALL checksum(f_u)
371!       CALL checksum(f_q)
372!       CALL checksum(f_hfluxt)
373!       CALL checksum(f_wfluxt)
374!       CALL checksum(f_u)
375!       CALL checksum(f_rhodz)
376
377       CALL enter_profile(id_adv)
378       IF(MOD(it,itau_adv)==0) THEN
379          CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz, & ! update q and rhodz after RK step
380               adv_over_out, f_masst,f_qmasst,f_massfluxt, f_qfluxt)  ! accumulate mass and fluxes if diagflux_on
381          fluxt_zero=.TRUE. ! restart accumulation of hfluxt and qfluxt at next call to explicit_scheme / HEVI_scheme
382          ! At this point advect_tracer has obtained the halos of u and rhodz,
383          ! needed for correct computation of kinetic energy
384          IF(diagflux_on) CALL abort_acc("diagflux_on")
385          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)
386
387          DO ind=1,ndomain
388            IF (.NOT. assigned_domain(ind)) CYCLE
389            CALL swap_dimensions(ind)
390            CALL swap_geometry(ind)
391            rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
392            IF(caldyn_eta==eta_mass) THEN
393              CALL compute_rhodz(.TRUE., ps, rhodz, ondevice=.TRUE.) ! save rhodz for transport scheme before dynamics update ps
394            ELSE
395              DO l=ll_begin,ll_end
396                 rhodz(:,l)=mass(:,l)
397              ENDDO
398            END IF
399          END DO
400
401          IF(positive_theta) CALL abort_acc("positive_theta")
402          IF(positive_theta) CALL copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
403       END IF
404       CALL exit_profile(id_adv)
405!       CALL checksum(f_ps)
406!       CALL checksum(f_theta_rhodz)
407!       CALL checksum(f_u)
408!       CALL checksum(f_q)
409
410       CALL enter_profile(id_diags)
411!       IF (MOD(it,itau_physics)==0) THEN
412          CALL check_conserve_detailed(it, AAM_dyn, &
413            f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
414          CALL enter_profile(id_phys)
415          CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)       
416          CALL exit_profile(id_phys)
417          CALL check_conserve_detailed(it, AAM_phys, &
418               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
419          !$OMP MASTER
420          IF (first_physic) CALL SYSTEM_CLOCK(start_clock)
421          !$OMP END MASTER   
422          first_physic=.FALSE.
423!       END IF
424
425       IF (MOD(it,itau_check_conserv)==0) THEN
426          CALL update_host_field(f_ps)
427          CALL update_host_field(f_theta_rhodz)
428          CALL update_host_field(f_u)
429          CALL update_host_field(f_dps)
430          CALL update_host_field(f_q)
431          CALL check_conserve_detailed(it, AAM_dyn, &
432               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
433          CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,it) 
434       ENDIF
435
436       IF (mod(it,itau_out)==0 ) THEN
437          CALL transfert_request(f_u,req_e1_vect)
438          CALL update_host_field(f_ps)             
439          CALL update_host_field(f_mass)
440          CALL update_host_field(f_theta_rhodz)
441          CALL update_host_field(f_geopot)
442          CALL update_host_field(f_u)
443          CALL update_host_field(f_q)
444          CALL write_output_fields_basic(.FALSE.,f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
445       ENDIF
446       CALL exit_profile(id_diags)
447
448       CALL exit_profile(id_timeloop)
449    END DO
450   
451    CALL update_host_field(f_ps)
452    CALL update_host_field(f_theta_rhodz)
453    CALL update_host_field(f_u)
454    CALL update_host_field(f_q)
455    CALL update_host_field(f_geopot)
456
457!    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q)
458    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W) 
459
460    CALL update_host_field(f_dps)   
461    CALL check_conserve_detailed(it, AAM_dyn, &
462         f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
463    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,it) 
464   
465    !$OMP MASTER
466    CALL SYSTEM_CLOCK(stop_clock)
467   
468    IF (mpi_rank==0) THEN
469       PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock 
470    ENDIF
471    !$OMP END MASTER
472   
473    ! CONTAINS
474  END SUBROUTINE timeloop
475
476  SUBROUTINE print_iteration(it,itau0,itaumax,start_clock,rate_clock)
477    INTEGER :: it, itau0, itaumax, throughput
478    INTEGER(kind=8) :: start_clock, stop_clock, rate_clock
479    REAL :: per_step,total, elapsed
480   
481    IF(is_master) THEN
482       WRITE(*,'(A,I7,A,F14.1)') "It No :",it,"   t :",dt*it
483       IF(MOD(it,10)==0) THEN
484          CALL SYSTEM_CLOCK(stop_clock)
485          elapsed = (stop_clock-start_clock)*1./rate_clock
486          per_step = elapsed/(it-itau0)
487          throughput = INT(dt/per_step)
488          total = per_step*itaumax
489          WRITE(*,'(A,I5,A,F6.2,A,I6)') 'Time spent (s):',INT(elapsed), &
490               '  -- ms/step : ', 1000*per_step, &
491               '  -- Throughput :', throughput
492          WRITE(*,'(A,I5,A,I5)') 'Whole job (min) :', INT(total/60.), &
493               '  -- Completion in (min) : ', INT((total-elapsed)/60.)
494       END IF
495    END IF
496    IF(MOD(it,itau_prof)==0) CALL print_profile
497
498  END SUBROUTINE print_iteration
499
500  SUBROUTINE copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
501    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
502    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
503    INTEGER :: ind, iq
504    DO ind=1, ndomain
505       IF (.NOT. assigned_domain(ind)) CYCLE
506       CALL swap_dimensions(ind)
507       CALL swap_geometry(ind)
508       theta_rhodz=f_theta_rhodz(ind)
509       rhodz=f_rhodz(ind)
510       q=f_q(ind)
511       DO iq=1, nqdyn
512          q(:,:,iq)  = theta_rhodz(:,:,iq)/rhodz(:,:)
513       END DO
514    END DO
515  END SUBROUTINE copy_theta_to_q
516
517  SUBROUTINE copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
518    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
519    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
520    INTEGER :: ind, iq
521    DO ind=1, ndomain
522       IF (.NOT. assigned_domain(ind)) CYCLE
523       CALL swap_dimensions(ind)
524       CALL swap_geometry(ind)
525       theta_rhodz=f_theta_rhodz(ind)
526       rhodz=f_rhodz(ind)
527       q=f_q(ind)
528       DO iq=1,nqdyn
529          theta_rhodz(:,:,iq) = rhodz(:,:)*q(:,:,iq)
530       END DO
531    END DO
532  END SUBROUTINE copy_q_to_theta
533
534END MODULE timeloop_gcm_mod
Note: See TracBrowser for help on using the repository browser.