Skip to content

Commit aa8cb5c

Browse files
committed
convective_bdy_weight cleanup
This should restore the essential behavior of the previous control, except that it no longer applies to the "newly non-convective" regions. It only adds resolution at the current locations of convective boundaries.
1 parent a7b64d4 commit aa8cb5c

File tree

2 files changed

+27
-15
lines changed

2 files changed

+27
-15
lines changed

star/defaults/controls.defaults

Lines changed: 1 addition & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -6270,10 +6270,9 @@
62706270
!### convective_bdy_min_dt_yrs
62716271

62726272
! Mesh function to enhance resolution near convective boundaries
6273-
! including regions that are newly nonconvective because of moving boundary.
62746273

62756274
convective_bdy_weight = 0
6276-
convective_bdy_dq_limit = 1d-4
6275+
convective_bdy_dq_limit = 3d-5
62776276
convective_bdy_min_dt_yrs = 1d-3
62786277

62796278
! max_rel_delta_IE_for_mesh_total_energy_balance

star/private/mesh_functions.f90

Lines changed: 26 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -177,7 +177,7 @@ subroutine set_mesh_function_data( &
177177
i = i+1; names(i) = 'omega_function'
178178
end if
179179
if (s% convective_bdy_weight > 0) then
180-
i = i+1; names(i) = 'newly_nonconv'
180+
i = i+1; names(i) = 'convective_bdy'
181181
end if
182182
do k=1,num_xa_function
183183
if (do_mass_function(s, s% xa_function_species(k), s% xa_function_weight(k), j)) then
@@ -267,9 +267,8 @@ subroutine set_mesh_function_data( &
267267
vals(k,i) = s% omega_function_weight*log10(max(1d-99,abs(s% omega(k))))
268268
end do
269269

270-
else if (names(i) == 'newly_nonconv') then
271-
!call do_newly_nonconvective(i)
272-
call do_mix_type_change(i)
270+
else if (names(i) == 'convective_bdy') then
271+
call do_conv_bdy(i)
273272

274273
else
275274
do k=1,num_xa_function
@@ -283,21 +282,35 @@ subroutine set_mesh_function_data( &
283282
contains
284283

285284

286-
subroutine do_mix_type_change(i)
285+
subroutine do_conv_bdy(i)
287286
integer, intent(in) :: i
288287
integer :: k, j
288+
289+
vals(1:nz,i) = 0
290+
289291
if (s% dt < s% convective_bdy_min_dt_yrs*secyer) return
290-
do k=1,nz
292+
293+
!do k=1,nz
294+
! if (s% dq(k) < s% convective_bdy_dq_limit) cycle
295+
! if (s% cz_bdy_dq(k) /= 0) vals(k,i) = s% convective_bdy_weight
296+
!end do
297+
298+
do j = 1, s% num_conv_boundaries
299+
k = s% conv_bdy_loc(j)
291300
if (s% dq(k) < s% convective_bdy_dq_limit) cycle
292-
if (s% cz_bdy_dq(k) /= 0) vals(k,i) = s% convective_bdy_weight
293-
end do
294-
do j=1,5
295-
do k=2,nz
296-
vals(k,i) = vals(k,i) + vals(k-1,i)
297-
end do
301+
vals(k,i) = s% convective_bdy_weight
302+
if (k < nz .and. s% top_conv_bdy(i)) then
303+
vals(k+1,i) = s% convective_bdy_weight
304+
else if (k > 1 .and. .not. s% top_conv_bdy(i)) then
305+
vals(k-1,i) = s% convective_bdy_weight
306+
end if
298307
end do
299308

300-
end subroutine do_mix_type_change
309+
do k=2,nz
310+
vals(k,i) = vals(k,i) + vals(k-1,i)
311+
end do
312+
313+
end subroutine do_conv_bdy
301314

302315
subroutine do1_xa_function(k,i)
303316
integer, intent(in) :: k,i

0 commit comments

Comments
 (0)