Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2020 Adam Jermyn & 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 auto_diff_support
21 :
22 : use star_private_def
23 : use const_def, only: dp, sqrt_2_div_3
24 : use auto_diff
25 :
26 : implicit none
27 :
28 : ! current use of xtra's
29 : ! xtra1 is w_div_wc
30 :
31 : public
32 :
33 : contains
34 :
35 52348 : type(auto_diff_real_star_order1) function shift_p1(val_00) result(val_p1)
36 : type(auto_diff_real_star_order1), intent(in) :: val_00
37 : integer :: j
38 52348 : val_p1%val = val_00%val
39 52348 : do j=auto_diff_star_num_vars-2,1,-3 ! p1 gets 00, 00 gets m1, m1 gets 0d0
40 575828 : val_p1%d1Array(j+2) = val_00%d1Array(j+1)
41 575828 : val_p1%d1Array(j+1) = val_00%d1Array(j)
42 575828 : val_p1%d1Array(j) = 0d0
43 : end do
44 52348 : end function shift_p1
45 :
46 52304 : type(auto_diff_real_star_order1) function shift_m1(val_00) result(val_m1)
47 : type(auto_diff_real_star_order1), intent(in) :: val_00
48 : integer :: j
49 52304 : val_m1%val = val_00%val
50 52304 : do j=1,auto_diff_star_num_vars,3 ! m1 gets 00, 00 gets p1, p1 gets 0d0
51 575344 : val_m1%d1Array(j) = val_00%d1Array(j+1)
52 575344 : val_m1%d1Array(j+1) = val_00%d1Array(j+2)
53 575344 : val_m1%d1Array(j+2) = 0d0
54 : end do
55 52304 : end function shift_m1
56 :
57 209392 : subroutine unwrap(var, val, &
58 : dlnd_m1, dlnd_00, dlnd_p1, &
59 : dlnT_m1, dlnT_00, dlnT_p1, &
60 : dw_m1, dw_00, dw_p1, &
61 : dlnR_m1, dlnR_00, dlnR_p1, &
62 : dv_m1, dv_00, dv_p1, &
63 : dL_m1, dL_00, dL_p1, &
64 : dHp_m1, dHp_00, dHp_p1, &
65 : dw_div_wc_m1, dw_div_wc_00, dw_div_wc_p1, &
66 : djrot_m1, djrot_00, djrot_p1, &
67 : dxtra1_m1, dxtra1_00, dxtra1_p1, &
68 : dxtra2_m1, dxtra2_00, dxtra2_p1)
69 : type(auto_diff_real_star_order1), intent(in) :: var
70 : real(dp), intent(out) :: &
71 : val, dlnd_m1, dlnd_00, dlnd_p1, dlnT_m1, dlnT_00, dlnT_p1, &
72 : dw_m1, dw_00, dw_p1, dlnR_m1, dlnR_00, dlnR_p1, &
73 : dv_m1, dv_00, dv_p1, dL_m1, dL_00, dL_p1, &
74 : dHp_m1, dHp_00, dHp_p1, &
75 : dw_div_wc_m1, dw_div_wc_00, dw_div_wc_p1, &
76 : djrot_m1, djrot_00, djrot_p1, &
77 : dxtra1_m1, dxtra1_00, dxtra1_p1, &
78 : dxtra2_m1, dxtra2_00, dxtra2_p1
79 209392 : val = var%val
80 209392 : dlnd_m1 = var%d1Array(i_lnd_m1)
81 209392 : dlnd_00 = var%d1Array(i_lnd_00)
82 209392 : dlnd_p1 = var%d1Array(i_lnd_p1)
83 209392 : dlnT_m1 = var%d1Array(i_lnT_m1)
84 209392 : dlnT_00 = var%d1Array(i_lnT_00)
85 209392 : dlnT_p1 = var%d1Array(i_lnT_p1)
86 209392 : dw_m1 = var%d1Array(i_w_m1)
87 209392 : dw_00 = var%d1Array(i_w_00)
88 209392 : dw_p1 = var%d1Array(i_w_p1)
89 209392 : dlnR_m1 = var%d1Array(i_lnR_m1)
90 209392 : dlnR_00 = var%d1Array(i_lnR_00)
91 209392 : dlnR_p1 = var%d1Array(i_lnR_p1)
92 209392 : dv_m1 = var%d1Array(i_v_m1)
93 209392 : dv_00 = var%d1Array(i_v_00)
94 209392 : dv_p1 = var%d1Array(i_v_p1)
95 209392 : dL_m1 = var%d1Array(i_L_m1)
96 209392 : dL_00 = var%d1Array(i_L_00)
97 209392 : dL_p1 = var%d1Array(i_L_p1)
98 209392 : dHp_m1 = var%d1Array(i_Hp_m1)
99 209392 : dHp_00 = var%d1Array(i_Hp_00)
100 209392 : dHp_p1 = var%d1Array(i_Hp_p1)
101 209392 : dw_div_wc_m1 = var%d1Array(i_w_div_wc_m1)
102 209392 : dw_div_wc_00 = var%d1Array(i_w_div_wc_00)
103 209392 : dw_div_wc_p1 = var%d1Array(i_w_div_wc_p1)
104 209392 : djrot_m1 = var%d1Array(i_jrot_m1)
105 209392 : djrot_00 = var%d1Array(i_jrot_00)
106 209392 : djrot_p1 = var%d1Array(i_jrot_p1)
107 209392 : dxtra1_m1 = var%d1Array(i_xtra1_m1)
108 209392 : dxtra1_00 = var%d1Array(i_xtra1_00)
109 209392 : dxtra1_p1 = var%d1Array(i_xtra1_p1)
110 209392 : dxtra2_m1 = var%d1Array(i_xtra2_m1)
111 209392 : dxtra2_00 = var%d1Array(i_xtra2_00)
112 209392 : dxtra2_p1 = var%d1Array(i_xtra2_p1)
113 209392 : end subroutine unwrap
114 :
115 176 : subroutine wrap(var, val, &
116 : dlnd_m1, dlnd_00, dlnd_p1, &
117 : dlnT_m1, dlnT_00, dlnT_p1, &
118 : dw_m1, dw_00, dw_p1, &
119 : dlnR_m1, dlnR_00, dlnR_p1, &
120 : dv_m1, dv_00, dv_p1, &
121 : dL_m1, dL_00, dL_p1, &
122 : dHp_m1, dHp_00, dHp_p1, &
123 : dw_div_wc_m1, dw_div_wc_00, dw_div_wc_p1, &
124 : djrot_m1, djrot_00, djrot_p1, &
125 : dxtra1_m1, dxtra1_00, dxtra1_p1, &
126 : dxtra2_m1, dxtra2_00, dxtra2_p1)
127 : type(auto_diff_real_star_order1), intent(out) :: var
128 : real(dp), intent(in) :: &
129 : val, dlnd_m1, dlnd_00, dlnd_p1, dlnT_m1, dlnT_00, dlnT_p1, &
130 : dw_m1, dw_00, dw_p1, dlnR_m1, dlnR_00, dlnR_p1, &
131 : dv_m1, dv_00, dv_p1, dL_m1, dL_00, dL_p1, &
132 : dHp_m1, dHp_00, dHp_p1, &
133 : dw_div_wc_m1, dw_div_wc_00, dw_div_wc_p1, &
134 : djrot_m1, djrot_00, djrot_p1, &
135 : dxtra1_m1, dxtra1_00, dxtra1_p1, &
136 : dxtra2_m1, dxtra2_00, dxtra2_p1
137 176 : var%val = val
138 176 : var%d1Array(i_lnd_m1) = dlnd_m1
139 176 : var%d1Array(i_lnd_00) = dlnd_00
140 176 : var%d1Array(i_lnd_p1) = dlnd_p1
141 176 : var%d1Array(i_lnT_m1) = dlnT_m1
142 176 : var%d1Array(i_lnT_00) = dlnT_00
143 176 : var%d1Array(i_lnT_p1) = dlnT_p1
144 176 : var%d1Array(i_w_m1) = dw_m1
145 176 : var%d1Array(i_w_00) = dw_00
146 176 : var%d1Array(i_w_p1) = dw_p1
147 176 : var%d1Array(i_lnR_m1) = dlnR_m1
148 176 : var%d1Array(i_lnR_00) = dlnR_00
149 176 : var%d1Array(i_lnR_p1) = dlnR_p1
150 176 : var%d1Array(i_v_m1) = dv_m1
151 176 : var%d1Array(i_v_00) = dv_00
152 176 : var%d1Array(i_v_p1) = dv_p1
153 176 : var%d1Array(i_L_m1) = dL_m1
154 176 : var%d1Array(i_L_00) = dL_00
155 176 : var%d1Array(i_L_p1) = dL_p1
156 176 : var%d1Array(i_Hp_m1) = dHp_m1
157 176 : var%d1Array(i_Hp_00) = dHp_00
158 176 : var%d1Array(i_Hp_p1) = dHp_p1
159 176 : var%d1Array(i_w_div_wc_m1) = dw_div_wc_m1
160 176 : var%d1Array(i_w_div_wc_00) = dw_div_wc_00
161 176 : var%d1Array(i_w_div_wc_p1) = dw_div_wc_p1
162 176 : var%d1Array(i_jrot_m1) = djrot_m1
163 176 : var%d1Array(i_jrot_00) = djrot_00
164 176 : var%d1Array(i_jrot_p1) = djrot_p1
165 176 : var%d1Array(i_xtra1_m1) = dxtra1_m1
166 176 : var%d1Array(i_xtra1_00) = dxtra1_00
167 176 : var%d1Array(i_xtra1_p1) = dxtra1_p1
168 176 : var%d1Array(i_xtra2_m1) = dxtra2_m1
169 176 : var%d1Array(i_xtra2_00) = dxtra2_00
170 176 : var%d1Array(i_xtra2_p1) = dxtra2_p1
171 176 : end subroutine wrap
172 :
173 :
174 : ! The following routines turn regular star variables into auto_diff_real_star_order1 variables.
175 : ! For independent variables this is a straightforward wrapping. For dependent variables like eos and kap
176 : ! outputs we also pull in information about their partials from the relevant module.
177 :
178 : !----------------------------------------------------------------------------------------------------
179 : !
180 : ! We begin with quantities defined on cells (T, rho, e_turb, kap, pressure, energy, entropy).
181 : ! We need these on cells k-1, k and k+1.
182 : !
183 : !----------------------------------------------------------------------------------------------------
184 :
185 :
186 212160 : function wrap_T_m1(s, k) result(T_m1)
187 : type (star_info), pointer :: s
188 : type(auto_diff_real_star_order1) :: T_m1
189 : integer, intent(in) :: k
190 212160 : T_m1 = 0d0
191 212160 : if (k > 1) then
192 212160 : T_m1 % val = s%T(k-1)
193 212160 : T_m1 % d1Array(i_lnT_m1) = s%T(k-1)
194 : end if
195 212160 : end function wrap_T_m1
196 :
197 212292 : function wrap_T_00(s, k) result(T_00)
198 : type (star_info), pointer :: s
199 : type(auto_diff_real_star_order1) :: T_00
200 : integer, intent(in) :: k
201 212292 : T_00 = 0d0
202 212292 : T_00 % val = s%T(k)
203 212292 : T_00 % d1Array(i_lnT_00) = s%T(k)
204 212292 : end function wrap_T_00
205 :
206 0 : function wrap_T_p1(s, k) result(T_p1)
207 : type (star_info), pointer :: s
208 : type(auto_diff_real_star_order1) :: T_p1
209 : integer, intent(in) :: k
210 0 : T_p1 = 0d0
211 0 : if (k < s%nz) then
212 0 : T_p1 % val = s%T(k+1)
213 0 : T_p1 % d1Array(i_lnT_p1) = s%T(k+1)
214 : end if
215 0 : end function wrap_T_p1
216 :
217 0 : function wrap_lnT_m1(s, k) result(lnT_m1)
218 : type (star_info), pointer :: s
219 : type(auto_diff_real_star_order1) :: lnT_m1
220 : integer, intent(in) :: k
221 0 : lnT_m1 = 0d0
222 0 : if (k > 1) then
223 0 : lnT_m1 % val = s%lnT(k-1)
224 0 : lnT_m1 % d1Array(i_lnT_m1) = 1d0
225 : end if
226 0 : end function wrap_lnT_m1
227 :
228 44 : function wrap_lnT_00(s, k) result(lnT_00)
229 : type (star_info), pointer :: s
230 : type(auto_diff_real_star_order1) :: lnT_00
231 : integer, intent(in) :: k
232 44 : lnT_00 = 0d0
233 44 : lnT_00 % val = s%lnT(k)
234 44 : lnT_00 % d1Array(i_lnT_00) = 1d0
235 44 : end function wrap_lnT_00
236 :
237 0 : function wrap_lnT_p1(s, k) result(lnT_p1)
238 : type (star_info), pointer :: s
239 : type(auto_diff_real_star_order1) :: lnT_p1
240 : integer, intent(in) :: k
241 0 : lnT_p1 = 0d0
242 0 : if (k < s%nz) then
243 0 : lnT_p1 % val = s%lnT(k+1)
244 0 : lnT_p1 % d1Array(i_lnT_p1) = 1d0
245 : end if
246 0 : end function wrap_lnT_p1
247 :
248 0 : function wrap_dxh_lnT(s, k) result(dxh_lnT)
249 : type (star_info), pointer :: s
250 : type(auto_diff_real_star_order1) :: dxh_lnT
251 : integer, intent(in) :: k
252 0 : dxh_lnT = 0d0
253 0 : dxh_lnT % val = s%dxh_lnT(k)
254 0 : dxh_lnT % d1Array(i_lnT_00) = 1d0
255 0 : end function wrap_dxh_lnT
256 :
257 239784 : function wrap_d_m1(s, k) result(d_m1)
258 : type (star_info), pointer :: s
259 : type(auto_diff_real_star_order1) :: d_m1
260 : integer, intent(in) :: k
261 239784 : d_m1 = 0d0
262 239784 : if (k > 1) then
263 239784 : d_m1 % val = s%rho(k-1)
264 239784 : d_m1 % d1Array(i_lnd_m1) = s%rho(k-1)
265 : end if
266 239784 : end function wrap_d_m1
267 :
268 292330 : function wrap_d_00(s, k) result(d_00)
269 : type (star_info), pointer :: s
270 : type(auto_diff_real_star_order1) :: d_00
271 : integer, intent(in) :: k
272 292330 : d_00 = 0d0
273 292330 : d_00 % val = s%rho(k)
274 292330 : d_00 % d1Array(i_lnd_00) = s%rho(k)
275 292330 : end function wrap_d_00
276 :
277 0 : function wrap_d_p1(s, k) result(d_p1)
278 : type (star_info), pointer :: s
279 : type(auto_diff_real_star_order1) :: d_p1
280 : integer, intent(in) :: k
281 0 : d_p1 = 0d0
282 0 : if (k < s%nz) then
283 0 : d_p1 % val = s%rho(k+1)
284 0 : d_p1 % d1Array(i_lnd_p1) = s%rho(k+1)
285 : end if
286 0 : end function wrap_d_p1
287 :
288 79928 : function wrap_d_m1_start(s, k) result(d_m1)
289 : type (star_info), pointer :: s
290 : type(auto_diff_real_star_order1) :: d_m1
291 : integer, intent(in) :: k
292 79928 : d_m1 = 0d0
293 79928 : if (k > 1) then
294 79928 : d_m1 % val = s%rho_start(k-1)
295 79928 : d_m1 % d1Array(i_lnd_m1) = s%rho_start(k-1)
296 : end if
297 79928 : end function wrap_d_m1_start
298 :
299 79994 : function wrap_d_00_start(s, k) result(d_00)
300 : type (star_info), pointer :: s
301 : type(auto_diff_real_star_order1) :: d_00
302 : integer, intent(in) :: k
303 79994 : d_00 = 0d0
304 79994 : d_00 % val = s%rho_start(k)
305 79994 : d_00 % d1Array(i_lnd_00) = s%rho_start(k)
306 79994 : end function wrap_d_00_start
307 :
308 0 : function wrap_d_p1_start(s, k) result(d_p1)
309 : type (star_info), pointer :: s
310 : type(auto_diff_real_star_order1) :: d_p1
311 : integer, intent(in) :: k
312 0 : d_p1 = 0d0
313 0 : if (k < s%nz) then
314 0 : d_p1 % val = s%rho_start(k+1)
315 0 : d_p1 % d1Array(i_lnd_p1) = s%rho_start(k+1)
316 : end if
317 0 : end function wrap_d_p1_start
318 :
319 0 : function wrap_lnd_m1(s, k) result(lnd_m1)
320 : type (star_info), pointer :: s
321 : type(auto_diff_real_star_order1) :: lnd_m1
322 : integer, intent(in) :: k
323 0 : lnd_m1 = 0d0
324 0 : if (k > 1) then
325 0 : lnd_m1 % val = s%lnd(k-1)
326 0 : lnd_m1 % d1Array(i_lnd_m1) = 1d0
327 : end if
328 0 : end function wrap_lnd_m1
329 :
330 0 : function wrap_lnd_00(s, k) result(lnd_00)
331 : type (star_info), pointer :: s
332 : type(auto_diff_real_star_order1) :: lnd_00
333 : integer, intent(in) :: k
334 0 : lnd_00 = 0d0
335 0 : lnd_00 % val = s%lnd(k)
336 0 : lnd_00 % d1Array(i_lnd_00) = 1d0
337 0 : end function wrap_lnd_00
338 :
339 0 : function wrap_lnd_p1(s, k) result(lnd_p1)
340 : type (star_info), pointer :: s
341 : type(auto_diff_real_star_order1) :: lnd_p1
342 : integer, intent(in) :: k
343 0 : lnd_p1 = 0d0
344 0 : if (k < s%nz) then
345 0 : lnd_p1 % val = s%lnd(k+1)
346 0 : lnd_p1 % d1Array(i_lnd_p1) = 1d0
347 : end if
348 0 : end function wrap_lnd_p1
349 :
350 0 : function wrap_dxh_lnd(s, k) result(dxh_lnd)
351 : type (star_info), pointer :: s
352 : type(auto_diff_real_star_order1) :: dxh_lnd
353 : integer, intent(in) :: k
354 0 : dxh_lnd = 0d0
355 0 : dxh_lnd % val = s%dxh_lnd(k)
356 0 : dxh_lnd % d1Array(i_lnd_00) = 1d0
357 0 : end function wrap_dxh_lnd
358 :
359 0 : function wrap_w_m1(s, k) result(w_m1)
360 : type (star_info), pointer :: s
361 : type(auto_diff_real_star_order1) :: w_m1
362 : integer, intent(in) :: k
363 0 : w_m1 = 0d0
364 0 : if (k > 1) then
365 0 : w_m1 % val = s%w(k-1)
366 0 : w_m1 % d1Array(i_w_m1) = 1d0
367 : end if
368 0 : end function wrap_w_m1
369 :
370 0 : function wrap_w_00(s, k) result(w_00)
371 : type (star_info), pointer :: s
372 : type(auto_diff_real_star_order1) :: w_00
373 : integer, intent(in) :: k
374 0 : w_00 = 0d0
375 0 : w_00 % val = s%w(k)
376 0 : w_00 % d1Array(i_w_00) = 1d0
377 0 : end function wrap_w_00
378 :
379 0 : function wrap_w_p1(s, k) result(w_p1)
380 : type (star_info), pointer :: s
381 : type(auto_diff_real_star_order1) :: w_p1
382 : integer, intent(in) :: k
383 0 : w_p1 = 0d0
384 0 : if (k < s%nz) then
385 0 : w_p1 % val = s%w(k+1)
386 0 : w_p1 % d1Array(i_w_p1) = 1d0
387 : end if
388 0 : end function wrap_w_p1
389 :
390 0 : real(dp) function get_etrb(s,k)
391 : type (star_info), pointer :: s
392 : integer, intent(in) :: k
393 0 : get_etrb = pow2(s% w(k))
394 0 : end function get_etrb
395 :
396 0 : real(dp) function get_etrb_start(s,k)
397 : type (star_info), pointer :: s
398 : integer, intent(in) :: k
399 0 : get_etrb_start = pow2(s% w_start(k))
400 0 : end function get_etrb_start
401 :
402 0 : real(dp) function get_RSP2_conv_velocity(s,k) result (cv) ! at face k
403 : type (star_info), pointer :: s
404 : integer, intent(in) :: k
405 0 : real(dp) :: alfa, beta
406 0 : if (k == 1) then
407 : cv = 0d0
408 : else
409 0 : alfa = s% dq(k-1)/(s% dq(k-1) + s% dq(k))
410 0 : beta = 1d0 - alfa
411 0 : cv = sqrt_2_div_3*(alfa*s% w(k) + beta*s% w(k-1))
412 : end if
413 0 : end function get_RSP2_conv_velocity
414 :
415 0 : real(dp) function get_w(s,k)
416 : type (star_info), pointer :: s
417 : integer, intent(in) :: k
418 0 : get_w = s% w(k)
419 0 : end function get_w
420 :
421 0 : real(dp) function get_w_start(s,k)
422 : type (star_info), pointer :: s
423 : integer, intent(in) :: k
424 0 : get_w_start = s% w_start(k)
425 0 : end function get_w_start
426 :
427 0 : function wrap_etrb_m1(s, k) result(etrb_m1)
428 : type (star_info), pointer :: s
429 : type(auto_diff_real_star_order1) :: etrb_m1
430 : integer, intent(in) :: k
431 0 : etrb_m1 = pow2(wrap_w_m1(s,k))
432 0 : end function wrap_etrb_m1
433 :
434 0 : function wrap_etrb_00(s, k) result(etrb_00)
435 : type (star_info), pointer :: s
436 : type(auto_diff_real_star_order1) :: etrb_00
437 : integer, intent(in) :: k
438 0 : etrb_00 = pow2(wrap_w_00(s,k))
439 0 : end function wrap_etrb_00
440 :
441 0 : function wrap_etrb_p1(s, k) result(etrb_p1)
442 : type (star_info), pointer :: s
443 : type(auto_diff_real_star_order1) :: etrb_p1
444 : integer, intent(in) :: k
445 0 : etrb_p1 = pow2(wrap_w_p1(s,k))
446 0 : end function wrap_etrb_p1
447 :
448 159856 : function wrap_kap_m1(s, k) result(kap_m1)
449 : type (star_info), pointer :: s
450 : type(auto_diff_real_star_order1) :: kap_m1
451 : integer, intent(in) :: k
452 159856 : kap_m1 = 0d0
453 159856 : if (k > 1) then
454 159856 : kap_m1 % val = s%opacity(k-1)
455 159856 : kap_m1 % d1Array(i_lnd_m1) = s%d_opacity_dlnd(k-1)
456 159856 : kap_m1 % d1Array(i_lnT_m1) = s%d_opacity_dlnT(k-1)
457 : end if
458 159856 : end function wrap_kap_m1
459 :
460 159988 : function wrap_kap_00(s, k) result(kap_00)
461 : type (star_info), pointer :: s
462 : type(auto_diff_real_star_order1) :: kap_00
463 : integer, intent(in) :: k
464 159988 : kap_00 = 0d0
465 159988 : kap_00 % val = s%opacity(k)
466 159988 : kap_00 % d1Array(i_lnd_00) = s%d_opacity_dlnd(k)
467 159988 : kap_00 % d1Array(i_lnT_00) = s%d_opacity_dlnT(k)
468 159988 : end function wrap_kap_00
469 :
470 0 : function wrap_kap_p1(s, k) result(kap_p1)
471 : type (star_info), pointer :: s
472 : type(auto_diff_real_star_order1) :: kap_p1
473 : integer, intent(in) :: k
474 0 : kap_p1 = 0d0
475 0 : if (k < s%nz) then
476 0 : kap_p1 % val = s%opacity(k+1)
477 0 : kap_p1 % d1Array(i_lnd_p1) = s%d_opacity_dlnd(k+1)/s% rho(k+1)
478 0 : kap_p1 % d1Array(i_lnT_p1) = s%d_opacity_dlnT(k+1)
479 : end if
480 0 : end function wrap_kap_p1
481 :
482 0 : function wrap_s_m1(s, k) result(s_m1)
483 : type (star_info), pointer :: s
484 : type(auto_diff_real_star_order1) :: s_m1
485 : integer, intent(in) :: k
486 0 : s_m1 = 0d0
487 0 : if (k > 1) then
488 0 : s_m1%val = s% entropy(k-1)
489 0 : s_m1%d1Array(i_lnd_m1) = s% dS_dRho_for_partials(k-1)*s% rho(k-1)
490 0 : s_m1%d1Array(i_lnT_m1) = s% dS_dT_for_partials(k-1)*s% T(k-1)
491 : end if
492 0 : end function wrap_s_m1
493 :
494 0 : function wrap_s_00(s, k) result(s_00)
495 : type (star_info), pointer :: s
496 : type(auto_diff_real_star_order1) :: s_00
497 : integer, intent(in) :: k
498 0 : s_00 = 0d0
499 0 : s_00%val = s% entropy(k)
500 0 : s_00%d1Array(i_lnd_00) = s% dS_dRho_for_partials(k)*s% rho(k)
501 0 : s_00%d1Array(i_lnT_00) = s% dS_dT_for_partials(k)*s% T(k)
502 0 : end function wrap_s_00
503 :
504 0 : function wrap_s_p1(s, k) result(s_p1)
505 : type (star_info), pointer :: s
506 : type(auto_diff_real_star_order1) :: s_p1
507 : integer, intent(in) :: k
508 0 : s_p1 = 0d0
509 0 : if (k < s%nz) then
510 0 : s_p1%val = s% entropy(k+1)
511 0 : s_p1%d1Array(i_lnd_p1) = s% dS_dRho_for_partials(k+1)*s% rho(k+1)
512 0 : s_p1%d1Array(i_lnT_p1) = s% dS_dT_for_partials(k+1)*s% T(k+1)
513 : end if
514 0 : end function wrap_s_p1
515 :
516 0 : function wrap_e_m1(s, k) result(e_m1)
517 : type (star_info), pointer :: s
518 : type(auto_diff_real_star_order1) :: e_m1
519 : integer, intent(in) :: k
520 0 : e_m1 = 0d0
521 0 : if (k > 1) then
522 0 : e_m1%val = s% energy(k-1)
523 0 : e_m1%d1Array(i_lnd_m1) = s% dE_dRho_for_partials(k-1)*s% rho(k-1)
524 0 : e_m1%d1Array(i_lnT_m1) = s% Cv_for_partials(k-1)*s% T(k-1)
525 : end if
526 0 : end function wrap_e_m1
527 :
528 0 : function wrap_e_00(s, k) result(e_00)
529 : type (star_info), pointer :: s
530 : type(auto_diff_real_star_order1) :: e_00
531 : integer, intent(in) :: k
532 0 : e_00 = 0d0
533 0 : e_00%val = s% energy(k)
534 0 : e_00%d1Array(i_lnd_00) = s% dE_dRho_for_partials(k)*s% rho(k)
535 0 : e_00%d1Array(i_lnT_00) = s% Cv_for_partials(k)*s% T(k)
536 0 : end function wrap_e_00
537 :
538 0 : function wrap_e_p1(s, k) result(e_p1)
539 : type (star_info), pointer :: s
540 : type(auto_diff_real_star_order1) :: e_p1
541 : integer, intent(in) :: k
542 0 : e_p1 = 0d0
543 0 : if (k < s%nz) then
544 0 : e_p1%val = s% energy(k+1)
545 0 : e_p1%d1Array(i_lnd_p1) = s% dE_dRho_for_partials(k+1)*s% rho(k+1)
546 0 : e_p1%d1Array(i_lnT_p1) = s% Cv_for_partials(k+1)*s% T(k+1)
547 : end if
548 0 : end function wrap_e_p1
549 :
550 476624 : function wrap_Peos_m1(s, k) result(Peos_m1)
551 : type (star_info), pointer :: s
552 : type(auto_diff_real_star_order1) :: Peos_m1
553 : integer, intent(in) :: k
554 476624 : Peos_m1 = 0d0
555 476624 : if (k > 1) then
556 476624 : Peos_m1%val = s% Peos(k-1)
557 476624 : Peos_m1%d1Array(i_lnd_m1) = s%Peos(k-1) * s% chiRho_for_partials(k-1)
558 476624 : Peos_m1%d1Array(i_lnT_m1) = s%Peos(k-1) * s% chiT_for_partials(k-1)
559 : end if
560 476624 : end function wrap_Peos_m1
561 :
562 581540 : function wrap_Peos_00(s, k) result(Peos_00)
563 : type (star_info), pointer :: s
564 : type(auto_diff_real_star_order1) :: Peos_00
565 : integer, intent(in) :: k
566 581540 : Peos_00 = 0d0
567 581540 : Peos_00%val = s% Peos(k)
568 581540 : Peos_00%d1Array(i_lnd_00) = s%Peos(k) * s% chiRho_for_partials(k)
569 581540 : Peos_00%d1Array(i_lnT_00) = s%Peos(k) * s% chiT_for_partials(k)
570 581540 : end function wrap_Peos_00
571 :
572 0 : function wrap_Peos_p1(s, k) result(Peos_p1)
573 : type (star_info), pointer :: s
574 : type(auto_diff_real_star_order1) :: Peos_p1
575 : integer, intent(in) :: k
576 0 : Peos_p1 = 0d0
577 0 : if (k < s%nz) then
578 0 : Peos_p1%val = s% Peos(k+1)
579 0 : Peos_p1%d1Array(i_lnd_p1) = s%Peos(k+1) * s% chiRho_for_partials(k+1)
580 0 : Peos_p1%d1Array(i_lnT_p1) = s%Peos(k+1) * s% chiT_for_partials(k+1)
581 : end if
582 0 : end function wrap_Peos_p1
583 :
584 0 : function wrap_lnPeos_m1(s, k) result(lnPeos_m1)
585 : type (star_info), pointer :: s
586 : type(auto_diff_real_star_order1) :: lnPeos_m1
587 : integer, intent(in) :: k
588 0 : lnPeos_m1 = 0d0
589 0 : if (k > 1) then
590 0 : lnPeos_m1%val = s% lnPeos(k-1)
591 0 : lnPeos_m1%d1Array(i_lnd_m1) = s% chiRho_for_partials(k-1)
592 0 : lnPeos_m1%d1Array(i_lnT_m1) = s% chiT_for_partials(k-1)
593 : end if
594 0 : end function wrap_lnPeos_m1
595 :
596 44 : function wrap_lnPeos_00(s, k) result(lnPeos_00)
597 : type (star_info), pointer :: s
598 : type(auto_diff_real_star_order1) :: lnPeos_00
599 : integer, intent(in) :: k
600 44 : lnPeos_00 = 0d0
601 44 : lnPeos_00%val = s% lnPeos(k)
602 44 : lnPeos_00%d1Array(i_lnd_00) = s% chiRho_for_partials(k)
603 44 : lnPeos_00%d1Array(i_lnT_00) = s% chiT_for_partials(k)
604 44 : end function wrap_lnPeos_00
605 :
606 0 : function wrap_lnPeos_p1(s, k) result(lnPeos_p1)
607 : type (star_info), pointer :: s
608 : type(auto_diff_real_star_order1) :: lnPeos_p1
609 : integer, intent(in) :: k
610 0 : lnPeos_p1 = 0d0
611 0 : if (k < s%nz) then
612 0 : lnPeos_p1%val = s% lnPeos(k+1)
613 0 : lnPeos_p1%d1Array(i_lnd_p1) = s% chiRho_for_partials(k+1)
614 0 : lnPeos_p1%d1Array(i_lnT_p1) = s% chiT_for_partials(k+1)
615 : end if
616 0 : end function wrap_lnPeos_p1
617 :
618 79928 : function wrap_ChiRho_m1(s, k) result(ChiRho_m1)
619 : use eos_def, only: i_ChiRho
620 : type (star_info), pointer :: s
621 : type(auto_diff_real_star_order1) :: ChiRho_m1
622 : integer, intent(in) :: k
623 79928 : ChiRho_m1 = 0d0
624 79928 : if (k > 1) then
625 79928 : ChiRho_m1%val = s% ChiRho(k-1)
626 79928 : ChiRho_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_ChiRho,k-1)
627 79928 : ChiRho_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_ChiRho,k-1)
628 : end if
629 79928 : end function wrap_ChiRho_m1
630 :
631 79994 : function wrap_ChiRho_00(s, k) result(ChiRho_00)
632 79928 : use eos_def, only: i_ChiRho
633 : type (star_info), pointer :: s
634 : type(auto_diff_real_star_order1) :: ChiRho_00
635 : integer, intent(in) :: k
636 79994 : ChiRho_00 = 0d0
637 79994 : ChiRho_00%val = s% ChiRho(k)
638 79994 : ChiRho_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_ChiRho,k)
639 79994 : ChiRho_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_ChiRho,k)
640 79994 : end function wrap_ChiRho_00
641 :
642 0 : function wrap_ChiRho_p1(s, k) result(ChiRho_p1)
643 79994 : use eos_def, only: i_ChiRho
644 : type (star_info), pointer :: s
645 : type(auto_diff_real_star_order1) :: ChiRho_p1
646 : integer, intent(in) :: k
647 0 : ChiRho_p1 = 0d0
648 0 : if (k < s% nz) then
649 0 : ChiRho_p1%val = s% ChiRho(k+1)
650 0 : ChiRho_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_ChiRho,k+1)
651 0 : ChiRho_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_ChiRho,k+1)
652 : end if
653 0 : end function wrap_ChiRho_p1
654 :
655 79928 : function wrap_ChiT_m1(s, k) result(ChiT_m1)
656 0 : use eos_def, only: i_ChiT
657 : type (star_info), pointer :: s
658 : type(auto_diff_real_star_order1) :: ChiT_m1
659 : integer, intent(in) :: k
660 79928 : ChiT_m1 = 0d0
661 79928 : if (k > 1) then
662 79928 : ChiT_m1%val = s% ChiT(k-1)
663 79928 : ChiT_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_ChiT,k-1)
664 79928 : ChiT_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_ChiT,k-1)
665 : end if
666 79928 : end function wrap_ChiT_m1
667 :
668 79994 : function wrap_ChiT_00(s, k) result(ChiT_00)
669 79928 : use eos_def, only: i_ChiT
670 : type (star_info), pointer :: s
671 : type(auto_diff_real_star_order1) :: ChiT_00
672 : integer, intent(in) :: k
673 79994 : ChiT_00 = 0d0
674 79994 : ChiT_00%val = s% ChiT(k)
675 79994 : ChiT_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_ChiT,k)
676 79994 : ChiT_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_ChiT,k)
677 79994 : end function wrap_ChiT_00
678 :
679 0 : function wrap_ChiT_p1(s, k) result(ChiT_p1)
680 79994 : use eos_def, only: i_ChiT
681 : type (star_info), pointer :: s
682 : type(auto_diff_real_star_order1) :: ChiT_p1
683 : integer, intent(in) :: k
684 0 : ChiT_p1 = 0d0
685 0 : if (k < s% nz) then
686 0 : ChiT_p1%val = s% ChiT(k+1)
687 0 : ChiT_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_ChiT,k+1)
688 0 : ChiT_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_ChiT,k+1)
689 : end if
690 0 : end function wrap_ChiT_p1
691 :
692 79928 : function wrap_Cp_m1(s, k) result(Cp_m1)
693 0 : use eos_def, only: i_Cp
694 : type (star_info), pointer :: s
695 : type(auto_diff_real_star_order1) :: Cp_m1
696 : integer, intent(in) :: k
697 79928 : Cp_m1 = 0d0
698 79928 : if (k > 1) then
699 79928 : Cp_m1%val = s% Cp(k-1)
700 79928 : Cp_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_Cp,k-1)
701 79928 : Cp_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_Cp,k-1)
702 : end if
703 79928 : end function wrap_Cp_m1
704 :
705 79994 : function wrap_Cp_00(s, k) result(Cp_00)
706 79928 : use eos_def, only: i_Cp
707 : type (star_info), pointer :: s
708 : type(auto_diff_real_star_order1) :: Cp_00
709 : integer, intent(in) :: k
710 79994 : Cp_00 = 0d0
711 79994 : Cp_00%val = s% Cp(k)
712 79994 : Cp_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_Cp,k)
713 79994 : Cp_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_Cp,k)
714 79994 : end function wrap_Cp_00
715 :
716 0 : function wrap_Cp_p1(s, k) result(Cp_p1)
717 79994 : use eos_def, only: i_Cp
718 : type (star_info), pointer :: s
719 : type(auto_diff_real_star_order1) :: Cp_p1
720 : integer, intent(in) :: k
721 0 : Cp_p1 = 0d0
722 0 : if (k < s% nz) then
723 0 : Cp_p1%val = s% Cp(k+1)
724 0 : Cp_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_Cp,k+1)
725 0 : Cp_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_Cp,k+1)
726 : end if
727 0 : end function wrap_Cp_p1
728 :
729 79928 : function wrap_grad_ad_m1(s, k) result(grad_ad_m1)
730 0 : use eos_def, only: i_grad_ad
731 : type (star_info), pointer :: s
732 : type(auto_diff_real_star_order1) :: grad_ad_m1
733 : integer, intent(in) :: k
734 79928 : grad_ad_m1 = 0d0
735 79928 : if (k > 1) then
736 79928 : grad_ad_m1%val = s% grada(k-1)
737 79928 : grad_ad_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_grad_ad,k-1)
738 79928 : grad_ad_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_grad_ad,k-1)
739 : end if
740 79928 : end function wrap_grad_ad_m1
741 :
742 79994 : function wrap_grad_ad_00(s, k) result(grad_ad_00)
743 79928 : use eos_def, only: i_grad_ad
744 : type (star_info), pointer :: s
745 : type(auto_diff_real_star_order1) :: grad_ad_00
746 : integer, intent(in) :: k
747 79994 : grad_ad_00 = 0d0
748 79994 : grad_ad_00%val = s% grada(k)
749 79994 : grad_ad_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_grad_ad,k)
750 79994 : grad_ad_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_grad_ad,k)
751 79994 : end function wrap_grad_ad_00
752 :
753 0 : function wrap_grad_ad_p1(s, k) result(grad_ad_p1)
754 79994 : use eos_def, only: i_grad_ad
755 : type (star_info), pointer :: s
756 : type(auto_diff_real_star_order1) :: grad_ad_p1
757 : integer, intent(in) :: k
758 0 : grad_ad_p1 = 0d0
759 0 : if (k < s% nz) then
760 0 : grad_ad_p1%val = s% grada(k+1)
761 0 : grad_ad_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_grad_ad,k+1)
762 0 : grad_ad_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_grad_ad,k+1)
763 : end if
764 0 : end function wrap_grad_ad_p1
765 :
766 0 : function wrap_gamma1_m1(s, k) result(gamma1_m1)
767 0 : use eos_def, only: i_gamma1
768 : type (star_info), pointer :: s
769 : type(auto_diff_real_star_order1) :: gamma1_m1
770 : integer, intent(in) :: k
771 0 : gamma1_m1 = 0d0
772 0 : if (k > 1) then
773 0 : gamma1_m1%val = s% gamma1(k-1)
774 0 : gamma1_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_gamma1,k-1)
775 0 : gamma1_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_gamma1,k-1)
776 : end if
777 0 : end function wrap_gamma1_m1
778 :
779 0 : function wrap_gamma1_00(s, k) result(gamma1_00)
780 0 : use eos_def, only: i_gamma1
781 : type (star_info), pointer :: s
782 : type(auto_diff_real_star_order1) :: gamma1_00
783 : integer, intent(in) :: k
784 0 : gamma1_00 = 0d0
785 0 : gamma1_00%val = s% gamma1(k)
786 0 : gamma1_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_gamma1,k)
787 0 : gamma1_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_gamma1,k)
788 0 : end function wrap_gamma1_00
789 :
790 0 : function wrap_gamma1_p1(s, k) result(gamma1_p1)
791 0 : use eos_def, only: i_gamma1
792 : type (star_info), pointer :: s
793 : type(auto_diff_real_star_order1) :: gamma1_p1
794 : integer, intent(in) :: k
795 0 : gamma1_p1 = 0d0
796 0 : if (k < s% nz) then
797 0 : gamma1_p1%val = s% gamma1(k+1)
798 0 : gamma1_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_gamma1,k+1)
799 0 : gamma1_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_gamma1,k+1)
800 : end if
801 0 : end function wrap_gamma1_p1
802 :
803 0 : function wrap_latent_ddlnT_m1(s, k) result(latent_ddlnT_m1)
804 0 : use eos_def, only: i_latent_ddlnT
805 : type (star_info), pointer :: s
806 : type(auto_diff_real_star_order1) :: latent_ddlnT_m1
807 : integer, intent(in) :: k
808 0 : latent_ddlnT_m1 = 0d0
809 0 : if (k > 1) then
810 0 : latent_ddlnT_m1%val = s% latent_ddlnT(k-1)
811 0 : latent_ddlnT_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_latent_ddlnT,k-1)
812 0 : latent_ddlnT_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_latent_ddlnT,k-1)
813 : end if
814 0 : end function wrap_latent_ddlnT_m1
815 :
816 0 : function wrap_latent_ddlnT_00(s, k) result(latent_ddlnT_00)
817 0 : use eos_def, only: i_latent_ddlnT
818 : type (star_info), pointer :: s
819 : type(auto_diff_real_star_order1) :: latent_ddlnT_00
820 : integer, intent(in) :: k
821 0 : latent_ddlnT_00 = 0d0
822 0 : latent_ddlnT_00%val = s% latent_ddlnT(k)
823 0 : latent_ddlnT_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_latent_ddlnT,k)
824 0 : latent_ddlnT_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_latent_ddlnT,k)
825 0 : end function wrap_latent_ddlnT_00
826 :
827 0 : function wrap_latent_ddlnT_p1(s, k) result(latent_ddlnT_p1)
828 0 : use eos_def, only: i_latent_ddlnT
829 : type (star_info), pointer :: s
830 : type(auto_diff_real_star_order1) :: latent_ddlnT_p1
831 : integer, intent(in) :: k
832 0 : latent_ddlnT_p1 = 0d0
833 0 : if (k < s% nz) then
834 0 : latent_ddlnT_p1%val = s% latent_ddlnT(k+1)
835 0 : latent_ddlnT_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_latent_ddlnT,k+1)
836 0 : latent_ddlnT_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_latent_ddlnT,k+1)
837 : end if
838 0 : end function wrap_latent_ddlnT_p1
839 :
840 0 : function wrap_latent_ddlnRho_m1(s, k) result(latent_ddlnRho_m1)
841 0 : use eos_def, only: i_latent_ddlnRho
842 : type (star_info), pointer :: s
843 : type(auto_diff_real_star_order1) :: latent_ddlnRho_m1
844 : integer, intent(in) :: k
845 0 : latent_ddlnRho_m1 = 0d0
846 0 : if (k > 1) then
847 0 : latent_ddlnRho_m1%val = s% latent_ddlnRho(k-1)
848 0 : latent_ddlnRho_m1%d1Array(i_lnd_m1) = s% d_eos_dlnd(i_latent_ddlnRho,k-1)
849 0 : latent_ddlnRho_m1%d1Array(i_lnT_m1) = s% d_eos_dlnT(i_latent_ddlnRho,k-1)
850 : end if
851 0 : end function wrap_latent_ddlnRho_m1
852 :
853 0 : function wrap_latent_ddlnRho_00(s, k) result(latent_ddlnRho_00)
854 0 : use eos_def, only: i_latent_ddlnRho
855 : type (star_info), pointer :: s
856 : type(auto_diff_real_star_order1) :: latent_ddlnRho_00
857 : integer, intent(in) :: k
858 0 : latent_ddlnRho_00 = 0d0
859 0 : latent_ddlnRho_00%val = s% latent_ddlnRho(k)
860 0 : latent_ddlnRho_00%d1Array(i_lnd_00) = s% d_eos_dlnd(i_latent_ddlnRho,k)
861 0 : latent_ddlnRho_00%d1Array(i_lnT_00) = s% d_eos_dlnT(i_latent_ddlnRho,k)
862 0 : end function wrap_latent_ddlnRho_00
863 :
864 0 : function wrap_latent_ddlnRho_p1(s, k) result(latent_ddlnRho_p1)
865 0 : use eos_def, only: i_latent_ddlnRho
866 : type (star_info), pointer :: s
867 : type(auto_diff_real_star_order1) :: latent_ddlnRho_p1
868 : integer, intent(in) :: k
869 0 : latent_ddlnRho_p1 = 0d0
870 0 : if (k < s% nz) then
871 0 : latent_ddlnRho_p1%val = s% latent_ddlnRho(k+1)
872 0 : latent_ddlnRho_p1%d1Array(i_lnd_p1) = s% d_eos_dlnd(i_latent_ddlnRho,k+1)
873 0 : latent_ddlnRho_p1%d1Array(i_lnT_p1) = s% d_eos_dlnT(i_latent_ddlnRho,k+1)
874 : end if
875 0 : end function wrap_latent_ddlnRho_p1
876 :
877 0 : function wrap_L_m1(s, k) result(L_m1)
878 : type (star_info), pointer :: s
879 : type(auto_diff_real_star_order1) :: L_m1
880 : integer, intent(in) :: k
881 0 : L_m1 = 0d0
882 0 : if (k > 1) then
883 0 : L_m1 % val = s%L(k-1)
884 0 : L_m1 % d1Array(i_L_m1) = 1d0
885 : end if
886 0 : end function wrap_L_m1
887 :
888 212336 : function wrap_L_00(s, k) result(L_00)
889 : type (star_info), pointer :: s
890 : type(auto_diff_real_star_order1) :: L_00
891 : integer, intent(in) :: k
892 212336 : L_00 = 0d0
893 212336 : L_00 % val = s%L(k)
894 212336 : L_00 % d1Array(i_L_00) = 1d0
895 212336 : end function wrap_L_00
896 :
897 52348 : function wrap_L_p1(s, k) result(L_p1)
898 : type (star_info), pointer :: s
899 : type(auto_diff_real_star_order1) :: L_p1
900 : integer, intent(in) :: k
901 52348 : L_p1 = 0d0
902 52348 : if (k < s%nz) then
903 52304 : L_p1 % val = s%L(k+1)
904 52304 : L_p1 % d1Array(i_L_p1) = 1d0
905 : else
906 44 : L_p1 %val = s%L_center
907 : ! L_center is a constant
908 : end if
909 52348 : end function wrap_L_p1
910 :
911 0 : function wrap_r_m1(s, k) result(r_m1)
912 : type (star_info), pointer :: s
913 : type(auto_diff_real_star_order1) :: r_m1
914 : integer, intent(in) :: k
915 0 : r_m1 = 0d0
916 0 : if (k > 1) then
917 0 : r_m1 % val = s%r(k-1)
918 0 : r_m1 % d1Array(i_lnR_m1) = s%r(k-1)
919 : end if
920 0 : end function wrap_r_m1
921 :
922 449242 : function wrap_r_00(s, k) result(r_00)
923 : type (star_info), pointer :: s
924 : type(auto_diff_real_star_order1) :: r_00
925 : integer, intent(in) :: k
926 449242 : r_00 = 0d0
927 449242 : r_00 % val = s%r(k)
928 449242 : r_00 % d1Array(i_lnR_00) = s%r(k)
929 449242 : end function wrap_r_00
930 :
931 52348 : function wrap_r_p1(s, k) result(r_p1)
932 : type (star_info), pointer :: s
933 : type(auto_diff_real_star_order1) :: r_p1
934 : integer, intent(in) :: k
935 52348 : r_p1 = 0d0
936 52348 : if (k < s%nz) then
937 52304 : r_p1 % val = s%r(k+1)
938 52304 : r_p1 % d1Array(i_lnR_p1) = s%r(k+1)
939 : else
940 44 : r_p1 % val = s%r_center
941 : end if
942 52348 : end function wrap_r_p1
943 :
944 0 : function wrap_lnR_m1(s, k) result(lnR_m1)
945 : type (star_info), pointer :: s
946 : type(auto_diff_real_star_order1) :: lnR_m1
947 : integer, intent(in) :: k
948 0 : lnR_m1 = 0d0
949 0 : if (k > 1) then
950 0 : lnR_m1 % val = s%lnR(k-1)
951 0 : lnR_m1 % d1Array(i_lnR_m1) = 1d0
952 : end if
953 0 : end function wrap_lnR_m1
954 :
955 52348 : function wrap_lnR_00(s, k) result(lnR_00)
956 : type (star_info), pointer :: s
957 : type(auto_diff_real_star_order1) :: lnR_00
958 : integer, intent(in) :: k
959 52348 : lnR_00 = 0d0
960 52348 : lnR_00 % val = s%lnR(k)
961 52348 : lnR_00 % d1Array(i_lnR_00) = 1d0
962 52348 : end function wrap_lnR_00
963 :
964 0 : function wrap_lnR_p1(s, k) result(lnR_p1)
965 : type (star_info), pointer :: s
966 : type(auto_diff_real_star_order1) :: lnR_p1
967 : integer, intent(in) :: k
968 0 : lnR_p1 = 0d0
969 0 : if (k < s%nz) then
970 0 : lnR_p1 % val = s%lnR(k+1)
971 0 : lnR_p1 % d1Array(i_lnR_p1) = 1d0
972 : else
973 0 : lnR_p1 % val = log(max(1d0,s%r_center))
974 : end if
975 0 : end function wrap_lnR_p1
976 :
977 0 : function wrap_dxh_lnR(s, k) result(dxh_lnR)
978 : type (star_info), pointer :: s
979 : type(auto_diff_real_star_order1) :: dxh_lnR
980 : integer, intent(in) :: k
981 0 : dxh_lnR = 0d0
982 0 : dxh_lnR % val = s%dxh_lnR(k)
983 0 : dxh_lnR % d1Array(i_lnR_00) = 1d0
984 0 : end function wrap_dxh_lnR
985 :
986 0 : function wrap_dxh_v_face(s, k) result(dxh_v) ! ! wrap_dxh_v_face_00 technically
987 : type (star_info), pointer :: s
988 : type(auto_diff_real_star_order1) :: dxh_v
989 : integer, intent(in) :: k
990 0 : if (s% u_flag) then
991 0 : dxh_v = wrap_dxh_u_face(s,k) ! wrap_dxh_u_face_00 technically
992 0 : return
993 : end if
994 0 : dxh_v = 0d0
995 0 : dxh_v % val = s% dxh_v(k)
996 0 : dxh_v % d1Array(i_v_00) = 1d0
997 0 : end function wrap_dxh_v_face
998 :
999 :
1000 0 : function wrap_geff_face(s, k) result(geff)
1001 : type (star_info), pointer :: s
1002 : type(auto_diff_real_star_order1) :: geff, r2
1003 : integer, intent(in) :: k
1004 0 : geff = 0d0
1005 0 : r2 = 0d0
1006 0 : if (s% include_mlt_in_velocity_time_centering) then
1007 0 : r2 = pow2(wrap_opt_time_center_r_00(s,k))
1008 : else
1009 0 : r2 = pow2(wrap_r_00(s,k))
1010 : end if
1011 :
1012 0 : if (s% rotation_flag .and. s% use_gravity_rotation_correction) then
1013 0 : geff = s%fp_rot(k)*s%cgrav(k)*s%m_grav(k)/r2
1014 : else
1015 0 : geff = s%cgrav(k)*s%m_grav(k)/r2
1016 : end if
1017 :
1018 0 : end function wrap_geff_face
1019 :
1020 0 : subroutine get_face_weights_ad_support(s, k, alfa, beta)
1021 : type (star_info), pointer :: s
1022 : integer, intent(in) :: k
1023 : real(dp), intent(out) :: alfa, beta
1024 : ! face_value(k) = alfa*cell_value(k) + beta*cell_value(k-1)
1025 0 : if (k == 1) call mesa_error(__FILE__,__LINE__,'bad k==1 for get_face_weights')
1026 0 : alfa = s% dq(k-1)/(s% dq(k-1) + s% dq(k))
1027 0 : beta = 1d0 - alfa
1028 0 : end subroutine get_face_weights_ad_support
1029 :
1030 0 : function wrap_dxh_u_face(s, k) result(dxh_u_face) ! wrap_dxh_u_face_00 technically
1031 : type (star_info), pointer :: s
1032 : integer, intent(in) :: k
1033 : type(auto_diff_real_star_order1) :: dxh_u_face
1034 0 : real(dp) :: alpha, beta
1035 :
1036 0 : dxh_u_face = 0d0
1037 0 : if (k == 1) then
1038 : ! at the inner boundary just take the first cell‐update
1039 0 : dxh_u_face%val = s%dxh_u(1)
1040 0 : dxh_u_face%d1Array(i_v_00) = 1d0
1041 : else
1042 : ! get mass weighted face interpolation coefficients
1043 0 : call get_face_weights_ad_support(s, k, alpha, beta)
1044 0 : dxh_u_face%val = alpha*s%dxh_u(k) + beta*s%dxh_u(k-1)
1045 : ! derivatives wrt the two neighbouring dxh_u's
1046 0 : dxh_u_face%d1Array(i_v_00) = alpha
1047 0 : dxh_u_face%d1Array(i_v_m1) = beta
1048 : end if
1049 0 : end function wrap_dxh_u_face
1050 :
1051 :
1052 0 : function wrap_u_face_m1(s, k) result(v_m1)
1053 : type (star_info), pointer :: s
1054 : type(auto_diff_real_star_order1) :: v_m1
1055 : integer, intent(in) :: k
1056 0 : v_m1 = 0
1057 0 : if (k > 1) v_m1 = shift_m1(s% u_face_ad(k-1))
1058 0 : end function wrap_u_face_m1
1059 :
1060 0 : function wrap_u_face_00(s, k) result(v_00)
1061 : type (star_info), pointer :: s
1062 : type(auto_diff_real_star_order1) :: v_00
1063 : integer, intent(in) :: k
1064 0 : v_00 = s% u_face_ad(k)
1065 0 : end function wrap_u_face_00
1066 :
1067 0 : function wrap_u_face_p1(s, k) result(v_p1)
1068 : type (star_info), pointer :: s
1069 : type(auto_diff_real_star_order1) :: v_p1
1070 : integer, intent(in) :: k
1071 0 : v_p1 = 0
1072 0 : if (k < s% nz) v_p1 = shift_p1(s% u_face_ad(k+1))
1073 0 : end function wrap_u_face_p1
1074 :
1075 0 : function wrap_v_m1(s, k) result(v_m1)
1076 : type (star_info), pointer :: s
1077 : type(auto_diff_real_star_order1) :: v_m1
1078 : integer, intent(in) :: k
1079 0 : if (s% u_flag) then
1080 0 : v_m1 = wrap_u_face_m1(s,k)
1081 0 : return
1082 : end if
1083 0 : v_m1 = 0d0
1084 0 : if (k > 1) then
1085 0 : v_m1 % val = s%v(k-1)
1086 0 : v_m1 % d1Array(i_v_m1) = 1d0
1087 : end if
1088 0 : end function wrap_v_m1
1089 :
1090 0 : function wrap_v_00(s, k) result(v_00)
1091 : type (star_info), pointer :: s
1092 : type(auto_diff_real_star_order1) :: v_00
1093 : integer, intent(in) :: k
1094 0 : if (s% u_flag) then
1095 0 : v_00 = wrap_u_face_00(s,k)
1096 0 : return
1097 : end if
1098 0 : v_00 = 0d0
1099 0 : v_00 % val = s%v(k)
1100 0 : v_00 % d1Array(i_v_00) = 1d0
1101 0 : end function wrap_v_00
1102 :
1103 0 : function wrap_v_p1(s, k) result(v_p1)
1104 : type (star_info), pointer :: s
1105 : type(auto_diff_real_star_order1) :: v_p1
1106 : integer, intent(in) :: k
1107 0 : if (s% u_flag) then
1108 0 : v_p1 = wrap_u_face_p1(s,k)
1109 0 : return
1110 : end if
1111 0 : v_p1 = 0d0
1112 0 : if (k < s%nz) then
1113 0 : v_p1 % val = s%v(k+1)
1114 0 : v_p1 % d1Array(i_v_p1) = 1d0
1115 : else
1116 0 : v_p1 % val = s%v_center
1117 : ! v_center is a constant
1118 : end if
1119 0 : end function wrap_v_p1
1120 :
1121 0 : function wrap_opt_time_center_r_m1(s, k) result(r_tc)
1122 : type (star_info), pointer :: s
1123 : integer, intent(in) :: k
1124 : type(auto_diff_real_star_order1) :: r_tc
1125 0 : r_tc = wrap_r_m1(s,k)
1126 0 : if (s% using_velocity_time_centering) then
1127 0 : if (k > 1) r_tc = 0.5d0*(r_tc + s% r_start(k-1))
1128 : end if
1129 0 : end function wrap_opt_time_center_r_m1
1130 :
1131 0 : function wrap_opt_time_center_r_00(s, k) result(r_tc)
1132 : type (star_info), pointer :: s
1133 : integer, intent(in) :: k
1134 : type(auto_diff_real_star_order1) :: r_tc
1135 0 : r_tc = wrap_r_00(s,k)
1136 0 : if (s% using_velocity_time_centering) &
1137 0 : r_tc = 0.5d0*(r_tc + s% r_start(k))
1138 0 : end function wrap_opt_time_center_r_00
1139 :
1140 0 : function wrap_opt_time_center_r_p1(s, k) result(r_tc)
1141 : type (star_info), pointer :: s
1142 : integer, intent(in) :: k
1143 : type(auto_diff_real_star_order1) :: r_tc
1144 0 : r_tc = wrap_r_p1(s,k)
1145 0 : if (s% using_velocity_time_centering) then
1146 0 : if (k < s% nz) then
1147 0 : r_tc = 0.5d0*(r_tc + s% r_start(k+1))
1148 : else
1149 0 : r_tc = 0.5d0*(r_tc + s% r_center)
1150 : end if
1151 : end if
1152 0 : end function wrap_opt_time_center_r_p1
1153 :
1154 0 : function wrap_opt_time_center_v_m1(s, k) result(v_tc)
1155 : type (star_info), pointer :: s
1156 : integer, intent(in) :: k
1157 : type(auto_diff_real_star_order1) :: v_tc
1158 0 : v_tc = 0
1159 0 : if (k == 1) return
1160 0 : if (s% v_flag) then
1161 0 : v_tc = wrap_v_m1(s,k)
1162 0 : if (s% using_velocity_time_centering) &
1163 0 : v_tc = 0.5d0*(v_tc + s% v_start(k-1))
1164 0 : else if (s% u_flag) then
1165 0 : v_tc = wrap_u_face_m1(s,k)
1166 0 : if (s% using_velocity_time_centering) &
1167 0 : v_tc = 0.5d0*(v_tc + s% u_face_start(k-1))
1168 : end if
1169 0 : end function wrap_opt_time_center_v_m1
1170 :
1171 0 : function wrap_opt_time_center_v_00(s, k) result(v_tc)
1172 : type (star_info), pointer :: s
1173 : integer, intent(in) :: k
1174 : type(auto_diff_real_star_order1) :: v_tc
1175 0 : v_tc = 0
1176 0 : if (s% v_flag) then
1177 0 : v_tc = wrap_v_00(s,k)
1178 0 : if (s% using_velocity_time_centering) &
1179 0 : v_tc = 0.5d0*(v_tc + s% v_start(k))
1180 0 : else if (s% u_flag) then
1181 0 : v_tc = wrap_u_face_00(s,k)
1182 0 : if (s% using_velocity_time_centering) &
1183 0 : v_tc = 0.5d0*(v_tc + s% u_face_start(k))
1184 : end if
1185 0 : end function wrap_opt_time_center_v_00
1186 :
1187 0 : function wrap_opt_time_center_v_p1(s, k) result(v_tc)
1188 : type (star_info), pointer :: s
1189 : integer, intent(in) :: k
1190 : type(auto_diff_real_star_order1) :: v_tc
1191 0 : v_tc = 0
1192 0 : if (k == s% nz) return
1193 0 : if (s% v_flag) then
1194 0 : v_tc = wrap_v_p1(s,k)
1195 0 : if (s% using_velocity_time_centering) &
1196 0 : v_tc = 0.5d0*(v_tc + s% v_start(k+1))
1197 0 : else if (s% u_flag) then
1198 0 : v_tc = wrap_u_face_p1(s,k)
1199 0 : if (s% using_velocity_time_centering) &
1200 0 : v_tc = 0.5d0*(v_tc + s% u_face_start(k+1))
1201 : end if
1202 0 : end function wrap_opt_time_center_v_p1
1203 :
1204 0 : function wrap_u_m1(s, k) result(v_m1)
1205 : type (star_info), pointer :: s
1206 : type(auto_diff_real_star_order1) :: v_m1
1207 : integer, intent(in) :: k
1208 0 : v_m1 = 0d0
1209 0 : if (k > 1) then
1210 0 : v_m1 % val = s%u(k-1)
1211 0 : v_m1 % d1Array(i_v_m1) = 1d0
1212 : end if
1213 0 : end function wrap_u_m1
1214 :
1215 0 : function wrap_u_00(s, k) result(v_00)
1216 : type (star_info), pointer :: s
1217 : type(auto_diff_real_star_order1) :: v_00
1218 : integer, intent(in) :: k
1219 0 : v_00 = 0d0
1220 0 : v_00 % val = s%u(k)
1221 0 : v_00 % d1Array(i_v_00) = 1d0
1222 0 : end function wrap_u_00
1223 :
1224 0 : function wrap_u_p1(s, k) result(v_p1)
1225 : type (star_info), pointer :: s
1226 : type(auto_diff_real_star_order1) :: v_p1
1227 : integer, intent(in) :: k
1228 0 : v_p1 = 0d0
1229 0 : if (k < s%nz) then
1230 0 : v_p1 % val = s%u(k+1)
1231 0 : v_p1 % d1Array(i_v_p1) = 1d0
1232 : else
1233 0 : v_p1 % val = 0d0
1234 : ! v_center is a constant
1235 : end if
1236 0 : end function wrap_u_p1
1237 :
1238 0 : function wrap_Hp_m1(s, k) result(Hp_m1)
1239 : type (star_info), pointer :: s
1240 : type(auto_diff_real_star_order1) :: Hp_m1
1241 : integer, intent(in) :: k
1242 0 : Hp_m1 = 0d0
1243 0 : if (k > 1) then
1244 0 : Hp_m1 % val = s%Hp_face(k-1)
1245 0 : Hp_m1 % d1Array(i_Hp_m1) = 1d0
1246 : end if
1247 0 : end function wrap_Hp_m1
1248 :
1249 0 : function wrap_Hp_00(s, k) result(Hp_00)
1250 : type (star_info), pointer :: s
1251 : type(auto_diff_real_star_order1) :: Hp_00
1252 : integer, intent(in) :: k
1253 0 : Hp_00 = 0d0
1254 0 : Hp_00 % val = s%Hp_face(k)
1255 0 : Hp_00 % d1Array(i_Hp_00) = 1d0
1256 0 : end function wrap_Hp_00
1257 :
1258 0 : function wrap_Hp_p1(s, k) result(Hp_p1)
1259 : type (star_info), pointer :: s
1260 : type(auto_diff_real_star_order1) :: Hp_p1
1261 : integer, intent(in) :: k
1262 0 : Hp_p1 = 0d0
1263 0 : if (k < s%nz) then
1264 0 : Hp_p1 % val = s%Hp_face(k+1)
1265 0 : Hp_p1 % d1Array(i_Hp_p1) = 1d0
1266 : end if
1267 0 : end function wrap_Hp_p1
1268 :
1269 :
1270 0 : function wrap_w_div_wc_m1(s, k) result(w_div_wc_m1)
1271 : type (star_info), pointer :: s
1272 : type(auto_diff_real_star_order1) :: w_div_wc_m1
1273 : integer, intent(in) :: k
1274 0 : w_div_wc_m1 = 0d0
1275 0 : if (k > 1) then
1276 0 : w_div_wc_m1 % val = s% w_div_w_crit_roche(k-1)
1277 0 : w_div_wc_m1 % d1Array(i_w_div_wc_m1) = 1d0
1278 : end if
1279 0 : end function wrap_w_div_wc_m1
1280 :
1281 0 : function wrap_w_div_wc_00(s, k) result(w_div_wc_00)
1282 : type (star_info), pointer :: s
1283 : type(auto_diff_real_star_order1) :: w_div_wc_00
1284 : integer, intent(in) :: k
1285 0 : w_div_wc_00 = 0d0
1286 0 : w_div_wc_00 % val = s% w_div_w_crit_roche(k)
1287 0 : w_div_wc_00 % d1Array(i_w_div_wc_00) = 1d0
1288 0 : end function wrap_w_div_wc_00
1289 :
1290 0 : function wrap_w_div_wc_p1(s, k) result(w_div_wc_p1)
1291 : type (star_info), pointer :: s
1292 : type(auto_diff_real_star_order1) :: w_div_wc_p1
1293 : integer, intent(in) :: k
1294 0 : w_div_wc_p1 = 0d0
1295 0 : if (k < s%nz) then
1296 0 : w_div_wc_p1 % val = s% w_div_w_crit_roche(k+1)
1297 0 : w_div_wc_p1 % d1Array(i_w_div_wc_p1) = 1d0
1298 : end if
1299 0 : end function wrap_w_div_wc_p1
1300 :
1301 :
1302 0 : function wrap_jrot_m1(s, k) result(jrot_m1)
1303 : type (star_info), pointer :: s
1304 : type(auto_diff_real_star_order1) :: jrot_m1
1305 : integer, intent(in) :: k
1306 0 : jrot_m1 = 0d0
1307 0 : if (k > 1) then
1308 0 : jrot_m1 % val = s% j_rot(k-1)
1309 0 : jrot_m1 % d1Array(i_jrot_m1) = 1d0
1310 : end if
1311 0 : end function wrap_jrot_m1
1312 :
1313 0 : function wrap_jrot_00(s, k) result(jrot_00)
1314 : type (star_info), pointer :: s
1315 : type(auto_diff_real_star_order1) :: jrot_00
1316 : integer, intent(in) :: k
1317 0 : jrot_00 = 0d0
1318 0 : jrot_00 % val = s% j_rot(k)
1319 0 : jrot_00 % d1Array(i_jrot_00) = 1d0
1320 0 : end function wrap_jrot_00
1321 :
1322 0 : function wrap_jrot_p1(s, k) result(jrot_p1)
1323 : type (star_info), pointer :: s
1324 : type(auto_diff_real_star_order1) :: jrot_p1
1325 : integer, intent(in) :: k
1326 0 : jrot_p1 = 0d0
1327 0 : if (k < s%nz) then
1328 0 : jrot_p1 % val = s% j_rot(k+1)
1329 0 : jrot_p1 % d1Array(i_jrot_p1) = 1d0
1330 : end if
1331 0 : end function wrap_jrot_p1
1332 :
1333 :
1334 0 : function wrap_omega_m1(s, k) result(omega_m1)
1335 : type (star_info), pointer :: s
1336 : type(auto_diff_real_star_order1) :: omega_m1, jrot_m1
1337 : integer, intent(in) :: k
1338 0 : omega_m1 = 0d0
1339 0 : jrot_m1 = wrap_jrot_m1(s,k)
1340 0 : if (k > 1) then
1341 0 : omega_m1 = jrot_m1/shift_m1(s% i_rot(k-1))
1342 : end if
1343 0 : end function wrap_omega_m1
1344 :
1345 0 : function wrap_omega_00(s, k) result(omega_00)
1346 : type (star_info), pointer :: s
1347 : type(auto_diff_real_star_order1) :: omega_00, jrot_00
1348 : integer, intent(in) :: k
1349 0 : jrot_00 = wrap_jrot_00(s,k)
1350 0 : omega_00 = jrot_00/s% i_rot(k)
1351 0 : end function wrap_omega_00
1352 :
1353 0 : function wrap_omega_p1(s, k) result(omega_p1)
1354 : type (star_info), pointer :: s
1355 : type(auto_diff_real_star_order1) :: omega_p1, jrot_p1
1356 : integer, intent(in) :: k
1357 0 : omega_p1 = 0d0
1358 0 : jrot_p1 = wrap_jrot_p1(s,k)
1359 0 : if (k < s%nz) then
1360 0 : omega_p1 = jrot_p1/shift_p1(s% i_rot(k+1))
1361 : end if
1362 0 : end function wrap_omega_p1
1363 :
1364 :
1365 0 : function wrap_xtra1_m1(s, k) result(xtra1_m1)
1366 : type (star_info), pointer :: s
1367 : type(auto_diff_real_star_order1) :: xtra1_m1
1368 : integer, intent(in) :: k
1369 0 : xtra1_m1 = 0d0
1370 0 : if (k > 1) then
1371 0 : xtra1_m1 % val = 0d0
1372 0 : xtra1_m1 % d1Array(i_xtra1_m1) = 1d0
1373 : end if
1374 0 : end function wrap_xtra1_m1
1375 :
1376 0 : function wrap_xtra1_00(s, k) result(xtra1_00)
1377 : type (star_info), pointer :: s
1378 : type(auto_diff_real_star_order1) :: xtra1_00
1379 : integer, intent(in) :: k
1380 0 : xtra1_00 = 0d0
1381 0 : xtra1_00 % val = 0d0
1382 0 : xtra1_00 % d1Array(i_xtra1_00) = 1d0
1383 0 : end function wrap_xtra1_00
1384 :
1385 0 : function wrap_xtra1_p1(s, k) result(xtra1_p1)
1386 : type (star_info), pointer :: s
1387 : type(auto_diff_real_star_order1) :: xtra1_p1
1388 : integer, intent(in) :: k
1389 0 : xtra1_p1 = 0d0
1390 0 : if (k < s%nz) then
1391 0 : xtra1_p1 % val = 0d0
1392 0 : xtra1_p1 % d1Array(i_xtra1_p1) = 1d0
1393 : end if
1394 0 : end function wrap_xtra1_p1
1395 :
1396 :
1397 0 : function wrap_xtra2_m1(s, k) result(xtra2_m1)
1398 : type (star_info), pointer :: s
1399 : type(auto_diff_real_star_order1) :: xtra2_m1
1400 : integer, intent(in) :: k
1401 0 : xtra2_m1 = 0d0
1402 0 : if (k > 1) then
1403 0 : xtra2_m1 % val = 0d0
1404 0 : xtra2_m1 % d1Array(i_xtra2_m1) = 1d0
1405 : end if
1406 0 : end function wrap_xtra2_m1
1407 :
1408 0 : function wrap_xtra2_00(s, k) result(xtra2_00)
1409 : type (star_info), pointer :: s
1410 : type(auto_diff_real_star_order1) :: xtra2_00
1411 : integer, intent(in) :: k
1412 0 : xtra2_00 = 0d0
1413 0 : xtra2_00 % val = 0d0
1414 0 : xtra2_00 % d1Array(i_xtra2_00) = 1d0
1415 0 : end function wrap_xtra2_00
1416 :
1417 0 : function wrap_xtra2_p1(s, k) result(xtra2_p1)
1418 : type (star_info), pointer :: s
1419 : type(auto_diff_real_star_order1) :: xtra2_p1
1420 : integer, intent(in) :: k
1421 0 : xtra2_p1 = 0d0
1422 0 : if (k < s%nz) then
1423 0 : xtra2_p1 % val = 0d0
1424 0 : xtra2_p1 % d1Array(i_xtra2_p1) = 1d0
1425 : end if
1426 0 : end function wrap_xtra2_p1
1427 :
1428 : end module auto_diff_support
|