LCOV - code coverage report
Current view: top level - star/private - turb_support.f90 (source / functions) Coverage Total Hit
Test: coverage.info Lines: 44.4 % 205 91
Test Date: 2025-11-11 15:19:17 Functions: 50.0 % 4 2

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

Generated by: LCOV version 2.0-1