Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010-2021 The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 : module turb_support
21 :
22 : use star_private_def
23 : use const_def, only: dp, crad, no_mixing
24 : use num_lib
25 : use utils_lib
26 : use auto_diff_support
27 : use star_utils
28 : use turb
29 :
30 : implicit none
31 :
32 : private
33 : public :: get_gradT
34 : public :: do1_mlt_eval
35 : public :: Get_results
36 :
37 : contains
38 :
39 : !> Determines if it is safe (physically) to use TDC instead of MLT.
40 : !!
41 : !! Currently we only know we have to fall back to MLT in cells that get touched
42 : !! by adjust_mass, because there the convection speeds at the start of the
43 : !! step can be badly out of whack.
44 : !!
45 : !! @param s star pointer
46 : !! @param k face index
47 : !! @param fallback False if we can use TDC, True if we can fall back to MLT.
48 65210 : logical function check_if_must_fall_back_to_MLT(s, k) result(fallback)
49 : type (star_info), pointer :: s
50 : integer, intent(in) :: k
51 :
52 65210 : fallback = .false.
53 0 : if (abs(s%mstar_dot) > 1d-99 .and. k < s% k_const_mass) then
54 : fallback = .true.
55 : end if
56 : end function check_if_must_fall_back_to_MLT
57 :
58 0 : subroutine get_gradT(s, MLT_option, & ! used to create models
59 : r, L, T, P, opacity, rho, chiRho, chiT, Cp, gradr, grada, scale_height, &
60 : iso, XH1, cgrav, m, gradL_composition_term, mixing_length_alpha, &
61 : mixing_type, gradT, Y_face, conv_vel, D, Gamma, ierr)
62 : type (star_info), pointer :: s
63 : character (len=*), intent(in) :: MLT_option
64 : real(dp), intent(in) :: &
65 : r, L, T, P, opacity, rho, chiRho, chiT, Cp, gradr, grada, scale_height, &
66 : XH1, cgrav, m, gradL_composition_term, mixing_length_alpha
67 : integer, intent(in) :: iso
68 : real(dp), intent(out) :: gradT, Y_face, conv_vel, D, Gamma
69 : integer, intent(out) :: mixing_type, ierr
70 : type(auto_diff_real_star_order1) :: &
71 : gradr_ad, grada_ad, scale_height_ad, gradT_ad, Y_face_ad, mlt_vc_ad, D_ad, &
72 : Gamma_ad, r_ad, L_ad, T_ad, P_ad, opacity_ad, rho_ad, dV_ad, chiRho_ad, chiT_ad, Cp_ad
73 0 : ierr = 0
74 0 : r_ad = r
75 0 : L_ad = L
76 0 : T_ad = T
77 0 : P_ad = P
78 0 : opacity_ad = opacity
79 0 : rho_ad = rho
80 0 : dV_ad = 0d0
81 0 : chiRho_ad = chiRho
82 0 : chiT_ad = chiT
83 0 : Cp_ad = Cp
84 0 : gradr_ad = gradr
85 0 : grada_ad = grada
86 0 : scale_height_ad = scale_height
87 : call Get_results(s, 0, MLT_option, &
88 : r_ad, L_ad, T_ad, P_ad, opacity_ad, rho_ad, dV_ad, chiRho_ad, chiT_ad, Cp_ad, &
89 : gradr_ad, grada_ad, scale_height_ad, &
90 : iso, XH1, cgrav, m, gradL_composition_term, mixing_length_alpha, &
91 : s% alpha_semiconvection, s% thermohaline_coeff, &
92 0 : mixing_type, gradT_ad, Y_face_ad, mlt_vc_ad, D_ad, Gamma_ad, ierr)
93 0 : gradT = gradT_ad%val
94 0 : Y_face = Y_face_ad%val
95 0 : conv_vel = mlt_vc_ad%val
96 0 : D = D_ad%val
97 0 : Gamma = Gamma_ad%val
98 0 : end subroutine get_gradT
99 :
100 :
101 79994 : subroutine do1_mlt_eval( &
102 : s, k, MLT_option, gradL_composition_term, &
103 : gradr_in, grada, scale_height, mixing_length_alpha, &
104 : mixing_type, gradT, Y_face, mlt_vc, D, Gamma, ierr)
105 : use chem_def, only: ih1
106 : use const_def, only: ln10
107 : use starspots, only: starspot_tweak_gradr
108 : type (star_info), pointer :: s
109 : integer, intent(in) :: k
110 : character (len=*), intent(in) :: MLT_option
111 : type(auto_diff_real_star_order1), intent(in) :: gradr_in, grada, scale_height
112 : real(dp), intent(in) :: gradL_composition_term, mixing_length_alpha
113 : integer, intent(out) :: mixing_type
114 : type(auto_diff_real_star_order1), intent(out) :: &
115 : gradT, Y_face, mlt_vc, D, Gamma
116 : integer, intent(out) :: ierr
117 :
118 79994 : real(dp) :: cgrav, m, XH1, P_theta, L_theta
119 : integer :: iso
120 : type(auto_diff_real_star_order1) :: gradr, r, L, T, P, opacity, rho, dV, chiRho, chiT, Cp, rho_start
121 : include 'formats'
122 79994 : ierr = 0
123 :
124 :
125 79994 : P = get_Peos_face(s,k) ! if u_flag, should this be P_face_ad? (time centered in riemann)
126 79994 : if (s% include_mlt_in_velocity_time_centering) then
127 : ! could be cleaner with a wrapper for time_centered P and L
128 : if (s% using_velocity_time_centering .and. &
129 0 : s% include_P_in_velocity_time_centering .and. &
130 : s% lnT(k)/ln10 <= s% max_logT_for_include_P_and_L_in_velocity_time_centering) then
131 0 : P_theta = s% P_theta_for_velocity_time_centering
132 : else
133 0 : P_theta = 1d0
134 : end if
135 : ! consder building a wrapper : wrap_opt_time_center_L_00(s,k)
136 : if (s% using_velocity_time_centering .and. &
137 0 : s% include_L_in_velocity_time_centering .and. &
138 : s% lnT(k)/ln10 <= s% max_logT_for_include_P_and_L_in_velocity_time_centering) then
139 0 : L_theta = s% L_theta_for_velocity_time_centering
140 : else
141 0 : L_theta = 1d0
142 : end if
143 0 : L = L_theta*wrap_L_00(s, k) + (1d0 - L_theta)*s% L_start(k)
144 0 : P = P_theta*P + (1d0-P_theta)*s% Peos_face_start(k)
145 0 : r = wrap_opt_time_center_r_00(s,k)
146 : else
147 79994 : L = wrap_L_00(s,k)
148 79994 : r = wrap_r_00(s,k)
149 : end if
150 79994 : gradr = gradr_in
151 79994 : cgrav = s% cgrav(k)
152 79994 : m = s% m_grav(k)
153 79994 : T = get_T_face(s,k)
154 79994 : opacity = get_kap_face(s,k)
155 79994 : rho = get_rho_face(s,k)
156 79994 : rho_start = get_rho_start_face(s,k)
157 79994 : dV = 1d0/rho - 1d0/rho_start ! both variables are face wrapped.
158 79994 : chiRho = get_ChiRho_face(s,k)
159 79994 : chiT = get_ChiT_face(s,k)
160 79994 : Cp = get_Cp_face(s,k)
161 79994 : iso = s% dominant_iso_for_thermohaline(k)
162 79994 : XH1 = s% xa(s% net_iso(ih1),k)
163 :
164 79994 : if (s% use_other_mlt_results) then
165 : call s% other_mlt_results(s% id, k, MLT_option, &
166 : r, L, T, P, opacity, rho, chiRho, chiT, Cp, gradr, grada, scale_height, &
167 : iso, XH1, cgrav, m, gradL_composition_term, mixing_length_alpha, &
168 : s% alpha_semiconvection, s% thermohaline_coeff, &
169 0 : mixing_type, gradT, Y_face, mlt_vc, D, Gamma, ierr)
170 : else
171 : ! starspot YREC routine
172 79994 : if (s% do_starspots) then
173 : !dV = 0d0 ! dV = 1/rho - 1/rho_start and we assume rho = rho_start.
174 0 : call starspot_tweak_gradr(s, P, gradr_in, gradr)
175 : end if
176 : call Get_results(s, k, MLT_option, &
177 : r, L, T, P, opacity, rho, dV, chiRho, chiT, Cp, gradr, grada, scale_height, &
178 : iso, XH1, cgrav, m, gradL_composition_term, mixing_length_alpha, &
179 : s% alpha_semiconvection, s% thermohaline_coeff, &
180 79994 : mixing_type, gradT, Y_face, mlt_vc, D, Gamma, ierr)
181 : end if
182 :
183 79994 : end subroutine do1_mlt_eval
184 :
185 :
186 79994 : subroutine Get_results(s, k, MLT_option, & ! NOTE: k=0 is a valid arg
187 : r, L, T, P, opacity, rho, dV, chiRho, chiT, Cp, gradr, grada, scale_height, &
188 : iso, XH1, cgrav, m, gradL_composition_term, mixing_length_alpha, &
189 : alpha_semiconvection, thermohaline_coeff, &
190 : mixing_type, gradT, Y_face, conv_vel, D, Gamma, ierr)
191 79994 : use star_utils
192 : use tdc_hydro, only: compute_tdc_Eq_cell, compute_tdc_Eq_div_w_face
193 : type (star_info), pointer :: s
194 : integer, intent(in) :: k
195 : character (len=*), intent(in) :: MLT_option
196 : type(auto_diff_real_star_order1), intent(in) :: &
197 : r, L, T, P, opacity, rho, dV, chiRho, chiT, Cp, gradr, grada, scale_height
198 : integer, intent(in) :: iso
199 : real(dp), intent(in) :: &
200 : XH1, cgrav, m, gradL_composition_term, &
201 : mixing_length_alpha, alpha_semiconvection, thermohaline_coeff
202 : integer, intent(out) :: mixing_type
203 : type(auto_diff_real_star_order1), intent(out) :: gradT, Y_face, conv_vel, D, Gamma
204 : integer, intent(out) :: ierr
205 :
206 : type(auto_diff_real_star_order1) :: Pr, Pg, grav, Lambda, gradL, beta
207 79994 : real(dp) :: conv_vel_start, scale, max_conv_vel
208 :
209 : ! these are used by use_superad_reduction
210 79994 : real(dp) :: Gamma_limit, scale_value1, scale_value2, diff_grads_limit, reduction_limit, lambda_limit
211 : type(auto_diff_real_star_order1) :: Lrad_div_Ledd, Gamma_inv_threshold, Gamma_factor, alfa0, &
212 : diff_grads_factor, Gamma_term, exp_limit, grad_scale, gradr_scaled, Eq_div_w, check_Eq, mlt_Pturb, Ptot
213 : logical :: test_partials, using_TDC
214 : logical, parameter :: report = .false.
215 : include 'formats'
216 :
217 : ! check if this particular k can be done with TDC
218 79994 : using_TDC = .false.
219 79994 : if (s% MLT_option == 'TDC') using_TDC = .true.
220 79994 : if (.not. s% have_mlt_vc) using_TDC = .false.
221 79994 : if (k <= 0 .or. s%dt <= 0d0) using_TDC = .false.
222 140276 : if (using_TDC) using_TDC = .not. check_if_must_fall_back_to_MLT(s, k)
223 :
224 : ! Pre-calculate some things.
225 79994 : Eq_div_w = 0d0
226 79994 : if ((s% v_flag .or. s% u_flag) .and. k > 0 ) then ! only include Eq_div_w if v_flag or u_flag is true.
227 0 : if (using_TDC .and. s% alpha_TDC_DampM > 0) then
228 0 : if (s% mlt_vc(k) > 0) then ! calculate using mlt_vc from current timestep.
229 0 : check_Eq = compute_tdc_Eq_div_w_face(s, k, ierr)
230 0 : Eq_div_w = check_Eq
231 : end if
232 : end if
233 : end if
234 :
235 : ! Wrap Pturb into P
236 79994 : if (s% okay_to_set_mlt_vc .and. s% include_mlt_Pturb_in_thermodynamic_gradients .and. k > 0) then
237 0 : mlt_Pturb = s% mlt_Pturb_factor*pow2(s% mlt_vc_old(k))*get_rho_face(s,k)/3d0
238 0 : Ptot = P + mlt_Pturb
239 : else
240 79994 : Ptot = P
241 : end if
242 :
243 79994 : Pr = crad*pow4(T)/3d0
244 79994 : Pg = Ptot - Pr
245 79994 : beta = Pg / Ptot
246 79994 : Lambda = mixing_length_alpha*scale_height
247 :
248 79994 : if (k == 0) then
249 0 : grav = cgrav*m/pow2(r)
250 : else
251 79994 : grav = cgrav*m/pow2(r) !try replacing with wrap_geff_face(s,k)
252 : end if
253 :
254 79994 : if (s% use_Ledoux_criterion) then
255 0 : gradL = grada + gradL_composition_term ! Ledoux temperature gradient
256 : else
257 : gradL = grada
258 : end if
259 :
260 : ! maximum convection velocity.
261 79994 : if (k > 0) then
262 79994 : if (s% q(k) <= s% max_conv_vel_div_csound_maxq) then
263 0 : max_conv_vel = s% csound_face(k) * s% max_conv_vel_div_csound
264 : else
265 79994 : max_conv_vel = 1d99
266 : end if
267 : else ! if k == 0
268 0 : max_conv_vel = 1d99
269 : end if
270 :
271 :
272 : ! Initialize with no mixing
273 79994 : mixing_type = no_mixing
274 79994 : gradT = gradr
275 79994 : Y_face = gradT - gradL
276 79994 : conv_vel = 0d0
277 79994 : D = 0d0
278 79994 : Gamma = 0d0
279 79994 : if (k /= 0) s% superad_reduction_factor(k) = 1d0
280 :
281 : ! Bail if we asked for no mixing, or if parameters are bad.
282 : if (MLT_option == 'none' .or. beta < 1d-10 .or. mixing_length_alpha <= 0d0 .or. &
283 : opacity%val < 1d-10 .or. P%val < 1d-20 .or. T%val < 1d-10 .or. Rho%val < 1d-20 &
284 79994 : .or. m < 1d-10 .or. r%val < 1d-10 .or. cgrav < 1d-10) return
285 :
286 : !test_partials = (k == s% solver_test_partials_k)
287 79994 : test_partials = .false.
288 79994 : ierr = 0
289 79994 : if (k > 0) then
290 79994 : s% tdc_num_iters(k) = 0
291 : end if
292 :
293 : if (report) then
294 : write(*,'(A)')
295 : write(*,4) 'enter Get_results k slvr_itr model gradr grada scale_height ' // trim(MLT_option), &
296 : k, s% solver_iter, s% model_number, gradr%val, grada%val, scale_height%val
297 : end if
298 :
299 79994 : if (k >= 1) then
300 79994 : s% dvc_dt_TDC(k) = 0d0
301 : end if
302 79994 : if (using_TDC) then
303 : if (report) write(*,3) 'call set_TDC', k, s% solver_iter
304 65210 : if (s% okay_to_set_mlt_vc) then
305 42492 : conv_vel_start = s% mlt_vc_old(k)
306 : else
307 22718 : conv_vel_start = s% mlt_vc(k)
308 : end if
309 :
310 : ! Set scale for judging the TDC luminosity equation Q(Y)=0.
311 : ! Q has units of a luminosity, so the scale should be a luminosity.
312 65210 : if (s% solver_iter == 0) then
313 39866797 : scale = max(abs(s% L(k)), 1d-3*maxval(s% L(1:s% nz)))
314 : else
315 34776621 : scale = max(abs(s% L_start(k)), 1d-3*maxval(s% L_start(1:s% nz)))
316 : end if
317 :
318 : call set_TDC(&
319 : conv_vel_start, mixing_length_alpha, s% alpha_TDC_DAMP, s%alpha_TDC_DAMPR, s%alpha_TDC_PtdVdt, &
320 : s%dt, cgrav, m, report, &
321 : mixing_type, scale, chiT, chiRho, gradr, r, Ptot, T, rho, dV, Cp, opacity, &
322 : scale_height, gradL, grada, conv_vel, D, Y_face, gradT, s%tdc_num_iters(k), max_conv_vel, &
323 65210 : Eq_div_w, grav, s% include_mlt_corr_to_TDC, s% alpha_TDC_C, s% alpha_TDC_S, ierr)
324 65210 : s% dvc_dt_TDC(k) = (conv_vel%val - conv_vel_start) / s%dt
325 :
326 65210 : if (ierr /= 0) then
327 0 : if (s% report_ierr) write(*,*) 'ierr from set_TDC'
328 0 : return
329 : end if
330 :
331 : ! Experimental method to lower superadiabaticity. Call TDC again with an artificially reduced
332 : ! gradr if the resulting gradT would lead to the radiative luminosity approaching the Eddington
333 : ! limit, or when a density inversion is expected to happen.
334 : ! This is meant as an implicit alternative to okay_to_reduce_gradT_excess
335 65210 : if (s% use_superad_reduction) then
336 0 : call set_superad_reduction
337 0 : if (Gamma_factor > 1d0) then
338 : call set_TDC(&
339 : conv_vel_start, mixing_length_alpha, s% alpha_TDC_DAMP, s%alpha_TDC_DAMPR, s%alpha_TDC_PtdVdt, &
340 : s%dt, cgrav, m, report, &
341 : mixing_type, scale, chiT, chiRho, gradr_scaled, r, Ptot, T, rho, dV, Cp, opacity, &
342 : scale_height, gradL, grada, conv_vel, D, Y_face, gradT, s%tdc_num_iters(k), max_conv_vel, &
343 0 : Eq_div_w, grav, s% include_mlt_corr_to_TDC, s% alpha_TDC_C, s% alpha_TDC_S, ierr)
344 0 : s% dvc_dt_TDC(k) = (conv_vel%val - conv_vel_start) / s%dt
345 0 : if (ierr /= 0) then
346 0 : if (s% report_ierr) write(*,*) 'ierr from set_TDC when using superad_reduction'
347 0 : return
348 : end if
349 : end if
350 : end if
351 :
352 14784 : else if (gradr > gradL) then
353 : if (report) write(*,3) 'call set_MLT', k, s% solver_iter
354 : call set_MLT(MLT_option, mixing_length_alpha, s% Henyey_MLT_nu_param, s% Henyey_MLT_y_param, &
355 : chiT, chiRho, Cp, grav, Lambda, rho, Ptot, T, opacity, &
356 : gradr, grada, gradL, &
357 2729 : Gamma, gradT, Y_face, conv_vel, D, mixing_type, max_conv_vel, ierr)
358 :
359 :
360 2729 : if (ierr /= 0) then
361 0 : if (s% report_ierr) write(*,*) 'ierr from set_MLT'
362 0 : return
363 : end if
364 :
365 : ! Experimental method to lower superadiabaticity. Call MLT again with an artificially reduced
366 : ! gradr if the resulting gradT would lead to the radiative luminosity approaching the Eddington
367 : ! limit, or when a density inversion is expected to happen.
368 : ! This is meant as an implicit alternative to okay_to_reduce_gradT_excess
369 2729 : if (s% use_superad_reduction) then
370 0 : call set_superad_reduction
371 0 : if (Gamma_factor > 1d0) then
372 : call set_MLT(MLT_option, mixing_length_alpha, s% Henyey_MLT_nu_param, s% Henyey_MLT_y_param, &
373 : chiT, chiRho, Cp, grav, Lambda, rho, Ptot, T, opacity, &
374 : gradr_scaled, grada, gradL, &
375 0 : Gamma, gradT, Y_face, conv_vel, D, mixing_type, max_conv_vel, ierr)
376 :
377 0 : if (ierr /= 0) then
378 0 : if (s% report_ierr) write(*,*) 'ierr from set_MLT when using superad_reduction'
379 0 : return
380 : end if
381 : end if
382 : end if
383 : end if
384 :
385 : ! If we're not convecting, try thermohaline and semiconvection.
386 79994 : if (mixing_type == no_mixing) then
387 68285 : if (gradL_composition_term < 0) then
388 : if (report) write(*,3) 'call set_thermohaline', k, s% solver_iter
389 : call set_thermohaline(s%thermohaline_option, Lambda, grada, gradr, T, opacity, rho, Cp, gradL_composition_term, &
390 : iso, XH1, thermohaline_coeff, &
391 0 : D, gradT, Y_face, conv_vel, mixing_type, ierr)
392 0 : if (ierr /= 0) then
393 0 : if (s% report_ierr) write(*,*) 'ierr from set_thermohaline'
394 0 : return
395 : end if
396 68285 : else if (gradr > grada) then
397 : if (report) write(*,3) 'call set_semiconvection', k, s% solver_iter
398 : call set_semiconvection(L, Lambda, m, T, Ptot, Pr, beta, opacity, rho, alpha_semiconvection, &
399 : s% semiconvection_option, cgrav, Cp, gradr, grada, gradL, &
400 : gradL_composition_term, &
401 0 : gradT, Y_face, conv_vel, D, mixing_type, ierr)
402 0 : if (ierr /= 0) then
403 0 : if (s% report_ierr) write(*,*) 'ierr from set_semiconvection'
404 0 : return
405 : end if
406 : end if
407 : end if
408 :
409 : ! If there's too-little mixing to bother, or we hit a bad value, fall back on no mixing.
410 79994 : if (D%val < s% remove_small_D_limit .or. is_bad(D%val)) then
411 : if (report) write(*,2) 'D < s% remove_small_D_limit', k, D%val, s% remove_small_D_limit
412 68285 : mixing_type = no_mixing
413 68285 : gradT = gradr
414 68285 : Y_face = gradT - gradL
415 68285 : conv_vel = 0d0
416 68285 : D = 0d0
417 68285 : Gamma = 0d0
418 : end if
419 :
420 : ! Prevent convection near center of model for MLT or TDC pulsations
421 : ! We don't check for the using_TDC flag, because mlt is sometimes called when using TDC
422 79994 : if ( s% TDC_num_innermost_cells_forced_nonturbulent > 0 .and. &
423 79994 : k > s% nz - s% TDC_num_innermost_cells_forced_nonturbulent) then
424 : if (report) write(*,2) 'make TDC center cells non-turbulent', k
425 0 : mixing_type = no_mixing
426 0 : gradT = gradr
427 0 : Y_face = gradT - gradL
428 0 : conv_vel = 0d0
429 0 : D = 0d0
430 0 : Gamma = 0d0
431 : end if
432 :
433 :
434 : contains
435 :
436 0 : subroutine set_superad_reduction()
437 0 : Gamma_limit = s% superad_reduction_Gamma_limit
438 0 : scale_value1 = s% superad_reduction_Gamma_limit_scale
439 0 : scale_value2 = s% superad_reduction_Gamma_inv_scale
440 0 : diff_grads_limit = s% superad_reduction_diff_grads_limit
441 0 : reduction_limit = s% superad_reduction_limit
442 0 : Lrad_div_Ledd = 4d0*crad/3d0*pow4(T)/P*gradT
443 0 : Gamma_inv_threshold = 4d0*(1d0-beta)/(4d0-3*beta)
444 :
445 0 : Gamma_factor = 1d0
446 0 : if (gradT > gradL) then
447 0 : if (Lrad_div_Ledd > Gamma_limit .or. Lrad_div_Ledd > Gamma_inv_threshold) then
448 0 : alfa0 = (gradT-gradL)/diff_grads_limit
449 0 : if (alfa0 < 1d0) then
450 0 : diff_grads_factor = -alfa0*alfa0*alfa0*(-10d0 + alfa0*(15d0 - 6d0*alfa0))
451 : else
452 0 : diff_grads_factor = 1d0
453 : end if
454 :
455 0 : Gamma_term = 0d0
456 : !if (Lrad_div_Ledd > Gamma_limit) then
457 : ! Gamma_term = Gamma_term + scale_value1*pow2(Lrad_div_Ledd/Gamma_limit-1d0)
458 : !end if
459 : !if (Lrad_div_Ledd% val > Gamma_inv_threshold) then
460 : ! Gamma_term = Gamma_term + scale_value2*pow2(Lrad_div_Ledd/Gamma_inv_threshold-1d0)
461 : !end if
462 0 : if (Lrad_div_Ledd > Gamma_limit) then
463 0 : alfa0 = Lrad_div_Ledd/Gamma_limit-1d0
464 0 : if (alfa0 < 1d0) then
465 0 : Gamma_term = Gamma_term + scale_value1*(0.5d0*alfa0*alfa0)
466 : else
467 0 : Gamma_term = Gamma_term + scale_value1*(alfa0-0.5d0)
468 : end if
469 : !Gamma_term = Gamma_term + scale_value1*pow2(Lrad_div_Ledd/Gamma_limit-1d0)
470 : end if
471 0 : if (Lrad_div_Ledd% val > Gamma_inv_threshold) then
472 0 : alfa0 = Lrad_div_Ledd/Gamma_inv_threshold-1d0
473 0 : if (alfa0 < 1d0) then
474 0 : Gamma_term = Gamma_term + scale_value1*(0.5d0*alfa0*alfa0)
475 : else
476 0 : Gamma_term = Gamma_term + scale_value1*(alfa0-0.5d0)
477 : end if
478 : !Gamma_term = Gamma_term + scale_value2*pow2(Lrad_div_Ledd/Gamma_inv_threshold-1d0)
479 : end if
480 :
481 0 : if (Gamma_term > 0d0) then
482 0 : Gamma_factor = Gamma_term/pow(beta,0.5d0)*diff_grads_factor
483 0 : Gamma_factor = Gamma_factor + 1d0
484 0 : if (reduction_limit > 1d0) then
485 0 : lambda_limit = 2d0/(reduction_limit-1d0)
486 0 : exp_limit = exp(-lambda_limit*(Gamma_factor-1d0))
487 0 : Gamma_factor = 2d0*(reduction_limit-1d0)*(1d0/(1d0+exp_limit)-0.5d0)+1d0
488 : end if
489 : end if
490 : end if
491 : end if
492 0 : if (k /= 0) s% superad_reduction_factor(k) = Gamma_factor% val
493 0 : if (Gamma_factor > 1d0) then
494 0 : grad_scale = (gradr-gradL)/(Gamma_factor*gradr) + gradL/gradr
495 0 : gradr_scaled = grad_scale*gradr
496 : end if
497 79994 : end subroutine set_superad_reduction
498 : end subroutine Get_results
499 :
500 : end module turb_support
|