Skip to content

Commit 5606a78

Browse files
authored
Chem preprocess speedup (#851)
1 parent 3e6a84f commit 5606a78

File tree

2 files changed

+29
-5
lines changed

2 files changed

+29
-5
lines changed

src/common/m_chemistry.fpp

Lines changed: 28 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -9,7 +9,8 @@
99
module m_chemistry
1010

1111
use m_thermochem, only: &
12-
num_species, molecular_weights, get_temperature, get_net_production_rates
12+
num_species, molecular_weights, get_temperature, get_net_production_rates, &
13+
gas_constant, get_mixture_molecular_weight
1314

1415
use m_global_parameters
1516

@@ -58,6 +59,32 @@ contains
5859

5960
end subroutine s_compute_q_T_sf
6061

62+
subroutine s_compute_T_from_primitives(q_T_sf, q_prim_vf, bounds)
63+
64+
type(scalar_field), intent(inout) :: q_T_sf
65+
type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
66+
type(int_bounds_info), dimension(1:3), intent(in) :: bounds
67+
68+
integer :: x, y, z, i
69+
real(wp), dimension(num_species) :: Ys
70+
real(wp) :: mix_mol_weight
71+
72+
do z = bounds(3)%beg, bounds(3)%end
73+
do y = bounds(2)%beg, bounds(2)%end
74+
do x = bounds(1)%beg, bounds(1)%end
75+
!$acc loop seq
76+
do i = chemxb, chemxe
77+
Ys(i - chemxb + 1) = q_prim_vf(i)%sf(x, y, z)
78+
end do
79+
80+
call get_mixture_molecular_weight(Ys, mix_mol_weight)
81+
q_T_sf%sf(x, y, z) = q_prim_vf(E_idx)%sf(x, y, z)*mix_mol_weight/(gas_constant*q_prim_vf(1)%sf(x, y, z))
82+
end do
83+
end do
84+
end do
85+
86+
end subroutine s_compute_T_from_primitives
87+
6188
subroutine s_compute_chemistry_reaction_flux(rhs_vf, q_cons_qp, q_T_sf, q_prim_qp, bounds)
6289

6390
type(scalar_field), dimension(sys_size), intent(inout) :: rhs_vf

src/pre_process/m_initial_condition.fpp

Lines changed: 1 addition & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -175,9 +175,6 @@ contains
175175
!! primitive variables are converted to conservative ones.
176176
subroutine s_generate_initial_condition
177177

178-
! First, compute the temperature field from the conservative variables.
179-
if (chemistry) call s_compute_q_T_sf(q_T_sf, q_cons_vf, idwint)
180-
181178
! Converting the conservative variables to the primitive ones given
182179
! preexisting initial condition data files were read in on start-up
183180
if (old_ic) then
@@ -198,7 +195,7 @@ contains
198195
! Converting the primitive variables to the conservative ones
199196
call s_convert_primitive_to_conservative_variables(q_prim_vf, q_cons_vf)
200197

201-
if (chemistry) call s_compute_q_T_sf(q_T_sf, q_cons_vf, idwint)
198+
if (chemistry) call s_compute_T_from_primitives(q_T_sf, q_prim_vf, idwint)
202199

203200
if (qbmm .and. .not. polytropic) then
204201
!Initialize pb and mv

0 commit comments

Comments
 (0)