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

Last change on this file since 396 was 387, checked in by dubos, 8 years ago

Infrastructure for multiple dynamical tracers - tested with JW06 and moist baroclinic wave

File size: 13.5 KB
RevLine 
[12]1MODULE timeloop_gcm_mod
[151]2  USE icosa
[360]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
[350]9  IMPLICIT NONE
[133]10  PRIVATE
[12]11 
[186]12  INTEGER, PARAMETER :: itau_sync=10
[377]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
[151]15
[360]16  PUBLIC :: init_timeloop, timeloop
[151]17
[12]18CONTAINS
19 
[151]20  SUBROUTINE init_timeloop
[360]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
[12]31
[360]32    CHARACTER(len=255) :: def
33   
34    IF (xios_output) itau_out=1
35    IF (.NOT. enable_io) itau_out=HUGE(itau_out)
[377]36
37    positive_theta=.FALSE.
38    CALL getin('positive_theta',positive_theta)
39    IF(positive_theta .AND. nqtot<1) THEN
40       PRINT *, 'nqtot must be >0 if positive_theta is .TRUE.'
41       STOP
42    END IF
43
[387]44    nqdyn = 1 ! one dynamical tracer = theta for the moment
45
[376]46    def='ARK2.3'
[360]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
[159]86    CALL allocate_field(f_phis,field_t,type_real,name='phis')
[360]87    ! Model state at current time step
[159]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')
[360]90    CALL allocate_field(f_rhodz,field_t,type_real,llm,name='rhodz')
[387]91    CALL allocate_field(f_theta_rhodz,field_t,type_real,llm,nqdyn,name='theta_rhodz')
[159]92    CALL allocate_field(f_u,field_u,type_real,llm,name='u')
[366]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')
[266]95    CALL allocate_field(f_q,field_t,type_real,llm,nqtot,'q')
[360]96    ! Mass fluxes
[159]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
[360]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
[159]105       CALL allocate_field(f_dps,field_t,type_real,name='dps')
[360]106       CALL allocate_field(f_dmass,field_t,type_real,llm, name='dmass')
[387]107       CALL allocate_field(f_dtheta_rhodz,field_t,type_real,llm,nqdyn,name='dtheta_rhodz')
[360]108       CALL allocate_field(f_du,field_u,type_real,llm,name='du')
109       ! Model state at previous time step (RK/MLF)
[159]110       CALL allocate_field(f_psm1,field_t,type_real,name='psm1')
[162]111       CALL allocate_field(f_massm1,field_t,type_real,llm, name='massm1')
[387]112       CALL allocate_field(f_theta_rhodzm1,field_t,type_real,llm,nqdyn,name='theta_rhodzm1')
[360]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')
[387]118       CALL allocate_fields(nb_stage,f_dtheta_rhodz_slow, field_t,type_real,llm,nqdyn,name='dtheta_rhodz_fast')
[360]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')
[366]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')
[360]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
[151]129
[360]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)
[387]135       CALL allocate_field(f_theta_rhodzm2,field_t,type_real,llm,nqdyn)
[360]136       CALL allocate_field(f_um2,field_u,type_real,llm)
[151]137    END SELECT
138
[295]139    CALL init_theta2theta_rhodz
[151]140    CALL init_dissip
[327]141    CALL init_sponge
[354]142    CALL init_observable
[151]143    CALL init_caldyn
144    CALL init_guided
145    CALL init_advect_tracer
146    CALL init_check_conserve
[186]147
[366]148    CALL etat0(f_ps,f_mass,f_phis,f_theta_rhodz,f_u, f_geopot,f_W, f_q)
[151]149
150    CALL transfert_request(f_phis,req_i0) 
151    CALL transfert_request(f_phis,req_i1) 
152    CALL writefield("phis",f_phis,once=.TRUE.)
153
154    CALL init_message(f_ps,req_i0,req_ps0)
[162]155    CALL init_message(f_mass,req_i0,req_mass0)
[151]156    CALL init_message(f_theta_rhodz,req_i0,req_theta_rhodz0)
157    CALL init_message(f_u,req_e0_vect,req_u0)
158    CALL init_message(f_q,req_i0,req_q0)
[377]159    CALL init_message(f_geopot,req_i0,req_geopot0)
160    CALL init_message(f_W,req_i0,req_W0)
[151]161
162  END SUBROUTINE init_timeloop
[12]163 
[151]164  SUBROUTINE timeloop
[360]165    USE dissip_gcm_mod
166    USE sponge_mod
167    USE observable_mod
168    USE etat0_mod
169    USE guided_mod
170    USE caldyn_mod
171    USE advect_tracer_mod
172    USE physics_mod
173    USE mpipara
174    USE transfert_mod
175    USE check_conserve_mod
176    USE xios_mod
177    USE output_field_mod
178    USE write_etat0_mod
179    USE checksum_mod
180    USE explicit_scheme_mod
181    USE hevi_scheme_mod
182    REAL(rstd),POINTER :: rhodz(:,:), mass(:,:), ps(:)
[12]183
[151]184    INTEGER :: ind
[360]185    INTEGER :: it,i,j,l,n,  stage
[151]186    LOGICAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating fluxes in time
[360]187    LOGICAL, PARAMETER :: check_rhodz=.FALSE.
188    INTEGER :: start_clock, stop_clock, rate_clock
[327]189    LOGICAL,SAVE :: first_physic=.TRUE.
[360]190    !$OMP THREADPRIVATE(first_physic)   
191
[295]192    CALL switch_omp_distrib_level
[350]193    CALL caldyn_BC(f_phis, f_geopot, f_wflux) ! set constant values in first/last interfaces
[132]194
[360]195    !$OMP BARRIER
196    DO ind=1,ndomain
197       IF (.NOT. assigned_domain(ind)) CYCLE
198       CALL swap_dimensions(ind)
199       CALL swap_geometry(ind)
200       rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
201       IF(caldyn_eta==eta_mass) THEN
202          CALL compute_rhodz(.TRUE., ps, rhodz) ! save rhodz for transport scheme before dynamics update ps
203       ELSE
204          DO l=ll_begin,ll_end
205             rhodz(:,l)=mass(:,l)
206          ENDDO
207       END IF
208    END DO
209    !$OMP BARRIER
210    fluxt_zero=.TRUE.
[186]211
[387]212    IF(positive_theta) CALL copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
[377]213
[360]214    !$OMP MASTER
215    CALL SYSTEM_CLOCK(start_clock)
[364]216    CALL SYSTEM_CLOCK(count_rate=rate_clock)
[360]217    !$OMP END MASTER   
[266]218
[360]219    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,itau0) 
[326]220
[360]221    CALL trace_on
[186]222
[360]223    DO it=itau0+1,itau0+itaumax
[186]224
[364]225       IF (is_master) CALL print_iteration(it, itau0, itaumax, start_clock, rate_clock)
226       IF (xios_output) THEN
227          CALL xios_update_calendar(it)
228       ELSE
229          CALL update_time_counter(dt*it)
230       END IF
[295]231
[360]232       IF (it==itau0+1 .OR. MOD(it,itau_sync)==0) THEN
233          CALL send_message(f_ps,req_ps0)
234          CALL wait_message(req_ps0)
235          CALL send_message(f_mass,req_mass0)
236          CALL wait_message(req_mass0)
237          CALL send_message(f_theta_rhodz,req_theta_rhodz0) 
238          CALL wait_message(req_theta_rhodz0) 
239          CALL send_message(f_u,req_u0)
240          CALL wait_message(req_u0)
241          CALL send_message(f_q,req_q0) 
[377]242          CALL wait_message(req_q0)
243          IF(.NOT. hydrostatic) THEN
244!             CALL send_message(f_geopot,req_geopot0)
245!             CALL wait_message(req_geopot0)
246!             CALL send_message(f_W,req_W0)
247!             CALL wait_message(req_W0)
248          END IF
[167]249       ENDIF
[327]250
[360]251       IF (mod(it,itau_out)==0 ) THEN
[365]252          CALL transfert_request(f_u,req_e1_vect)
[374]253          CALL write_output_fields_basic(f_ps, f_mass, f_geopot, f_u, f_W, f_q)
[327]254       ENDIF
[130]255
[360]256       CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q)
[326]257
[360]258       SELECT CASE(scheme_family)
259       CASE(explicit)
260          CALL explicit_scheme(it, fluxt_zero)
261       CASE(hevi)
262          CALL HEVI_scheme(it, fluxt_zero)
263       END SELECT
[364]264       
[360]265       IF (MOD(it,itau_dissip)==0) THEN
266         
267          IF(caldyn_eta==eta_mass) THEN
268             !ym flush ps
269             !$OMP BARRIER
270             DO ind=1,ndomain
271                IF (.NOT. assigned_domain(ind)) CYCLE
272                CALL swap_dimensions(ind)
273                CALL swap_geometry(ind)
274                mass=f_mass(ind); ps=f_ps(ind);
275                CALL compute_rhodz(.TRUE., ps, mass)
[162]276             END DO
[360]277          ENDIF
278         
[365]279          CALL check_conserve_detailed(it, AAM_dyn, &
280               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
281       
[360]282          CALL dissip(f_u,f_du,f_mass,f_phis, f_theta_rhodz,f_dtheta_rhodz)
283         
284          CALL euler_scheme(.FALSE.)  ! update only u, theta
285          IF (iflag_sponge > 0) THEN
286             CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz)
287             CALL euler_scheme(.FALSE.)  ! update only u, theta
[162]288          END IF
[365]289
290          CALL check_conserve_detailed(it, AAM_dissip, &
291               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
[360]292       END IF
[148]293       
[360]294       IF(MOD(it,itau_adv)==0) THEN
[364]295          CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz)  ! update q and rhodz after RK step
[360]296          fluxt_zero=.TRUE.
297          ! FIXME : check that rhodz is consistent with ps
298          IF (check_rhodz) THEN
299             DO ind=1,ndomain
300                IF (.NOT. assigned_domain(ind)) CYCLE
301                CALL swap_dimensions(ind)
302                CALL swap_geometry(ind)
303                rhodz=f_rhodz(ind); ps=f_ps(ind);
304                CALL compute_rhodz(.FALSE., ps, rhodz)   
305             END DO
306          ENDIF
[387]307          IF(positive_theta) CALL copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
[159]308       END IF
[360]309       
310       IF (MOD(it,itau_physics)==0) THEN
[365]311          CALL check_conserve_detailed(it, AAM_dyn, &
312            f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
[360]313          CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)       
[365]314          CALL check_conserve_detailed(it, AAM_phys, &
315               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
[360]316          !$OMP MASTER
317          IF (first_physic) CALL SYSTEM_CLOCK(start_clock)
318          !$OMP END MASTER   
319          first_physic=.FALSE.
[159]320       END IF
[365]321
322       IF (MOD(it,itau_check_conserv)==0) THEN
323          CALL check_conserve_detailed(it, AAM_dyn, &
324               f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
325          CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
326       ENDIF       
[360]327    END DO
328   
329    CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q) 
330   
[365]331    CALL check_conserve_detailed(it, AAM_dyn, &
332         f_ps,f_dps,f_u,f_theta_rhodz,f_phis)
[360]333    CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
334   
335    !$OMP MASTER
336    CALL SYSTEM_CLOCK(stop_clock)
337   
338    IF (mpi_rank==0) THEN
339       PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock 
340    ENDIF
341    !$OMP END MASTER
342   
343    ! CONTAINS
344  END SUBROUTINE timeloop
[364]345
346  SUBROUTINE print_iteration(it,itau0,itaumax,start_clock,rate_clock)
347    INTEGER :: it, itau0, itaumax, start_clock, stop_clock, rate_clock, throughput
348    REAL :: per_step,total, elapsed
[365]349    WRITE(*,'(A,I7,A,F14.1)') "It No :",it,"   t :",dt*it
[364]350    IF(MOD(it,10)==0) THEN
351       CALL SYSTEM_CLOCK(stop_clock)
352       elapsed = (stop_clock-start_clock)*1./rate_clock
353       per_step = elapsed/(it-itau0)
354       throughput = dt/per_step
355       total = per_step*(itaumax-itau0)
356       WRITE(*,'(A,I5,A,F6.2,A,I6)') 'Time spent (s):',INT(elapsed), &
357            '  -- ms/step : ', 1000*per_step, &
358            '  -- Throughput :', throughput
359       WRITE(*,'(A,I5,A,I5)') 'Whole job (min) :', INT(total/60.), &
360            '  -- Completion in (min) : ', INT((total-elapsed)/60.)
361    END IF
362  END SUBROUTINE print_iteration
363
[387]364  SUBROUTINE copy_theta_to_q(f_theta_rhodz,f_rhodz,f_q)
[377]365    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
[387]366    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
367    INTEGER :: ind, iq
[377]368    DO ind=1, ndomain
369       IF (.NOT. assigned_domain(ind)) CYCLE
370       CALL swap_dimensions(ind)
371       CALL swap_geometry(ind)
372       theta_rhodz=f_theta_rhodz(ind)
373       rhodz=f_rhodz(ind)
374       q=f_q(ind)
[387]375       DO iq=1, nqdyn
376          q(:,:,iq)  = theta_rhodz(:,:,iq)/rhodz(:,:)
377       END DO
[377]378    END DO
[387]379  END SUBROUTINE copy_theta_to_q
[377]380
[387]381  SUBROUTINE copy_q_to_theta(f_theta_rhodz,f_rhodz,f_q)
[377]382    TYPE(t_field),POINTER :: f_theta_rhodz(:),f_rhodz(:), f_q(:)
[387]383    REAL(rstd), POINTER :: theta_rhodz(:,:,:), rhodz(:,:), q(:,:,:)
384    INTEGER :: ind, iq
[377]385    DO ind=1, ndomain
386       IF (.NOT. assigned_domain(ind)) CYCLE
387       CALL swap_dimensions(ind)
388       CALL swap_geometry(ind)
389       theta_rhodz=f_theta_rhodz(ind)
390       rhodz=f_rhodz(ind)
391       q=f_q(ind)
[387]392       DO iq=1,nqdyn
393          theta_rhodz(:,:,iq) = rhodz(:,:)*q(:,:,iq)
394       END DO
[377]395    END DO
[387]396  END SUBROUTINE copy_q_to_theta
[377]397
[12]398END MODULE timeloop_gcm_mod
Note: See TracBrowser for help on using the repository browser.