Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010-2019 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 overshoot_step
21 :
22 : use star_private_def
23 : use overshoot_utils
24 :
25 : implicit none
26 :
27 : private
28 : public :: eval_overshoot_step
29 :
30 : contains
31 :
32 0 : subroutine eval_overshoot_step (s, i, j, k_a, k_b, D, vc, ierr)
33 :
34 : type(star_info), pointer :: s
35 : integer, intent(in) :: i
36 : integer, intent(in) :: j
37 : integer, intent(out) :: k_a
38 : integer, intent(out) :: k_b
39 : real(dp), intent(out) :: D(:)
40 : real(dp), intent(out) :: vc(:)
41 : integer, intent(out) :: ierr
42 :
43 0 : real(dp) :: f
44 : real(dp) :: f0
45 0 : real(dp) :: D0
46 0 : real(dp) :: Delta0
47 0 : real(dp) :: w
48 0 : real(dp) :: factor
49 0 : real(dp) :: Hp_cb
50 : integer :: k_ob
51 0 : real(dp) :: r_ob
52 0 : real(dp) :: D_ob
53 0 : real(dp) :: vc_ob
54 : logical :: outward
55 : integer :: dk
56 : integer :: k
57 0 : real(dp) :: r
58 0 : real(dp) :: dr
59 :
60 : ! Evaluate the overshoot diffusion coefficient D(k_a:k_b) and
61 : ! mixing velocity vc(k_a:k_b) at the i'th convective boundary,
62 : ! using the j'th set of overshoot parameters. The overshoot
63 : ! follows a simple step scheme
64 :
65 0 : ierr = 0
66 :
67 : ! Extract parameters
68 :
69 0 : f = s%overshoot_f(j)
70 0 : f0 = s%overshoot_f0(j)
71 :
72 0 : D0 = s%overshoot_D0(j)
73 0 : Delta0 = s%overshoot_Delta0(j)
74 :
75 0 : if (f <= 0._dp .OR. f0 <= 0._dp) then
76 0 : write(*,*) 'ERROR: for step overshooting, must set f and f0 > 0'
77 0 : write(*,*) 'see description of overshooting in star/defaults/control.defaults'
78 0 : ierr = -1
79 0 : return
80 : end if
81 :
82 : ! Apply mass limits
83 :
84 0 : if (s%star_mass < s%overshoot_mass_full_on(j)) then
85 0 : if (s%star_mass > s%overshoot_mass_full_off(j)) then
86 : w = (s%star_mass - s%overshoot_mass_full_off(j)) / &
87 0 : (s%overshoot_mass_full_on(j) - s%overshoot_mass_full_off(j))
88 0 : factor = 0.5_dp*(1._dp - cospi(w))
89 0 : f = f*factor
90 0 : f0 = f0*factor
91 : else
92 0 : f = 0._dp
93 0 : f0 = 0._dp
94 : end if
95 : end if
96 :
97 : ! Evaluate convective boundary (_cb) parameters
98 :
99 0 : call eval_conv_bdy_Hp(s, i, Hp_cb, ierr)
100 0 : if (ierr /= 0) return
101 :
102 : ! Evaluate overshoot boundary (_ob) parameters
103 :
104 0 : call eval_over_bdy_params(s, i, f0, k_ob, r_ob, D_ob, vc_ob, ierr)
105 0 : if (ierr /= 0) return
106 :
107 : ! Loop over cell faces, adding overshoot until D <= overshoot_D_min
108 :
109 0 : outward = s%top_conv_bdy(i)
110 :
111 0 : if (outward) then
112 0 : k_a = k_ob
113 0 : k_b = 1
114 0 : dk = -1
115 : else
116 0 : k_a = k_ob+1
117 0 : k_b = s%nz
118 0 : dk = 1
119 : end if
120 :
121 0 : face_loop : do k = k_a, k_b, dk
122 :
123 : ! Evaluate the step factor
124 :
125 0 : r = s%r(k)
126 :
127 0 : if (outward) then
128 0 : dr = r - r_ob
129 : else
130 0 : dr = r_ob - r
131 : end if
132 :
133 0 : if (dr < f*Hp_cb) then
134 : factor = 1._dp
135 : else
136 0 : factor = 0._dp
137 : end if
138 :
139 : ! Store the diffusion coefficient and velocity
140 :
141 0 : D(k) = (D0 + Delta0*D_ob)*factor
142 0 : if(D_ob /= 0d0) then
143 0 : vc(k) = (D0/D_ob + Delta0)*vc_ob*factor
144 : else
145 0 : vc(k) = 0d0
146 : end if
147 : ! Check for early overshoot completion
148 :
149 0 : if (D(k) < s%overshoot_D_min) then
150 0 : k_b = k
151 0 : exit face_loop
152 : end if
153 :
154 : end do face_loop
155 :
156 0 : ierr = 0
157 :
158 0 : return
159 :
160 : end subroutine eval_overshoot_step
161 :
162 : end module overshoot_step
|