source: codes/icosagcm/trunk/src/timeloop_gcm.f90 @ 415

Last change on this file since 415 was 415, checked in by ymipsl, 8 years ago

Writing restart file is not working well with intel mpi library. Call is descativated for now.

YM

File size: 13.7 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 checksum_mod
178    USE explicit_scheme_mod
179    USE hevi_scheme_mod
180    REAL(rstd),POINTER :: rhodz(:,:), mass(:,:), ps(:)
181
182    INTEGER :: ind
183    INTEGER :: it,i,j,l,n,  stage
184    LOGICAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating fluxes in time
185    LOGICAL, PARAMETER :: check_rhodz=.FALSE.
186    INTEGER :: start_clock, stop_clock, rate_clock
187    LOGICAL,SAVE :: first_physic=.TRUE.
188    !$OMP THREADPRIVATE(first_physic)   
189
190    CALL switch_omp_distrib_level
191    CALL caldyn_BC(f_phis, f_geopot, f_wflux) ! set constant values in first/last interfaces
192
193    !$OMP BARRIER
194    DO ind=1,ndomain
195       IF (.NOT. assigned_domain(ind)) CYCLE
196       CALL swap_dimensions(ind)
197       CALL swap_geometry(ind)
198       rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
199       IF(caldyn_eta==eta_mass) THEN
200          CALL compute_rhodz(.TRUE., ps, rhodz) ! save rhodz for transport scheme before dynamics update ps
201       ELSE
202          DO l=ll_begin,ll_end
203             rhodz(:,l)=mass(:,l)
204          ENDDO
205       END IF
206    END DO
207    !$OMP BARRIER
208    fluxt_zero=.TRUE.
209
210    IF(positive_theta) CALL copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
211
212    !$OMP MASTER
213    CALL SYSTEM_CLOCK(start_clock)
214    CALL SYSTEM_CLOCK(count_rate=rate_clock)
215    !$OMP END MASTER   
216
217    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,itau0) 
218
219    CALL trace_on
220
221    IF (xios_output) THEN ! we must call update_calendar before any XIOS output
222       CALL xios_update_calendar(1)
223    END IF
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    DO it=itau0+1,itau0+itaumax
227
228       IF (is_master) CALL print_iteration(it, itau0, itaumax, start_clock, rate_clock)
229
230       IF (xios_output) THEN
231          IF(it>itau0+1) CALL xios_update_calendar(it-itau0)
232       ELSE
233          CALL update_time_counter(dt*it)
234       END IF
235
236       IF (it==itau0+1 .OR. MOD(it,itau_sync)==0) THEN
237          CALL send_message(f_ps,req_ps0)
238          CALL wait_message(req_ps0)
239          CALL send_message(f_mass,req_mass0)
240          CALL wait_message(req_mass0)
241          CALL send_message(f_theta_rhodz,req_theta_rhodz0) 
242          CALL wait_message(req_theta_rhodz0) 
243          CALL send_message(f_u,req_u0)
244          CALL wait_message(req_u0)
245          CALL send_message(f_q,req_q0) 
246          CALL wait_message(req_q0)
247          IF(.NOT. hydrostatic) THEN
248!             CALL send_message(f_geopot,req_geopot0)
249!             CALL wait_message(req_geopot0)
250!             CALL send_message(f_W,req_W0)
251!             CALL wait_message(req_W0)
252          END IF
253       ENDIF
254
255       CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q)
256
257       SELECT CASE(scheme_family)
258       CASE(explicit)
259          CALL explicit_scheme(it, fluxt_zero)
260       CASE(hevi)
261          CALL HEVI_scheme(it, fluxt_zero)
262       END SELECT
263       
264       IF (MOD(it,itau_dissip)==0) THEN
265         
266          IF(caldyn_eta==eta_mass) THEN
267             !ym flush ps
268             !$OMP BARRIER
269             DO ind=1,ndomain
270                IF (.NOT. assigned_domain(ind)) CYCLE
271                CALL swap_dimensions(ind)
272                CALL swap_geometry(ind)
273                mass=f_mass(ind); ps=f_ps(ind);
274                CALL compute_rhodz(.TRUE., ps, mass)
275             END DO
276          ENDIF
277         
278          CALL check_conserve_detailed(it, AAM_dyn, &
279               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
280       
281          CALL dissip(f_u,f_du,f_mass,f_phis, f_theta_rhodz,f_dtheta_rhodz)
282         
283          CALL euler_scheme(.FALSE.)  ! update only u, theta
284          IF (iflag_sponge > 0) THEN
285             CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz)
286             CALL euler_scheme(.FALSE.)  ! update only u, theta
287          END IF
288
289          CALL check_conserve_detailed(it, AAM_dissip, &
290               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
291       END IF
292       
293       IF(MOD(it,itau_adv)==0) THEN
294          CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz)  ! update q and rhodz after RK step
295          fluxt_zero=.TRUE.
296          ! FIXME : check that rhodz is consistent with ps
297          IF (check_rhodz) THEN
298             DO ind=1,ndomain
299                IF (.NOT. assigned_domain(ind)) CYCLE
300                CALL swap_dimensions(ind)
301                CALL swap_geometry(ind)
302                rhodz=f_rhodz(ind); ps=f_ps(ind);
303                CALL compute_rhodz(.FALSE., ps, rhodz)   
304             END DO
305          ENDIF
306          IF(positive_theta) CALL copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
307       END IF
308       
309       IF (MOD(it,itau_physics)==0) THEN
310          CALL check_conserve_detailed(it, AAM_dyn, &
311            f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
312          CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)       
313          CALL check_conserve_detailed(it, AAM_phys, &
314               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
315          !$OMP MASTER
316          IF (first_physic) CALL SYSTEM_CLOCK(start_clock)
317          !$OMP END MASTER   
318          first_physic=.FALSE.
319       END IF
320
321       IF (MOD(it,itau_check_conserv)==0) THEN
322          CALL check_conserve_detailed(it, AAM_dyn, &
323               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
324          CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
325       ENDIF
326
327       IF (mod(it,itau_out)==0 ) THEN
328          CALL transfert_request(f_u,req_e1_vect)
329          CALL write_output_fields_basic(.FALSE.,f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
330       ENDIF
331
332    END DO
333   
334!    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q)
335   
336    CALL check_conserve_detailed(it, AAM_dyn, &
337         f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
338    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
339   
340    !$OMP MASTER
341    CALL SYSTEM_CLOCK(stop_clock)
342   
343    IF (mpi_rank==0) THEN
344       PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock 
345    ENDIF
346    !$OMP END MASTER
347   
348    ! CONTAINS
349  END SUBROUTINE timeloop
350
351  SUBROUTINE print_iteration(it,itau0,itaumax,start_clock,rate_clock)
352    INTEGER :: it, itau0, itaumax, start_clock, stop_clock, rate_clock, throughput
353    REAL :: per_step,total, elapsed
354    WRITE(*,'(A,I7,A,F14.1)') "It No :",it,"   t :",dt*it
355    IF(MOD(it,10)==0) THEN
356       CALL SYSTEM_CLOCK(stop_clock)
357       elapsed = (stop_clock-start_clock)*1./rate_clock
358       per_step = elapsed/(it-itau0)
359       throughput = dt/per_step
360       total = per_step*(itaumax-itau0)
361       WRITE(*,'(A,I5,A,F6.2,A,I6)') 'Time spent (s):',INT(elapsed), &
362            '  -- ms/step : ', 1000*per_step, &
363            '  -- Throughput :', throughput
364       WRITE(*,'(A,I5,A,I5)') 'Whole job (min) :', INT(total/60.), &
365            '  -- Completion in (min) : ', INT((total-elapsed)/60.)
366    END IF
367  END SUBROUTINE print_iteration
368
369  SUBROUTINE copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
370    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
371    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
372    INTEGER :: ind, iq
373    DO ind=1, ndomain
374       IF (.NOT. assigned_domain(ind)) CYCLE
375       CALL swap_dimensions(ind)
376       CALL swap_geometry(ind)
377       theta_rhodz=f_theta_rhodz(ind)
378       rhodz=f_rhodz(ind)
379       q=f_q(ind)
380       DO iq=1, nqdyn
381          q(:,:,iq)  = theta_rhodz(:,:,iq)/rhodz(:,:)
382       END DO
383    END DO
384  END SUBROUTINE copy_theta_to_q
385
386  SUBROUTINE copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
387    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
388    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
389    INTEGER :: ind, iq
390    DO ind=1, ndomain
391       IF (.NOT. assigned_domain(ind)) CYCLE
392       CALL swap_dimensions(ind)
393       CALL swap_geometry(ind)
394       theta_rhodz=f_theta_rhodz(ind)
395       rhodz=f_rhodz(ind)
396       q=f_q(ind)
397       DO iq=1,nqdyn
398          theta_rhodz(:,:,iq) = rhodz(:,:)*q(:,:,iq)
399       END DO
400    END DO
401  END SUBROUTINE copy_q_to_theta
402
403END MODULE timeloop_gcm_mod
Note: See TracBrowser for help on using the repository browser.