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