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

Last change on this file since 953 was 953, checked in by adurocher, 5 years ago

trunk : GPU implementation with OpenACC ( merge from glcp.idris.fr )

File size: 17.1 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 :: 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
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=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    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    LOGICAL,SAVE :: first_physic=.TRUE.
206    !$OMP THREADPRIVATE(first_physic)   
207
208    CALL switch_omp_distrib_level
209    CALL caldyn_BC(f_phis, f_geopot, f_wflux) ! set constant values in first/last interfaces
210
211    !$OMP BARRIER
212    DO ind=1,ndomain
213       IF (.NOT. assigned_domain(ind)) CYCLE
214       CALL swap_dimensions(ind)
215       CALL swap_geometry(ind)
216       rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
217       IF(caldyn_eta==eta_mass) THEN
218          CALL compute_rhodz(.TRUE., ps, rhodz, ondevice=.FALSE.) ! save rhodz for transport scheme before dynamics update ps
219       ELSE
220          DO l=ll_begin,ll_end
221             rhodz(:,l)=mass(:,l)
222          ENDDO
223       END IF
224    END DO
225    !$OMP BARRIER
226    fluxt_zero=.TRUE.
227
228    IF(positive_theta) CALL copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
229    IF(diagflux_on) THEN
230       adv_over_out = itau_adv*(1./itau_out)
231    ELSE
232       adv_over_out = 0.
233    END IF
234
235    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,itau0) 
236
237    Call trace_on
238
239    IF (xios_output) THEN ! we must call update_calendar before any XIOS output
240      IF (is_omp_master) CALL xios_update_calendar(1)
241    END IF
242   
243!    IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q)
244    IF (write_start) CALL write_etat0(itau0,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W)
245   
246    CALL write_output_fields_basic(.TRUE., f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
247
248    !$OMP MASTER
249    CALL SYSTEM_CLOCK(start_clock, rate_clock)
250    !$OMP END MASTER   
251    call update_device_field(f_ps)
252    call update_device_field(f_mass)
253    CALL update_device_field(f_theta_rhodz)
254    CALL update_device_field(f_u)
255    CALL update_device_field(f_q)
256    CALL update_device_field(f_geopot)
257    CALL update_device_field(f_wflux)
258    CALL update_device_field(f_rhodz)
259
260
261    DO it=itau0+1,itau0+itaumax
262       IF (is_master) CALL print_iteration(it, itau0, itaumax, start_clock, rate_clock)
263
264       CALL enter_profile(id_timeloop)
265
266       IF (xios_output) THEN
267
268          IF(it>itau0+1 .AND. is_omp_master) CALL xios_update_calendar(it-itau0)
269       ELSE
270          CALL update_time_counter(dt*it)
271       END IF
272
273       IF (it==itau0+1 .OR. MOD(it,itau_sync)==0) THEN
274          CALL send_message(f_ps,req_ps0)
275          CALL wait_message(req_ps0)
276          CALL send_message(f_mass,req_mass0)
277          CALL wait_message(req_mass0)
278          CALL send_message(f_theta_rhodz,req_theta_rhodz0) 
279          CALL wait_message(req_theta_rhodz0)
280          CALL send_message(f_u,req_u0)
281          CALL wait_message(req_u0)
282          CALL send_message(f_q,req_q0) 
283          CALL wait_message(req_q0)
284          IF(.NOT. hydrostatic) THEN
285             CALL send_message(f_geopot,req_geopot0)
286             CALL wait_message(req_geopot0)
287             CALL send_message(f_W,req_W0)
288             CALL wait_message(req_W0)
289          END IF
290       ENDIF
291
292       CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q)
293
294       CALL enter_profile(id_dyn)
295       SELECT CASE(scheme_family)
296       CASE(explicit)
297          CALL abort_acc("explicit_scheme")
298          CALL explicit_scheme(it, fluxt_zero)
299       CASE(hevi)
300          CALL HEVI_scheme(it, fluxt_zero)
301       END SELECT
302       CALL exit_profile(id_dyn)
303
304       CALL enter_profile(id_dissip)
305       IF (MOD(it,itau_dissip)==0) THEN
306         
307          IF(caldyn_eta==eta_mass) THEN
308             !ym flush ps
309             !$OMP BARRIER
310             DO ind=1,ndomain
311                IF (.NOT. assigned_domain(ind)) CYCLE
312                CALL swap_dimensions(ind)
313                CALL swap_geometry(ind)
314                mass=f_mass(ind); ps=f_ps(ind);
315                CALL compute_rhodz(.TRUE., ps, mass, ondevice=.TRUE.)
316             END DO
317          ENDIF
318         
319          CALL enter_profile(id_diags)
320          CALL check_conserve_detailed(it, AAM_dyn, &
321               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
322          CALL exit_profile(id_diags)
323       
324          CALL dissip(f_ps,f_mass,f_phis,f_geopot,f_theta_rhodz,f_u, f_dtheta_rhodz,f_du)
325         
326          CALL euler_scheme(.FALSE.)  ! update only u, theta
327          IF (iflag_sponge > 0) THEN
328             CALL abort_acc("iflag_sponge>0")
329             CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz)
330             CALL euler_scheme(.FALSE.)  ! update only u, theta
331          END IF
332
333          CALL enter_profile(id_diags)
334          CALL check_conserve_detailed(it, AAM_dissip, &
335               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
336          CALL exit_profile(id_diags)
337       END IF
338       CALL exit_profile(id_dissip)
339
340       CALL enter_profile(id_adv)
341       IF(MOD(it,itau_adv)==0) THEN
342          CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz, & ! update q and rhodz after RK step
343               adv_over_out, f_masst,f_qmasst,f_massfluxt, f_qfluxt)  ! accumulate mass and fluxes if diagflux_on
344          fluxt_zero=.TRUE. ! restart accumulation of hfluxt and qfluxt at next call to explicit_scheme / HEVI_scheme
345          ! At this point advect_tracer has obtained the halos of u and rhodz,
346          ! needed for correct computation of kinetic energy
347          IF(diagflux_on) CALL abort_acc("diagflux_on")
348          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)
349
350          IF (check_rhodz) THEN ! check that rhodz is consistent with ps
351             DO ind=1,ndomain
352                IF (.NOT. assigned_domain(ind)) CYCLE
353                CALL swap_dimensions(ind)
354                CALL swap_geometry(ind)
355                rhodz=f_rhodz(ind); ps=f_ps(ind);
356                CALL compute_rhodz(.FALSE., ps, rhodz)   
357             END DO
358          ENDIF
359          IF(positive_theta) CALL abort_acc("positive_theta")
360          IF(positive_theta) CALL copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
361       END IF
362       CALL exit_profile(id_adv)
363
364       CALL enter_profile(id_diags)
365!       IF (MOD(it,itau_physics)==0) THEN
366          CALL check_conserve_detailed(it, AAM_dyn, &
367            f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
368          CALL enter_profile(id_phys)
369          CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)       
370          CALL exit_profile(id_phys)
371          CALL check_conserve_detailed(it, AAM_phys, &
372               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
373          !$OMP MASTER
374          IF (first_physic) CALL SYSTEM_CLOCK(start_clock)
375          !$OMP END MASTER   
376          first_physic=.FALSE.
377!       END IF
378
379       IF (MOD(it,itau_check_conserv)==0) THEN
380          CALL update_host_field(f_ps)
381          CALL update_host_field(f_theta_rhodz)
382          CALL update_host_field(f_u)
383          CALL update_host_field(f_dps)
384          CALL update_host_field(f_q)
385          CALL check_conserve_detailed(it, AAM_dyn, &
386               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
387          CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,it) 
388       ENDIF
389
390       IF (mod(it,itau_out)==0 ) THEN
391          CALL transfert_request(f_u,req_e1_vect)
392          CALL update_host_field(f_ps)             
393          CALL update_host_field(f_mass)
394          CALL update_host_field(f_theta_rhodz)
395          CALL update_host_field(f_geopot)
396          CALL update_host_field(f_u)
397          CALL update_host_field(f_q)
398          CALL write_output_fields_basic(.FALSE.,f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
399       ENDIF
400       CALL exit_profile(id_diags)
401
402       CALL exit_profile(id_timeloop)
403    END DO
404   
405    CALL update_host_field(f_ps)
406    CALL update_host_field(f_theta_rhodz)
407    CALL update_host_field(f_u)
408    CALL update_host_field(f_q)
409    CALL update_host_field(f_geopot)
410
411!    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q)
412    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q, f_geopot, f_W) 
413
414    CALL update_host_field(f_dps)   
415    CALL check_conserve_detailed(it, AAM_dyn, &
416         f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
417    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,f_q,it) 
418   
419    !$OMP MASTER
420    CALL SYSTEM_CLOCK(stop_clock)
421   
422    IF (mpi_rank==0) THEN
423       PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock 
424    ENDIF
425    !$OMP END MASTER
426   
427    ! CONTAINS
428  END SUBROUTINE timeloop
429
430  SUBROUTINE print_iteration(it,itau0,itaumax,start_clock,rate_clock)
431    INTEGER :: it, itau0, itaumax, throughput
432    INTEGER(kind=8) :: start_clock, stop_clock, rate_clock
433    REAL :: per_step,total, elapsed
434    WRITE(*,'(A,I7,A,F14.1)') "It No :",it,"   t :",dt*it
435    IF(MOD(it,10)==0) THEN
436       CALL SYSTEM_CLOCK(stop_clock)
437       elapsed = (stop_clock-start_clock)*1./rate_clock
438       per_step = elapsed/(it-itau0)
439       throughput = INT(dt/per_step)
440       total = per_step*itaumax
441       WRITE(*,'(A,I5,A,F6.2,A,I6)') 'Time spent (s):',INT(elapsed), &
442            '  -- ms/step : ', 1000*per_step, &
443            '  -- Throughput :', throughput
444       WRITE(*,'(A,I5,A,I5)') 'Whole job (min) :', INT(total/60.), &
445            '  -- Completion in (min) : ', INT((total-elapsed)/60.)
446    END IF
447    IF(MOD(it,itau_prof)==0) CALL print_profile
448
449  END SUBROUTINE print_iteration
450
451  SUBROUTINE copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
452    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
453    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
454    INTEGER :: ind, iq
455    DO ind=1, ndomain
456       IF (.NOT. assigned_domain(ind)) CYCLE
457       CALL swap_dimensions(ind)
458       CALL swap_geometry(ind)
459       theta_rhodz=f_theta_rhodz(ind)
460       rhodz=f_rhodz(ind)
461       q=f_q(ind)
462       DO iq=1, nqdyn
463          q(:,:,iq)  = theta_rhodz(:,:,iq)/rhodz(:,:)
464       END DO
465    END DO
466  END SUBROUTINE copy_theta_to_q
467
468  SUBROUTINE copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
469    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
470    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
471    INTEGER :: ind, iq
472    DO ind=1, ndomain
473       IF (.NOT. assigned_domain(ind)) CYCLE
474       CALL swap_dimensions(ind)
475       CALL swap_geometry(ind)
476       theta_rhodz=f_theta_rhodz(ind)
477       rhodz=f_rhodz(ind)
478       q=f_q(ind)
479       DO iq=1,nqdyn
480          theta_rhodz(:,:,iq) = rhodz(:,:)*q(:,:,iq)
481       END DO
482    END DO
483  END SUBROUTINE copy_q_to_theta
484
485END MODULE timeloop_gcm_mod
Note: See TracBrowser for help on using the repository browser.