source: codes/icosagcm/trunk/src/time/timeloop_gcm.f90 @ 581

Last change on this file since 581 was 581, checked in by dubos, 7 years ago

trunk : upgrading to devel

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