Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010-2019 Meridith Joyce & 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 for handling starspot calculations
21 : ! contains parameterized YREC routines
22 : ! note: MESA models a pressure contrast rather than temperature contrast
23 :
24 : module starspots
25 :
26 : use star_private_def
27 : use const_def, only: dp, pi4, boltz_sigma, cgas
28 : use utils_lib
29 :
30 : implicit none
31 :
32 : private
33 : public :: starspot_tweak_gradr
34 : public :: starspot_tweak_PT
35 : public :: starspot_restore_PT
36 :
37 : real(dp) :: L_init
38 :
39 : contains
40 :
41 0 : subroutine starspot_tweak_gradr(s, P, gradr, gradr_spot)
42 : ! Purpose: adjusts the gradient of the radius to account for starspots.
43 : !
44 : ! This subroutine is called at the beginning of Get_results()
45 : ! in turb_support.f90
46 :
47 : use auto_diff_support
48 : type(star_info), pointer :: s
49 : type(auto_diff_real_star_order1), intent(in) :: P
50 : type(auto_diff_real_star_order1), intent(in) :: gradr
51 : type(auto_diff_real_star_order1), intent(out) :: gradr_spot
52 :
53 0 : real(dp) :: mu_ideal_gas, R2, Teff_local, PB_i
54 : type(auto_diff_real_star_order1) :: xspot_of_r ! xspot4
55 :
56 0 : if (.not. s%doing_relax .and. .not. s%doing_first_model_of_run) then
57 0 : mu_ideal_gas = s%mu(1) !1.00794d0 ! for hydrogen, 1 gram per mole
58 0 : R2 = pow2(s%R(1))
59 0 : Teff_local = pow(s%L(1)/(pi4*boltz_sigma*R2), 0.25d0)
60 0 : PB_i = (cgas*s%rho(1)/mu_ideal_gas)*(1.0d0 - s%xspot)*Teff_local
61 0 : xspot_of_r = (P - PB_i)/P
62 0 : gradr_spot = gradr/(s%fspot*pow4(xspot_of_r) + 1.0d0 - s%fspot)
63 : else
64 0 : gradr_spot = gradr
65 : end if
66 0 : end subroutine starspot_tweak_gradr
67 :
68 0 : subroutine starspot_tweak_PT(s)
69 : ! Purpose: saves the surface luminosity and adjusts it and the effective
70 : ! temperature to account for starspots.
71 : !
72 : ! This subroutine is called at the beginning of get_surf_PT()
73 : ! in hydro_vars.f90
74 :
75 : type(star_info), pointer :: s
76 :
77 0 : real(dp) :: alp
78 :
79 0 : alp = 1d0 - s%fspot + s%fspot*pow4(s%xspot)
80 :
81 : ! This is the surface-average value for luminosity
82 0 : L_init = s%L(1)
83 :
84 : ! Set the surface L to the unspotted, ambient L
85 0 : s%L(1) = s%L(1)/alp
86 :
87 : ! Now, set the Teff. Used in atm table lookup to set boundary conditions
88 0 : s%Teff = pow(s%L(1)/(pi4*pow2(s%r(1))*boltz_sigma), 0.25_dp)
89 :
90 0 : end subroutine starspot_tweak_PT
91 :
92 0 : subroutine starspot_restore_PT(s)
93 : ! Purpose: restores the surface luminosity and effective temperature.
94 : !
95 : ! Called at the end of get_surf_PT()
96 :
97 : type(star_info), pointer :: s
98 :
99 0 : s%Teff = pow(L_init/(pi4*pow2(s%r(1))*boltz_sigma), 0.25_dp)
100 0 : s%L(1) = L_init
101 :
102 0 : end subroutine starspot_restore_PT
103 :
104 : end module starspots
|