-
Notifications
You must be signed in to change notification settings - Fork 99
/
Copy pathm_time_steppers.fpp
1297 lines (1039 loc) · 45.4 KB
/
m_time_steppers.fpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!>
!! @file m_time_steppers.f90
!! @brief Contains module m_time_steppers
#:include 'macros.fpp'
!> @brief The following module features a variety of time-stepping schemes.
!! Currently, it includes the following Runge-Kutta (RK) algorithms:
!! 1) 1st Order TVD RK
!! 2) 2nd Order TVD RK
!! 3) 3rd Order TVD RK
!! where TVD designates a total-variation-diminishing time-stepper.
module m_time_steppers
use m_derived_types !< Definitions of the derived types
use m_global_parameters !< Definitions of the global parameters
use m_rhs !< Right-hane-side (RHS) evaluation procedures
use m_data_output !< Run-time info & solution data output procedures
use m_bubbles_EE !< Ensemble-averaged bubble dynamics routines
use m_bubbles_EL !< Lagrange bubble dynamics routines
use m_ibm
use m_hyperelastic
use m_mpi_proxy !< Message passing interface (MPI) module proxy
use m_boundary_conditions
use m_helper
use m_sim_helpers
use m_fftw
use m_nvtx
use m_thermochem, only: num_species
use m_body_forces
implicit none
type(vector_field), allocatable, dimension(:) :: q_cons_ts !<
!! Cell-average conservative variables at each time-stage (TS)
type(scalar_field), allocatable, dimension(:) :: q_prim_vf !<
!! Cell-average primitive variables at the current time-stage
type(scalar_field), allocatable, dimension(:) :: rhs_vf !<
!! Cell-average RHS variables at the current time-stage
type(vector_field), allocatable, dimension(:) :: rhs_ts_rkck
!! Cell-average RHS variables at each time-stage (TS)
!! Adaptive 4th/5th order Runge—Kutta–Cash–Karp (RKCK) time stepper
type(vector_field), allocatable, dimension(:) :: q_prim_ts !<
!! Cell-average primitive variables at consecutive TIMESTEPS
real(wp), allocatable, dimension(:, :, :, :, :) :: rhs_pb
type(scalar_field) :: q_T_sf !<
!! Cell-average temperature variables at the current time-stage
real(wp), allocatable, dimension(:, :, :, :, :) :: rhs_mv
real(wp), allocatable, dimension(:, :, :) :: max_dt
integer, private :: num_ts !<
!! Number of time stages in the time-stepping scheme
!$acc declare create(q_cons_ts, q_prim_vf, q_T_sf, rhs_vf, rhs_ts_rkck, q_prim_ts, rhs_mv, rhs_pb, max_dt)
contains
!> The computation of parameters, the allocation of memory,
!! the association of pointers and/or the execution of any
!! other procedures that are necessary to setup the module.
subroutine s_initialize_time_steppers_module
integer :: i, j !< Generic loop iterators
! Setting number of time-stages for selected time-stepping scheme
if (time_stepper == 1) then
num_ts = 1
elseif (any(time_stepper == (/2, 3, 4/))) then
num_ts = 2
end if
! Allocating the cell-average conservative variables
@:ALLOCATE(q_cons_ts(1:num_ts))
do i = 1, num_ts
@:ALLOCATE(q_cons_ts(i)%vf(1:sys_size))
end do
do i = 1, num_ts
do j = 1, sys_size
@:ALLOCATE(q_cons_ts(i)%vf(j)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
end do
@:ACC_SETUP_VFs(q_cons_ts(i))
end do
! Allocating the cell-average primitive ts variables
if (probe_wrt) then
@:ALLOCATE(q_prim_ts(0:3))
do i = 0, 3
@:ALLOCATE(q_prim_ts(i)%vf(1:sys_size))
end do
do i = 0, 3
do j = 1, sys_size
@:ALLOCATE(q_prim_ts(i)%vf(j)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
end do
end do
do i = 0, 3
@:ACC_SETUP_VFs(q_prim_ts(i))
end do
end if
! Allocating the cell-average primitive variables
@:ALLOCATE(q_prim_vf(1:sys_size))
do i = 1, adv_idx%end
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
if (bubbles_euler) then
do i = bub_idx%beg, bub_idx%end
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
if (adv_n) then
@:ALLOCATE(q_prim_vf(n_idx)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(n_idx))
end if
end if
if (mhd) then
do i = B_idx%beg, B_idx%end
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
end if
if (elasticity) then
do i = stress_idx%beg, stress_idx%end
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
end if
if (hyperelasticity) then
do i = xibeg, xiend + 1
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
end if
if (cont_damage) then
@:ALLOCATE(q_prim_vf(damage_idx)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(damage_idx))
end if
if (model_eqns == 3) then
do i = internalEnergies_idx%beg, internalEnergies_idx%end
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
end if
if (surface_tension) then
@:ALLOCATE(q_prim_vf(c_idx)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(c_idx))
end if
if (chemistry) then
do i = chemxb, chemxe
@:ALLOCATE(q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_prim_vf(i))
end do
@:ALLOCATE(q_T_sf%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end))
@:ACC_SETUP_SFs(q_T_sf)
end if
@:ALLOCATE(pb_ts(1:2))
!Initialize bubble variables pb and mv at all quadrature nodes for all R0 bins
if (qbmm .and. (.not. polytropic)) then
@:ALLOCATE(pb_ts(1)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(pb_ts(1))
@:ALLOCATE(pb_ts(2)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(pb_ts(2))
@:ALLOCATE(rhs_pb(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
else if (qbmm .and. polytropic) then
@:ALLOCATE(pb_ts(1)%sf(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(pb_ts(1))
@:ALLOCATE(pb_ts(2)%sf(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(pb_ts(2))
@:ALLOCATE(rhs_pb(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
end if
@:ALLOCATE(mv_ts(1:2))
if (qbmm .and. (.not. polytropic)) then
@:ALLOCATE(mv_ts(1)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(mv_ts(1))
@:ALLOCATE(mv_ts(2)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(mv_ts(2))
@:ALLOCATE(rhs_mv(idwbuff(1)%beg:idwbuff(1)%end, &
idwbuff(2)%beg:idwbuff(2)%end, &
idwbuff(3)%beg:idwbuff(3)%end, 1:nnode, 1:nb))
else if (qbmm .and. polytropic) then
@:ALLOCATE(mv_ts(1)%sf(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(mv_ts(1))
@:ALLOCATE(mv_ts(2)%sf(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
@:ACC_SETUP_SFs(mv_ts(2))
@:ALLOCATE(rhs_mv(idwbuff(1)%beg:idwbuff(1)%beg + 1, &
idwbuff(2)%beg:idwbuff(2)%beg + 1, &
idwbuff(3)%beg:idwbuff(3)%beg + 1, 1:nnode, 1:nb))
end if
! Allocating the cell-average RHS time-stages for adaptive RKCK stepper
if (bubbles_lagrange .and. time_stepper == 4) then
@:ALLOCATE(rhs_ts_rkck(1:num_ts_rkck))
do i = 1, num_ts_rkck
@:ALLOCATE(rhs_ts_rkck(i)%vf(1:sys_size))
end do
do i = 1, num_ts_rkck
do j = 1, sys_size
@:ALLOCATE(rhs_ts_rkck(i)%vf(j)%sf(0:m, 0:n, 0:p))
end do
@:ACC_SETUP_VFs(rhs_ts_rkck(i))
end do
else
! Allocating the cell-average RHS variables
@:ALLOCATE(rhs_vf(1:sys_size))
do i = 1, sys_size
@:ALLOCATE(rhs_vf(i)%sf(0:m, 0:n, 0:p))
@:ACC_SETUP_SFs(rhs_vf(i))
end do
end if
! Opening and writing the header of the run-time information file
if (proc_rank == 0 .and. run_time_info) then
call s_open_run_time_information_file()
end if
if (cfl_dt) then
@:ALLOCATE(max_dt(0:m, 0:n, 0:p))
end if
end subroutine s_initialize_time_steppers_module
!> 1st order TVD RK time-stepping algorithm
!! @param t_step Current time step
subroutine s_1st_order_tvd_rk(t_step, time_avg)
integer, intent(in) :: t_step
real(wp), intent(inout) :: time_avg
integer :: i, j, k, l, q !< Generic loop iterator
! Stage 1 of 1
call nvtxStartRange("TIMESTEP")
call s_compute_rhs(q_cons_ts(1)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(1)%sf, rhs_pb, mv_ts(1)%sf, rhs_mv, t_step, time_avg)
#ifdef DEBUG
print *, 'got rhs'
#endif
if (run_time_info) then
call s_write_run_time_information(q_prim_vf, t_step)
end if
#ifdef DEBUG
print *, 'wrote runtime info'
#endif
if (probe_wrt) then
call s_time_step_cycling(t_step)
end if
if (cfl_dt) then
if (mytime >= t_stop) return
else
if (t_step == t_step_stop) return
end if
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, stage=1)
call s_update_lagrange_tdv_rk(stage=1)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(1)%vf(i)%sf(j, k, l) = &
q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ dt*rhs_vf(i)%sf(j, k, l)
end do
end do
end do
end do
!Evolve pb and mv for non-polytropic qbmm
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(1)%sf(j, k, l, q, i) = &
pb_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_pb(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(1)%sf(j, k, l, q, i) = &
mv_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_mv(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, dt)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(1)%vf)
if (model_eqns == 3) call s_pressure_relaxation_procedure(q_cons_ts(1)%vf)
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(1)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf, pb_ts(1)%sf, mv_ts(1)%sf)
else
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf)
end if
end if
call nvtxEndRange
end subroutine s_1st_order_tvd_rk
!> 2nd order TVD RK time-stepping algorithm
!! @param t_step Current time-step
subroutine s_2nd_order_tvd_rk(t_step, time_avg)
integer, intent(in) :: t_step
real(wp), intent(inout) :: time_avg
integer :: i, j, k, l, q!< Generic loop iterator
real(wp) :: start, finish
! Stage 1 of 2
call cpu_time(start)
call nvtxStartRange("TIMESTEP")
call s_compute_rhs(q_cons_ts(1)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(1)%sf, rhs_pb, mv_ts(1)%sf, rhs_mv, t_step, time_avg)
if (run_time_info) then
call s_write_run_time_information(q_prim_vf, t_step)
end if
if (probe_wrt) then
call s_time_step_cycling(t_step)
end if
if (cfl_dt) then
if (mytime >= t_stop) return
else
if (t_step == t_step_stop) return
end if
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, stage=1)
call s_update_lagrange_tdv_rk(stage=1)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(2)%vf(i)%sf(j, k, l) = &
q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ dt*rhs_vf(i)%sf(j, k, l)
end do
end do
end do
end do
!Evolve pb and mv for non-polytropic qbmm
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(2)%sf(j, k, l, q, i) = &
pb_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_pb(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(2)%sf(j, k, l, q, i) = &
mv_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_mv(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, dt)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(2)%vf)
if (model_eqns == 3 .and. (.not. relax)) then
call s_pressure_relaxation_procedure(q_cons_ts(2)%vf)
end if
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(2)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf, pb_ts(2)%sf, mv_ts(2)%sf)
else
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf)
end if
end if
! Stage 2 of 2
call s_compute_rhs(q_cons_ts(2)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(2)%sf, rhs_pb, mv_ts(2)%sf, rhs_mv, t_step, time_avg)
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, stage=2)
call s_update_lagrange_tdv_rk(stage=2)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(1)%vf(i)%sf(j, k, l) = &
(q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ q_cons_ts(2)%vf(i)%sf(j, k, l) &
+ dt*rhs_vf(i)%sf(j, k, l))/2._wp
end do
end do
end do
end do
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(1)%sf(j, k, l, q, i) = &
(pb_ts(1)%sf(j, k, l, q, i) &
+ pb_ts(2)%sf(j, k, l, q, i) &
+ dt*rhs_pb(j, k, l, q, i))/2._wp
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(1)%sf(j, k, l, q, i) = &
(mv_ts(1)%sf(j, k, l, q, i) &
+ mv_ts(2)%sf(j, k, l, q, i) &
+ dt*rhs_mv(j, k, l, q, i))/2._wp
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, 2._wp*dt/3._wp)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(1)%vf)
if (model_eqns == 3 .and. (.not. relax)) then
call s_pressure_relaxation_procedure(q_cons_ts(1)%vf)
end if
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(1)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf, pb_ts(1)%sf, mv_ts(1)%sf)
else
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf)
end if
end if
call nvtxEndRange
call cpu_time(finish)
end subroutine s_2nd_order_tvd_rk
!> 3rd order TVD RK time-stepping algorithm
!! @param t_step Current time-step
subroutine s_3rd_order_tvd_rk(t_step, time_avg)
integer, intent(IN) :: t_step
real(wp), intent(INOUT) :: time_avg
integer :: i, j, k, l, q !< Generic loop iterator
real(wp) :: start, finish
! Stage 1 of 3
if (.not. adap_dt) then
call cpu_time(start)
call nvtxStartRange("TIMESTEP")
end if
call s_compute_rhs(q_cons_ts(1)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(1)%sf, rhs_pb, mv_ts(1)%sf, rhs_mv, t_step, time_avg)
if (run_time_info) then
call s_write_run_time_information(q_prim_vf, t_step)
end if
if (probe_wrt) then
call s_time_step_cycling(t_step)
end if
if (cfl_dt) then
if (mytime >= t_stop) return
else
if (t_step == t_step_stop) return
end if
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, stage=1)
call s_update_lagrange_tdv_rk(stage=1)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(2)%vf(i)%sf(j, k, l) = &
q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ dt*rhs_vf(i)%sf(j, k, l)
end do
end do
end do
end do
!Evolve pb and mv for non-polytropic qbmm
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(2)%sf(j, k, l, q, i) = &
pb_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_pb(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(2)%sf(j, k, l, q, i) = &
mv_ts(1)%sf(j, k, l, q, i) &
+ dt*rhs_mv(j, k, l, q, i)
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, dt)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(2)%vf)
if (model_eqns == 3 .and. (.not. relax)) then
call s_pressure_relaxation_procedure(q_cons_ts(2)%vf)
end if
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(2)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf, pb_ts(2)%sf, mv_ts(2)%sf)
else
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf)
end if
end if
! Stage 2 of 3
call s_compute_rhs(q_cons_ts(2)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(2)%sf, rhs_pb, mv_ts(2)%sf, rhs_mv, t_step, time_avg)
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, stage=2)
call s_update_lagrange_tdv_rk(stage=2)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(2)%vf(i)%sf(j, k, l) = &
(3._wp*q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ q_cons_ts(2)%vf(i)%sf(j, k, l) &
+ dt*rhs_vf(i)%sf(j, k, l))/4._wp
end do
end do
end do
end do
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(2)%sf(j, k, l, q, i) = &
(3._wp*pb_ts(1)%sf(j, k, l, q, i) &
+ pb_ts(2)%sf(j, k, l, q, i) &
+ dt*rhs_pb(j, k, l, q, i))/4._wp
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(2)%sf(j, k, l, q, i) = &
(3._wp*mv_ts(1)%sf(j, k, l, q, i) &
+ mv_ts(2)%sf(j, k, l, q, i) &
+ dt*rhs_mv(j, k, l, q, i))/4._wp
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, dt/4._wp)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(2)%vf)
if (model_eqns == 3 .and. (.not. relax)) then
call s_pressure_relaxation_procedure(q_cons_ts(2)%vf)
end if
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(2)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf, pb_ts(2)%sf, mv_ts(2)%sf)
else
call s_ibm_correct_state(q_cons_ts(2)%vf, q_prim_vf)
end if
end if
! Stage 3 of 3
call s_compute_rhs(q_cons_ts(2)%vf, q_T_sf, q_prim_vf, rhs_vf, pb_ts(2)%sf, rhs_pb, mv_ts(2)%sf, rhs_mv, t_step, time_avg)
if (bubbles_lagrange) then
call s_compute_EL_coupled_solver(q_cons_ts(2)%vf, q_prim_vf, rhs_vf, stage=3)
call s_update_lagrange_tdv_rk(stage=3)
end if
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, sys_size
do l = 0, p
do k = 0, n
do j = 0, m
q_cons_ts(1)%vf(i)%sf(j, k, l) = &
(q_cons_ts(1)%vf(i)%sf(j, k, l) &
+ 2._wp*q_cons_ts(2)%vf(i)%sf(j, k, l) &
+ 2._wp*dt*rhs_vf(i)%sf(j, k, l))/3._wp
end do
end do
end do
end do
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
pb_ts(1)%sf(j, k, l, q, i) = &
(pb_ts(1)%sf(j, k, l, q, i) &
+ 2._wp*pb_ts(2)%sf(j, k, l, q, i) &
+ 2._wp*dt*rhs_pb(j, k, l, q, i))/3._wp
end do
end do
end do
end do
end do
end if
if (qbmm .and. (.not. polytropic)) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
do l = 0, p
do k = 0, n
do j = 0, m
do q = 1, nnode
mv_ts(1)%sf(j, k, l, q, i) = &
(mv_ts(1)%sf(j, k, l, q, i) &
+ 2._wp*mv_ts(2)%sf(j, k, l, q, i) &
+ 2._wp*dt*rhs_mv(j, k, l, q, i))/3._wp
end do
end do
end do
end do
end do
end if
if (bodyForces) call s_apply_bodyforces(q_cons_ts(1)%vf, q_prim_vf, rhs_vf, 2._wp*dt/3._wp)
if (grid_geometry == 3) call s_apply_fourier_filter(q_cons_ts(1)%vf)
if (model_eqns == 3 .and. (.not. relax)) then
call s_pressure_relaxation_procedure(q_cons_ts(1)%vf)
end if
call nvtxStartRange("RHS-ELASTIC")
if (hyperelasticity) call s_hyperelastic_rmt_stress_update(q_cons_ts(1)%vf, q_prim_vf)
call nvtxEndRange
if (adv_n) call s_comp_alpha_from_n(q_cons_ts(1)%vf)
if (ib) then
if (qbmm .and. .not. polytropic) then
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf, pb_ts(1)%sf, mv_ts(1)%sf)
else
call s_ibm_correct_state(q_cons_ts(1)%vf, q_prim_vf)
end if
end if
if (.not. adap_dt) then
call nvtxEndRange
call cpu_time(finish)
time = time + (finish - start)
end if
end subroutine s_3rd_order_tvd_rk
!> Strang splitting scheme with 3rd order TVD RK time-stepping algorithm for
!! the flux term and adaptive time stepping algorithm for
!! the source term
!! @param t_step Current time-step
subroutine s_strang_splitting(t_step, time_avg)
integer, intent(in) :: t_step
real(wp), intent(inout) :: time_avg
real(wp) :: start, finish
call cpu_time(start)
call nvtxStartRange("TIMESTEP")
! Stage 1 of 3
call s_adaptive_dt_bubble(t_step)
! Stage 2 of 3
call s_3rd_order_tvd_rk(t_step, time_avg)
! Stage 3 of 3
call s_adaptive_dt_bubble(t_step)
call nvtxEndRange
call cpu_time(finish)
time = time + (finish - start)
end subroutine s_strang_splitting
!> Bubble source part in Strang operator splitting scheme
!! @param t_step Current time-step
subroutine s_adaptive_dt_bubble(t_step)
integer, intent(in) :: t_step
type(vector_field) :: gm_alpha_qp
call s_convert_conservative_to_primitive_variables( &
q_cons_ts(1)%vf, &
q_T_sf, &
q_prim_vf, &
idwint, &
gm_alpha_qp%vf)
call s_compute_bubble_EE_source(q_cons_ts(1)%vf, q_prim_vf, t_step, rhs_vf)
call s_comp_alpha_from_n(q_cons_ts(1)%vf)
end subroutine s_adaptive_dt_bubble
subroutine s_compute_dt()
real(wp) :: rho !< Cell-avg. density
real(wp), dimension(num_vels) :: vel !< Cell-avg. velocity
real(wp) :: vel_sum !< Cell-avg. velocity sum
real(wp) :: pres !< Cell-avg. pressure
real(wp), dimension(num_fluids) :: alpha !< Cell-avg. volume fraction
real(wp) :: gamma !< Cell-avg. sp. heat ratio
real(wp) :: pi_inf !< Cell-avg. liquid stiffness function
real(wp) :: c !< Cell-avg. sound speed
real(wp) :: H !< Cell-avg. enthalpy
real(wp), dimension(2) :: Re !< Cell-avg. Reynolds numbers
type(vector_field) :: gm_alpha_qp
real(wp) :: dt_local
integer :: j, k, l !< Generic loop iterators
call s_convert_conservative_to_primitive_variables( &
q_cons_ts(1)%vf, &
q_T_sf, &
q_prim_vf, &
idwint, &
gm_alpha_qp%vf)
!$acc parallel loop collapse(3) gang vector default(present) private(vel, alpha, Re)
do l = 0, p
do k = 0, n
do j = 0, m
call s_compute_enthalpy(q_prim_vf, pres, rho, gamma, pi_inf, Re, H, alpha, vel, vel_sum, j, k, l)
! Compute mixture sound speed
call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, H, alpha, vel_sum, 0._wp, c)
call s_compute_dt_from_cfl(vel, c, max_dt, rho, Re, j, k, l)
end do
end do
end do
!$acc kernels
dt_local = minval(max_dt)
!$acc end kernels
if (num_procs == 1) then
dt = dt_local
else
call s_mpi_allreduce_min(dt_local, dt)
end if
!$acc update device(dt)
end subroutine s_compute_dt
!> This subroutine applies the body forces source term at each
!! Runge-Kutta stage
subroutine s_apply_bodyforces(q_cons_vf, q_prim_vf, rhs_vf, ldt)
type(scalar_field), dimension(1:sys_size), intent(inout) :: q_cons_vf
type(scalar_field), dimension(1:sys_size), intent(in) :: q_prim_vf
type(scalar_field), dimension(1:sys_size), intent(inout) :: rhs_vf
real(wp), intent(in) :: ldt !< local dt
integer :: i, j, k, l