-
Notifications
You must be signed in to change notification settings - Fork 46
/
Copy pathhydro_rsp2.f90
1290 lines (1147 loc) · 49.6 KB
/
hydro_rsp2.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! ***********************************************************************
!
! Copyright (C) 2010-2020 The MESA Team
!
! MESA is free software; you can use it and/or modify
! it under the combined terms and restrictions of the MESA MANIFESTO
! and the GNU General Library Public License as published
! by the Free Software Foundation; either version 2 of the License,
! or (at your option) any later version.
!
! You should have received a copy of the MESA MANIFESTO along with
! this software; if not, it is available at the mesa website:
! http://mesa.sourceforge.net/
!
! MESA is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
! See the GNU Library General Public License for more details.
!
! You should have received a copy of the GNU Library General Public License
! along with this software; if not, write to the Free Software
! Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
!
! ***********************************************************************
module hydro_rsp2
use star_private_def
use const_def
use utils_lib, only: is_bad
use auto_diff
use auto_diff_support
use accurate_sum_auto_diff_star_order1
use star_utils
implicit none
private
public :: &
do1_rsp2_L_eqn, do1_turbulent_energy_eqn, do1_rsp2_Hp_eqn, &
compute_Eq_cell, compute_Uq_face, set_RSP2_vars, &
Hp_face_for_rsp2_val, Hp_face_for_rsp2_eqn, set_etrb_start_vars, &
RSP2_adjust_vars_before_call_solver, get_RSP2_alfa_beta_face_weights, &
set_viscosity_vars_TDC
real(dp), parameter :: &
x_ALFAP = 2.d0/3.d0, & ! Ptrb
x_ALFAS = (1.d0/2.d0)*sqrt_2_div_3, & ! PII_face and Lc
x_ALFAC = (1.d0/2.d0)*sqrt_2_div_3, & ! Lc
x_CEDE = (8.d0/3.d0)*sqrt_2_div_3, & ! DAMP
x_GAMMAR = 2.d0*sqrt(3.d0) ! DAMPR
contains
subroutine set_RSP2_vars(s,ierr)
type (star_info), pointer :: s
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: x
integer :: k, op_err
include 'formats'
ierr = 0
op_err = 0
!$OMP PARALLEL DO PRIVATE(k,op_err) SCHEDULE(dynamic,2)
do k=1,s%nz
! Hp_face(k) <= 0 means it needs to be set. e.g., after read file
if (s% Hp_face(k) <= 0) then
s% Hp_face(k) = get_scale_height_face_val(s,k)
s% xh(s% i_Hp,k) = s% Hp_face(k)
end if
x = compute_Y_face(s, k, op_err)
if (op_err /= 0) ierr = op_err
x = compute_PII_face(s, k, op_err)
if (op_err /= 0) ierr = op_err
!Pvsc skip?
end do
!$OMP END PARALLEL DO
if (ierr /= 0) then
if (s% report_ierr) write(*,2) 'failed in set_RSP2_vars loop 1', s% model_number
return
end if
!$OMP PARALLEL DO PRIVATE(k,op_err) SCHEDULE(dynamic,2)
do k=1,s% nz
x = compute_Chi_cell(s, k, op_err)
if (op_err /= 0) ierr = op_err
x = compute_Eq_cell(s, k, op_err)
if (op_err /= 0) ierr = op_err
x = compute_C(s, k, op_err) ! COUPL
if (op_err /= 0) ierr = op_err
x = compute_L_face(s, k, op_err) ! Lr, Lt, Lc
if (op_err /= 0) ierr = op_err
end do
!$OMP END PARALLEL DO
if (ierr /= 0) then
if (s% report_ierr) write(*,2) 'failed in set_RSP2_vars loop 2', s% model_number
return
end if
do k = 1, s% RSP2_num_outermost_cells_forced_nonturbulent
s% Eq(k) = 0d0; s% Eq_ad(k) = 0d0
s% Chi(k) = 0d0; s% Chi_ad(k) = 0d0
s% COUPL(k) = 0d0; s% COUPL_ad(k) = 0d0
!s% Ptrb(k) = 0d0;
s% Lc(k) = 0d0; s% Lc_ad(k) = 0d0
s% Lt(k) = 0d0; s% Lt_ad(k) = 0d0
end do
do k = s% nz + 1 - int(s% nz/s% RSP2_nz_div_IBOTOM) , s% nz
s% Eq(k) = 0d0; s% Eq_ad(k) = 0d0
s% Chi(k) = 0d0; s% Chi_ad(k) = 0d0
s% COUPL(k) = 0d0; s% COUPL_ad(k) = 0d0
!s% Ptrb(k) = 0d0;
s% Lc(k) = 0d0; s% Lc_ad(k) = 0d0
s% Lt(k) = 0d0; s% Lt_ad(k) = 0d0
end do
end subroutine set_RSP2_vars
! This routine is called to initialize eq and uq for TDC.
subroutine set_viscosity_vars_TDC(s,ierr)
type (star_info), pointer :: s
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: x
integer :: k, op_err
include 'formats'
ierr = 0
op_err = 0
!$OMP PARALLEL DO PRIVATE(k,op_err) SCHEDULE(dynamic,2)
do k=1,s%nz
! Hp_face(k) <= 0 means it needs to be set. e.g., after read file
if (s% Hp_face(k) <= 0) then
! this scale height for face is already calculated in TDC
s% Hp_face(k) = s% scale_height(k) !get_scale_height_face_val(s,k)
end if
end do
!$OMP END PARALLEL DO
if (ierr /= 0) then
if (s% report_ierr) write(*,2) 'failed in set_viscosity_vars_TDC loop 1', s% model_number
return
end if
!$OMP PARALLEL DO PRIVATE(k,op_err) SCHEDULE(dynamic,2)
do k=1,s% nz
x = compute_Chi_cell(s, k, op_err)
if (op_err /= 0) ierr = op_err
x = compute_Eq_cell(s, k, op_err)
if (op_err /= 0) ierr = op_err
x = compute_Uq_face(s, k, op_err)
if (op_err /= 0) ierr = op_err
end do
!$OMP END PARALLEL DO
if (ierr /= 0) then
if (s% report_ierr) write(*,2) 'failed in set_viscosity_vars_TDC loop 2', s% model_number
return
end if
if (.not. s% v_flag) then ! set values 0 if not using v_flag.
do k = 1, s% nz
s% Eq(k) = 0d0; s% Eq_ad(k) = 0d0
s% Chi(k) = 0d0; s% Chi_ad(k) = 0d0
s% Uq(k) = 0d0
end do
end if
end subroutine set_viscosity_vars_TDC
subroutine do1_rsp2_L_eqn(s, k, nvar, ierr)
use star_utils, only: save_eqn_residual_info
type (star_info), pointer :: s
integer, intent(in) :: k, nvar
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: &
L_expected, L_actual,resid
real(dp) :: scale, residual, L_start_max
logical :: test_partials
include 'formats'
!test_partials = (k == s% solver_test_partials_k)
test_partials = .false.
if (.not. s% RSP2_flag) then
ierr = -1
return
end if
ierr = 0
!L_expected = compute_L_face(s, k, ierr)
!if (ierr /= 0) return
L_expected = s% Lr_ad(k) + s% Lc_ad(k) + s% Lt_ad(k)
L_actual = wrap_L_00(s, k)
L_start_max = maxval(s% L_start(1:s% nz))
scale = 1d0/L_start_max
if (is_bad(scale)) then
write(*,2) 'do1_rsp2_L_eqn scale', k, scale
call mesa_error(__FILE__,__LINE__,'do1_rsp2_L_eqn')
end if
resid = (L_expected - L_actual)*scale
residual = resid%val
s% equ(s% i_equL, k) = residual
if (test_partials) then
s% solver_test_partials_val = residual
end if
call save_eqn_residual_info(s, k, nvar, s% i_equL, resid, 'do1_rsp2_L_eqn', ierr)
if (ierr /= 0) return
if (test_partials) then
s% solver_test_partials_var = s% i_lnR
s% solver_test_partials_dval_dx = resid%d1Array(i_lnR_00)
write(*,4) 'do1_rsp2_L_eqn', s% solver_test_partials_var
end if
end subroutine do1_rsp2_L_eqn
subroutine do1_rsp2_Hp_eqn(s, k, nvar, ierr)
use star_utils, only: save_eqn_residual_info
type (star_info), pointer :: s
integer, intent(in) :: k, nvar
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: &
Hp_expected, Hp_actual,resid
real(dp) :: residual, Hp_start
logical :: test_partials
include 'formats'
!test_partials = (k == s% solver_test_partials_k)
test_partials = .false.
if (.not. s% RSP2_flag) then
ierr = -1
return
end if
ierr = 0
Hp_expected = Hp_face_for_rsp2_eqn(s, k, ierr)
if (ierr /= 0) return
Hp_actual = wrap_Hp_00(s, k)
Hp_start = s% Hp_face_start(k)
resid = (Hp_expected - Hp_actual)/max(Hp_expected,Hp_actual)
residual = resid%val
s% equ(s% i_equ_Hp, k) = residual
if (test_partials) then
s% solver_test_partials_val = residual
end if
if (residual > 1d3) then
!$omp critical (hydro_rsp2_1)
write(*,2) 'residual', k, residual
write(*,2) 'Hp_expected', k, Hp_expected%val
write(*,2) 'Hp_actual', k, Hp_actual%val
call mesa_error(__FILE__,__LINE__,'do1_rsp2_Hp_eqn')
!$omp end critical (hydro_rsp2_1)
end if
call save_eqn_residual_info(s, k, nvar, s% i_equ_Hp, resid, 'do1_rsp2_Hp_eqn', ierr)
if (ierr /= 0) return
if (test_partials) then
s% solver_test_partials_var = s% i_lnR
s% solver_test_partials_dval_dx = resid%d1Array(i_lnR_00)
write(*,4) 'do1_rsp2_Hp_eqn', s% solver_test_partials_var
end if
end subroutine do1_rsp2_Hp_eqn
real(dp) function Hp_face_for_rsp2_val(s, k, ierr) result(Hp_face) ! cm
type (star_info), pointer :: s
integer, intent(in) :: k
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Hp_face_ad
ierr = 0
Hp_face_ad = Hp_face_for_rsp2_eqn(s, k, ierr)
if (ierr /= 0) return
Hp_face = Hp_face_ad%val
end function Hp_face_for_rsp2_val
function Hp_face_for_rsp2_eqn(s, k, ierr) result(Hp_face) ! cm
type (star_info), pointer :: s
integer, intent(in) :: k
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Hp_face
type(auto_diff_real_star_order1) :: &
rho_face, area, dlnPeos, &
r_00, Peos_00, d_00, Peos_m1, d_m1, Peos_div_rho, &
d_face, Peos_face, alt_Hp_face, A
real(dp) :: alfa, beta
include 'formats'
ierr = 0
if (k > s% nz) then
Hp_face = 1d0 ! not used
return
end if
if (k > 1 .and. .not. s% RSP2_assume_HSE) then
call get_RSP2_alfa_beta_face_weights(s, k, alfa, beta)
rho_face = alfa*wrap_d_00(s,k) + beta*wrap_d_m1(s,k)
area = 4d0*pi*pow2(wrap_r_00(s,k))
dlnPeos = wrap_lnPeos_m1(s,k) - wrap_lnPeos_00(s,k)
Hp_face = -s% dm_bar(k)/(area*rho_face*dlnPeos)
else
r_00 = wrap_r_00(s, k) ! not time-centered in RSP
d_00 = wrap_d_00(s, k)
Peos_00 = wrap_Peos_00(s, k)
if (k == 1) then
Peos_div_rho = Peos_00/d_00
Hp_face = pow2(r_00)*Peos_div_rho/(s% cgrav(k)*s% m(k))
else
d_m1 = wrap_d_m1(s, k)
Peos_m1 = wrap_Peos_m1(s, k)
call get_RSP2_alfa_beta_face_weights(s, k, alfa, beta)
Peos_div_rho = alfa*Peos_00/d_00 + beta*Peos_m1/d_m1
Hp_face = pow2(r_00)*Peos_div_rho/(s% cgrav(k)*s% m(k))
if (k==-104) then
write(*,3) 'RSP2 Hp P_div_rho Pdrho_00 Pdrho_m1', k, s% solver_iter, &
Hp_face%val, Peos_div_rho%val, Peos_00%val/d_00%val, Peos_m1%val/d_m1%val
!write(*,3) 'RSP2 Hp r2_div_Gm r_start r', k, s% solver_iter, &
! Hp_face%val, pow2(r_00%val)/(s% cgrav(k)*s% m(k)), &
! s% r_start(k), r_00%val
end if
if (s% alt_scale_height_flag) then
call mesa_error(__FILE__,__LINE__,'Hp_face_for_rsp2_eqn: cannot use alt_scale_height_flag')
! consider sound speed*hydro time scale as an alternative scale height
d_face = alfa*d_00 + beta*d_m1
Peos_face = alfa*Peos_00 + beta*Peos_m1
alt_Hp_face = sqrt(Peos_face/s% cgrav(k))/d_face
if (alt_Hp_face%val < Hp_face%val) then ! blend
A = pow2(alt_Hp_face/Hp_face) ! 0 <= A%val < 1
Hp_face = A*Hp_face + (1d0 - A)*alt_Hp_face
end if
end if
end if
end if
end function Hp_face_for_rsp2_eqn
subroutine do1_turbulent_energy_eqn(s, k, nvar, ierr)
use star_utils, only: set_energy_eqn_scal, save_eqn_residual_info
type (star_info), pointer :: s
integer, intent(in) :: k, nvar
integer, intent(out) :: ierr
! for OLD WAY
type(auto_diff_real_star_order1) :: &
d_turbulent_energy_ad, Ptrb_dV_ad, dt_C_ad, dt_Eq_ad
type(auto_diff_real_star_order1) :: w_00
type(auto_diff_real_star_order1) :: tst, resid_ad, dt_dLt_dm_ad
type(accurate_auto_diff_real_star_order1) :: esum_ad
logical :: non_turbulent_cell, test_partials
real(dp) :: residual, scal
include 'formats'
!test_partials = (k == s% solver_test_partials_k)
test_partials = .false.
ierr = 0
w_00 = wrap_w_00(s,k)
non_turbulent_cell = &
s% mixing_length_alpha == 0d0 .or. &
k <= s% RSP2_num_outermost_cells_forced_nonturbulent .or. &
k > s% nz - int(s% nz/s% RSP2_nz_div_IBOTOM)
if (.not. s% RSP2_flag) then
resid_ad = w_00 - s% w_start(k) ! just hold w constant when not using RSP2
else if (non_turbulent_cell) then
resid_ad = w_00/s% csound(k) ! make w = 0
else
call setup_d_turbulent_energy(ierr); if (ierr /= 0) return ! erg g^-1 = cm^2 s^-2
call setup_Ptrb_dV_ad(ierr); if (ierr /= 0) return ! erg g^-1
call setup_dt_dLt_dm_ad(ierr); if (ierr /= 0) return ! erg g^-1
call setup_dt_C_ad(ierr); if (ierr /= 0) return ! erg g^-1
call setup_dt_Eq_ad(ierr); if (ierr /= 0) return ! erg g^-1
call set_energy_eqn_scal(s, k, scal, ierr); if (ierr /= 0) return ! 1/(erg g^-1 s^-1)
! sum terms in esum_ad using accurate_auto_diff_real_star_order1
esum_ad = d_turbulent_energy_ad + Ptrb_dV_ad + dt_dLt_dm_ad - dt_C_ad - dt_Eq_ad ! erg g^-1
resid_ad = esum_ad
if (k==-35 .and. s% solver_iter == 1) then
write(*,3) 'RSP2 w dEt PdV dtC dtEq', k, s% solver_iter, &
w_00%val, d_turbulent_energy_ad%val, Ptrb_dV_ad%val, dt_C_ad%val, dt_Eq_ad%val
end if
resid_ad = resid_ad*scal/s%dt ! to make residual unitless, must cancel out the dt in scal
end if
residual = resid_ad%val
s% equ(s% i_detrb_dt, k) = residual
if (test_partials) then
tst = residual
s% solver_test_partials_val = tst%val
if (s% solver_iter == 12) &
write(*,*) 'do1_turbulent_energy_eqn', s% solver_test_partials_var, s% lnd(k), tst%val
end if
call save_eqn_residual_info(s, k, nvar, s% i_detrb_dt, resid_ad, 'do1_turbulent_energy_eqn', ierr)
if (ierr /= 0) return
if (test_partials) then
s% solver_test_partials_var = s% i_lnd
s% solver_test_partials_dval_dx = tst%d1Array(i_lnd_00) ! xi0 good , xi1 partial 0, xi2 good. Af horrible.'
write(*,*) 'do1_turbulent_energy_eqn', s% solver_test_partials_var, s% lnd(k)/ln10, tst%val
end if
contains
subroutine setup_d_turbulent_energy(ierr) ! erg g^-1
integer, intent(out) :: ierr
ierr = 0
d_turbulent_energy_ad = wrap_etrb_00(s,k) - get_etrb_start(s,k)
end subroutine setup_d_turbulent_energy
! Ptrb_dV_ad = Ptrb_ad*dV_ad
subroutine setup_Ptrb_dV_ad(ierr) ! erg g^-1
use star_utils, only: calc_Ptrb_ad_tw
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Ptrb_ad, PT0, dV_ad, d_00
call calc_Ptrb_ad_tw(s, k, Ptrb_ad, PT0, ierr)
if (ierr /= 0) return
d_00 = wrap_d_00(s,k)
dV_ad = 1d0/d_00 - 1d0/s% rho_start(k)
Ptrb_dV_ad = Ptrb_ad*dV_ad ! erg cm^-3 cm^-3 g^-1 = erg g^-1
end subroutine setup_Ptrb_dV_ad
subroutine setup_dt_dLt_dm_ad(ierr) ! erg g^-1
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Lt_00, Lt_p1
real(dp) :: L_theta
include 'formats'
ierr = 0
if (s% using_velocity_time_centering .and. &
s% include_L_in_velocity_time_centering) then
L_theta = s% L_theta_for_velocity_time_centering
else
L_theta = 1d0
end if
Lt_00 = L_theta*s% Lt_ad(k) + (1d0 - L_theta)*s% Lt_start(k)
if (k == s% nz) then
Lt_p1 = 0d0
else
Lt_p1 = L_theta*shift_p1(s% Lt_ad(k+1)) + (1d0 - L_theta)*s% Lt_start(k+1)
if (ierr /= 0) return
end if
dt_dLt_dm_ad = (Lt_00 - Lt_p1)*s%dt/s%dm(k)
end subroutine setup_dt_dLt_dm_ad
subroutine setup_dt_C_ad(ierr) ! erg g^-1
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: C
C = s% COUPL_ad(k) ! compute_C(s, k, ierr) ! erg g^-1 s^-1
if (ierr /= 0) return
dt_C_ad = s%dt*C
end subroutine setup_dt_C_ad
subroutine setup_dt_Eq_ad(ierr) ! erg g^-1
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Eq_cell
Eq_cell = s% Eq_ad(k) ! compute_Eq_cell(s, k, ierr) ! erg g^-1 s^-1
if (ierr /= 0) return
dt_Eq_ad = s%dt*Eq_cell
end subroutine setup_dt_Eq_ad
end subroutine do1_turbulent_energy_eqn
subroutine get_RSP2_alfa_beta_face_weights(s, k, alfa, beta)
type (star_info), pointer :: s
integer, intent(in) :: k
real(dp), intent(out) :: alfa, beta
! face_value = alfa*cell_value(k) + beta*cell_value(k-1)
if (k == 1) call mesa_error(__FILE__,__LINE__,'bad k==1 for get_RSP2_alfa_beta_face_weights')
if (s% RSP2_use_mass_interp_face_values) then
alfa = s% dq(k-1)/(s% dq(k-1) + s% dq(k))
beta = 1d0 - alfa
else
alfa = 0.5d0
beta = 0.5d0
end if
end subroutine get_RSP2_alfa_beta_face_weights
function compute_Y_face(s, k, ierr) result(Y_face) ! superadiabatic gradient [unitless]
type (star_info), pointer :: s
integer, intent(in) :: k
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Y_face
type(auto_diff_real_star_order1) :: Hp_face, Y1, Y2, QQ_div_Cp_face, &
r_00, d_00, Peos_00, Cp_00, T_00, chiT_00, chiRho_00, QQ_00, lnT_00, &
r_m1, d_m1, Peos_m1, Cp_m1, T_m1, chiT_m1, chiRho_m1, QQ_m1, lnT_m1, &
dlnT_dlnP, grad_ad_00, grad_ad_m1, grad_ad_face, dlnT, dlnP, alt_Y_face
real(dp) :: dm_bar, alfa, beta
include 'formats'
ierr = 0
if (k > s% nz) then
Y_face = 0d0
return
end if
if (k == 1 .or. s% mixing_length_alpha == 0d0) then
Y_face = 0d0
s% Y_face(k) = 0d0
s% Y_face_ad(k) = 0d0
return
end if
call get_RSP2_alfa_beta_face_weights(s, k, alfa, beta)
if (s% RSP2_use_RSP_eqn_for_Y_face) then
dm_bar = s% dm_bar(k)
Hp_face = wrap_Hp_00(s,k)
r_00 = wrap_r_00(s, k)
d_00 = wrap_d_00(s, k)
Peos_00 = wrap_Peos_00(s, k)
Cp_00 = wrap_Cp_00(s, k)
T_00 = wrap_T_00(s, k)
chiT_00 = wrap_chiT_00(s, k)
chiRho_00 = wrap_chiRho_00(s, k)
QQ_00 = chiT_00/(d_00*T_00*chiRho_00)
lnT_00 = wrap_lnT_00(s,k)
r_m1 = wrap_r_m1(s, k)
d_m1 = wrap_d_m1(s, k)
Peos_m1 = wrap_Peos_m1(s, k)
Cp_m1 = wrap_Cp_m1(s, k)
T_m1 = wrap_T_m1(s, k)
chiT_m1 = wrap_chiT_m1(s, k)
chiRho_m1 = wrap_chiRho_m1(s, k)
QQ_m1 = chiT_m1/(d_m1*T_m1*chiRho_m1)
lnT_m1 = wrap_lnT_m1(s,k)
QQ_div_Cp_face = alfa*QQ_00/Cp_00 + beta*QQ_m1/Cp_m1
! QQ units (g cm^-3 K)^-1 = g^-1 cm^3 K^-1
! Cp units erg g^-1 K^-1 = g cm^2 s^-2 g^-1 K^-1 = cm^2 s^-2 K^-1
! QQ/Cp units = (g^-1 cm^3 K^-1)/(cm^2 s^-2 K^-1)
! = g^-1 cm^3 K^-1 cm^-2 s^2 K
! = g^-1 cm s^2
! P units = erg cm^-3 = g cm^2 s^-2 cm^-3 = g cm^-1 s^-2
! QQ/Cp*P is unitless.
Y1 = QQ_div_Cp_face*(Peos_m1 - Peos_00) - (lnT_m1 - lnT_00)
! Y1 unitless
Y2 = 4d0*pi*pow2(r_00)*Hp_face*2d0/(1d0/d_00 + 1d0/d_m1)/dm_bar
! units = cm^2 cm / (cm^3 g^-1) / g
! = cm^2 cm cm^-3 g g^-1 = unitless
Y_face = Y1*Y2 ! unitless
if (k==-35) then
write(*,3) 'RSP2 Y_face Y1 Y2', k, s% solver_iter, s% Y_face(k), Y1%val, Y2%val
write(*,3) 'Peos', k, s% solver_iter, Peos_00%val
write(*,3) 'Peos', k-1, s% solver_iter, Peos_m1%val
write(*,3) 'QQ', k, s% solver_iter, QQ_00%val
write(*,3) 'QQ', k-1, s% solver_iter, QQ_m1%val
write(*,3) 'Cp', k, s% solver_iter, Cp_00%val
write(*,3) 'Cp', k-1, s% solver_iter, Cp_m1%val
write(*,3) 'lgT', k, s% solver_iter, lnT_00%val/ln10
write(*,3) 'lgT', k-1, s% solver_iter, lnT_m1%val/ln10
write(*,3) 'lgd', k, s% solver_iter, s% lnd(k)/ln10
write(*,3) 'lgd', k-1, s% solver_iter, s% lnd(k-1)/ln10
!call mesa_error(__FILE__,__LINE__,'compute_Y_face')
end if
else
grad_ad_00 = wrap_grad_ad_00(s,k)
grad_ad_m1 = wrap_grad_ad_m1(s,k)
grad_ad_face = alfa*grad_ad_00 + beta*grad_ad_m1
dlnT = wrap_lnT_m1(s,k) - wrap_lnT_00(s,k)
dlnP = wrap_lnPeos_m1(s,k) - wrap_lnPeos_00(s,k)
dlnT_dlnP = dlnT/dlnP
if (is_bad(dlnT_dlnP%val)) then
alt_Y_face = 0d0
else if (s% use_Ledoux_criterion .and. s% calculate_Brunt_B) then
! gradL = grada + gradL_composition_term
alt_Y_face = dlnT_dlnP - (grad_ad_face + s% gradL_composition_term(k))
else
alt_Y_face = dlnT_dlnP - grad_ad_face
end if
if (is_bad(alt_Y_face%val)) alt_Y_face = 0
Y_face = alt_Y_face
end if
s% Y_face_ad(k) = Y_face
s% Y_face(k) = Y_face%val
end function compute_Y_face
function compute_PII_face(s, k, ierr) result(PII_face) ! ergs g^-1 K^-1 (like Cp)
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: PII_face
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Cp_00, Cp_m1, Cp_face, Y_face
real(dp) :: ALFAS_ALFA, alfa, beta
include 'formats'
ierr = 0
if (k > s% nz) then
PII_face = 0d0
return
end if
if (k == 1 .or. s% mixing_length_alpha == 0d0 .or. &
k == s% nz) then ! just skip k == nz to be like RSP
PII_face = 0d0
s% PII(k) = 0d0
s% PII_ad(k) = 0d0
return
end if
Y_face = s% Y_face_ad(k) ! compute_Y_face(s, k, ierr)
if (ierr /= 0) return
Cp_00 = wrap_Cp_00(s, k)
Cp_m1 = wrap_Cp_m1(s, k)
call get_RSP2_alfa_beta_face_weights(s, k, alfa, beta)
Cp_face = alfa*Cp_00 + beta*Cp_m1 ! ergs g^-1 K^-1
ALFAS_ALFA = x_ALFAS*s% mixing_length_alpha
PII_face = ALFAS_ALFA*Cp_face*Y_face
s% PII(k) = PII_face%val
s% PII_ad(k) = PII_face
if (k == -2 .and. s% PII(k) < 0d0) then
write(*,2) 's% PII(k)', k, s% PII(k)
write(*,2) 'Cp_face', k, Cp_face%val
write(*,2) 'Y_face', k, Y_face%val
!write(*,2) 'PII_face%val', k, PII_face%val
!write(*,2) 'T_rho_face%val', k, T_rho_face%val
!write(*,2) '', k,
!write(*,2) '', k,
call mesa_error(__FILE__,__LINE__,'compute_PII_face')
end if
end function compute_PII_face
function compute_d_v_div_r(s, k, ierr) result(d_v_div_r) ! s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: d_v_div_r
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: v_00, v_p1, r_00, r_p1, term1, term2
logical :: dbg
include 'formats'
ierr = 0
dbg = .false.
v_00 = wrap_v_00(s,k)
v_p1 = wrap_v_p1(s,k)
r_00 = wrap_r_00(s,k)
r_p1 = wrap_r_p1(s,k)
if (r_p1%val == 0d0) r_p1 = 1d0
d_v_div_r = v_00/r_00 - v_p1/r_p1 ! units s^-1
! Debugging output to trace values
if (dbg .and. k == -63) then
write(*,*) 'test d_v_div_r, k:', k
write(*,*) 'v_00:', v_00%val, 'v_p1:', v_p1%val
write(*,*) 'r_00:', r_00%val, 'r_p1:', r_p1%val
write(*,*) 'd_v_div_r:', d_v_div_r %val
end if
end function compute_d_v_div_r
function compute_d_v_div_r_opt_time_center(s, k, ierr) result(d_v_div_r) ! s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: d_v_div_r
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: v_00, v_p1, r_00, r_p1
include 'formats'
ierr = 0
v_00 = wrap_opt_time_center_v_00(s,k)
v_p1 = wrap_opt_time_center_v_p1(s,k)
r_00 = wrap_opt_time_center_r_00(s,k)
r_p1 = wrap_opt_time_center_r_p1(s,k)
if (r_p1%val == 0d0) r_p1 = 1d0
d_v_div_r = v_00/r_00 - v_p1/r_p1 ! units s^-1
end function compute_d_v_div_r_opt_time_center
function wrap_Hp_cell(s, k) result(Hp_cell) ! cm
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Hp_cell
Hp_cell = 0.5d0*(wrap_Hp_00(s,k) + wrap_Hp_p1(s,k))
end function wrap_Hp_cell
function Hp_cell_for_Chi(s, k, ierr) result(Hp_cell) ! cm
type (star_info), pointer :: s
integer, intent(in) :: k
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Hp_cell
type(auto_diff_real_star_order1) :: d_00, Peos_00, rmid
real(dp) :: mmid, cgrav_mid
include 'formats'
ierr = 0
Hp_cell = wrap_Hp_cell(s, k)
return
d_00 = wrap_d_00(s, k)
Peos_00 = wrap_Peos_00(s, k)
if (k < s% nz) then
rmid = 0.5d0*(wrap_r_00(s,k) + wrap_r_p1(s,k))
mmid = 0.5d0*(s% m(k) + s% m(k+1))
cgrav_mid = 0.5d0*(s% cgrav(k) + s% cgrav(k+1))
else
rmid = 0.5d0*(wrap_r_00(s,k) + s% r_center)
mmid = 0.5d0*(s% m(k) + s% m_center)
cgrav_mid = s% cgrav(k)
end if
Hp_cell = pow2(rmid)*Peos_00/(d_00*cgrav_mid*mmid)
if (s% alt_scale_height_flag) then
call mesa_error(__FILE__,__LINE__,'Hp_cell_for_Chi: cannot use alt_scale_height_flag')
end if
end function Hp_cell_for_Chi
function compute_Chi_cell(s, k, ierr) result(Chi_cell)
! eddy viscosity energy (Kuhfuss 1986) [erg]
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Chi_cell
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: &
rho2, r6_cell, d_v_div_r, Hp_cell, w_00, d_00, r_00, r_p1
real(dp) :: f, ALFAM_ALFA
logical :: dbg
include 'formats'
ierr = 0
dbg = .false.
! check where we are getting alfam from.
if (s% MLT_option == 'TDC' .and. .not. s% RSP2_flag) then
ALFAM_ALFA = s% alpha_TDC_DampM * s% mixing_length_alpha
else if (s% RSP2_flag) then
ALFAM_ALFA = s% RSP2_alfam * s% mixing_length_alpha
else ! this is for safety, but probably is never called.
ALFAM_ALFA = 0d0
end if
if (ALFAM_ALFA == 0d0 .or. &
k <= s% RSP2_num_outermost_cells_forced_nonturbulent .or. &
k > s% nz - int(s% nz/s% RSP2_nz_div_IBOTOM)) then
Chi_cell = 0d0
if (k >= 1 .and. k <= s% nz) then
s% Chi(k) = 0d0
s% Chi_ad(k) = 0d0
end if
else
Hp_cell = Hp_cell_for_Chi(s, k, ierr)
if (ierr /= 0) return
d_v_div_r = compute_d_v_div_r(s, k, ierr)
if (ierr /= 0) return
! don't need to check if mlt_vc > 0 here.
if (s% MLT_option == 'TDC' .and. .not. s% RSP2_flag) then
if (s% have_mlt_vc .and. s% okay_to_set_mlt_vc) then
w_00 = s% mlt_vc_old(k)/sqrt_2_div_3! same as info%A0 from TDC
else
w_00 = s% mlt_vc(k)/sqrt_2_div_3! same as info%A0 from TDC
end if
else ! normal RSP2
w_00 = wrap_w_00(s,k)
end if
d_00 = wrap_d_00(s,k)
f = (16d0/3d0)*pi*ALFAM_ALFA/s% dm(k)
rho2 = pow2(d_00)
r_00 = wrap_r_00(s,k)
r_p1 = wrap_r_p1(s,k)
r6_cell = 0.5d0*(pow6(r_00) + pow6(r_p1))
Chi_cell = f*rho2*r6_cell*d_v_div_r*Hp_cell*w_00
! units = g^-1 cm s^-1 g^2 cm^-6 cm^6 s^-1 cm
! = g cm^2 s^-2
! = erg
end if
s% Chi(k) = Chi_cell%val
s% Chi_ad(k) = Chi_cell
if (dbg .and. k==-100) then
write(*,*) ' s% ALFAM_ALFA', ALFAM_ALFA
write(*,*) 'Hp_cell', Hp_cell %val
write(*,*) 'd_v_div_r', d_v_div_r %val
write(*,*) ' f', f
write(*,*) 'w_00',w_00 %val
write(*,*) 'd_00 ', d_00 %val
write(*,*) 'rho2 ', rho2 %val
write(*,*) 'r_00', r_00 %val
write(*,*) 'r_p1 ', r_p1 %val
write(*,*) 'r6_cell', r6_cell %val
end if
end function compute_Chi_cell
function compute_Eq_cell(s, k, ierr) result(Eq_cell) ! erg g^-1 s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Eq_cell
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: d_v_div_r, Chi_cell
include 'formats'
ierr = 0
if (s% mixing_length_alpha == 0d0 .or. &
k <= s% RSP2_num_outermost_cells_forced_nonturbulent .or. &
k > s% nz - int(s% nz/s% RSP2_nz_div_IBOTOM)) then
Eq_cell = 0d0
if (k >= 1 .and. k <= s% nz) s% Eq_ad(k) = 0d0
else
Chi_cell = s% Chi_ad(k) ! compute_Chi_cell(s,k,ierr)
if (ierr /= 0) return
d_v_div_r = compute_d_v_div_r_opt_time_center(s, k, ierr)
if (ierr /= 0) return
Eq_cell = 4d0*pi*Chi_cell*d_v_div_r/s% dm(k) ! erg s^-1 g^-1
end if
s% Eq(k) = Eq_cell%val
s% Eq_ad(k) = Eq_cell
end function compute_Eq_cell
function compute_Uq_face(s, k, ierr) result(Uq_face) ! cm s^-2, acceleration
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Uq_face
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Chi_00, Chi_m1, r_00
include 'formats'
ierr = 0
if (s% mixing_length_alpha == 0d0 .or. &
k <= s% RSP2_num_outermost_cells_forced_nonturbulent .or. &
k > s% nz - int(s% nz/s% RSP2_nz_div_IBOTOM)) then
Uq_face = 0d0
else
r_00 = wrap_opt_time_center_r_00(s,k)
! which do we adopt?
Chi_00 = compute_Chi_cell(s,k,ierr) ! s% Chi_ad(k) XXX
!Chi_00 = s% Chi_ad(k) ! compute_Chi_cell(s,k,ierr)
if (k > 1) then
Chi_m1 = shift_m1(compute_Chi_cell(s,k-1,ierr))
!Chi_m1 = shift_m1(s% Chi_ad(k-1)) XXX
if (ierr /= 0) return
else
Chi_m1 = 0d0
end if
Uq_face = 4d0*pi*(Chi_m1 - Chi_00)/(r_00*s% dm_bar(k))
if (k==-56) then
write(*,3) 'RSP2 Uq chi_m1 chi_00 r', k, s% solver_iter, &
Uq_face%val, Chi_m1%val, Chi_00%val, r_00%val
end if
end if
! erg g^-1 cm^-1 = g cm^2 s^-2 g^-1 cm^-1 = cm s^-2, acceleration
s% Uq(k) = Uq_face%val
end function compute_Uq_face
function compute_Source(s, k, ierr) result(Source) ! erg g^-1 s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Source
! source_div_w assumes RSP2_source_seed == 0
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: &
w_00, T_00, d_00, Peos_00, Cp_00, chiT_00, chiRho_00, QQ_00, &
Hp_face_00, Hp_face_p1, PII_face_00, PII_face_p1, PII_div_Hp_cell, &
P_QQ_div_Cp
include 'formats'
ierr = 0
w_00 = wrap_w_00(s, k)
T_00 = wrap_T_00(s, k)
d_00 = wrap_d_00(s, k)
Peos_00 = wrap_Peos_00(s, k)
Cp_00 = wrap_Cp_00(s, k)
chiT_00 = wrap_chiT_00(s, k)
chiRho_00 = wrap_chiRho_00(s, k)
QQ_00 = chiT_00/(d_00*T_00*chiRho_00)
Hp_face_00 = wrap_Hp_00(s,k)
PII_face_00 = s% PII_ad(k) ! compute_PII_face(s, k, ierr)
if (ierr /= 0) return
if (k == s% nz) then
PII_div_Hp_cell = PII_face_00/Hp_face_00
else
Hp_face_p1 = wrap_Hp_p1(s,k)
if (ierr /= 0) return
!PII_face_p1 = shift_p1(compute_PII_face(s, k+1, ierr))
PII_face_p1 = shift_p1(s% PII_ad(k+1))
if (ierr /= 0) return
PII_div_Hp_cell = 0.5d0*(PII_face_00/Hp_face_00 + PII_face_p1/Hp_face_p1)
end if
! Peos_00*QQ_00/Cp_00 = grad_ad if all perfect.
!grad_ad_00 = wrap_grad_ad_00(s, k)
P_QQ_div_Cp = Peos_00*QQ_00/Cp_00 ! use this to be same as RSP
Source = (w_00 + s% RSP2_source_seed)*PII_div_Hp_cell*T_00*P_QQ_div_Cp
! PII units same as Cp = erg g^-1 K^-1
! P*QQ/Cp is unitless (see Y_face)
! Source units = (erg g^-1 K^-1) cm^-1 cm s^-1 K
! = erg g^-1 s^-1
if (k==-109) then
write(*,3) 'RSP2 Source w PII_div_Hp T_P_QQ_div_Cp', k, s% solver_iter, &
Source%val, w_00%val, PII_div_Hp_cell%val, T_00%val*P_QQ_div_Cp% val
!write(*,3) 'RSP2 PII_00 PII_p1 Hp_00 Hp_p1', k, s% solver_iter, &
! PII_face_00%val, PII_face_p1%val, Hp_face_00%val, Hp_face_p1%val
end if
s% SOURCE(k) = Source%val
end function compute_Source
function compute_D(s, k, ierr) result(D) ! erg g^-1 s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: D
type(auto_diff_real_star_order1) :: dw3, w_00
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Hp_cell
include 'formats'
ierr = 0
if (s% mixing_length_alpha == 0d0) then
D = 0d0
else
Hp_cell = wrap_Hp_cell(s,k)
w_00 = wrap_w_00(s,k)
dw3 = pow3(w_00) - pow3(s% RSP2_w_min_for_damping)
D = (s% RSP2_alfad*x_CEDE/s% mixing_length_alpha)*dw3/Hp_cell
! units cm^3 s^-3 cm^-1 = cm^2 s^-3 = erg g^-1 s^-1
end if
if (k==-50) then
write(*,3) 'RSP2 DAMP w Hp_cell dw3', k, s% solver_iter, &
D%val, w_00%val, Hp_cell%val, dw3% val
end if
s% DAMP(k) = D%val
end function compute_D
function compute_Dr(s, k, ierr) result(Dr) ! erg g^-1 s^-1 = cm^2 s^-3
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: Dr
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: &
w_00, T_00, d_00, Cp_00, kap_00, Hp_cell, POM2
real(dp) :: gammar, alpha, POM
include 'formats'
ierr = 0
alpha = s% mixing_length_alpha
gammar = s% RSP2_alfar*x_GAMMAR
if (gammar == 0d0) then
Dr = 0d0
s% DAMPR(k) = 0d0
return
end if
w_00 = wrap_w_00(s,k)
T_00 = wrap_T_00(s,k)
d_00 = wrap_d_00(s,k)
Cp_00 = wrap_Cp_00(s,k)
kap_00 = wrap_kap_00(s,k)
Hp_cell = wrap_Hp_cell(s,k)
POM = 4d0*boltz_sigma*pow2(gammar/alpha) ! erg cm^-2 K^-4 s^-1
POM2 = pow3(T_00)/(pow2(d_00)*Cp_00*kap_00)
! K^3 / ((g cm^-3)^2 (erg g^-1 K^-1) (cm^2 g^-1))
! K^3 / (cm^-4 erg K^-1) = K^4 cm^4 erg^-1
Dr = get_etrb(s,k)*POM*POM2/pow2(Hp_cell)
! (erg cm^-2 K^-4 s^-1) (K^4 cm^4 erg^-1) cm^2 s^-2 cm^-2
! cm^2 s^-3 = erg g^-1 s^-1
s% DAMPR(k) = Dr%val
end function compute_Dr
function compute_C(s, k, ierr) result(C) ! erg g^-1 s^-1
type (star_info), pointer :: s
integer, intent(in) :: k
type(auto_diff_real_star_order1) :: C
integer, intent(out) :: ierr
type(auto_diff_real_star_order1) :: Source, D, Dr
if (s% mixing_length_alpha == 0d0 .or. &
k <= s% RSP2_num_outermost_cells_forced_nonturbulent .or. &
k > s% nz - int(s% nz/s% RSP2_nz_div_IBOTOM)) then
if (k >= 1 .and. k <= s% nz) then
s% SOURCE(k) = 0d0
s% DAMP(k) = 0d0
s% DAMPR(k) = 0d0
s% COUPL(k) = 0d0
s% COUPL_ad(k) = 0d0
end if
C = 0d0
return
end if
Source = compute_Source(s, k, ierr)
if (ierr /= 0) return
D = compute_D(s, k, ierr)
if (ierr /= 0) return
Dr = compute_Dr(s, k, ierr)
if (ierr /= 0) return
C = Source - D - Dr
s% COUPL(k) = C%val
s% COUPL_ad(k) = C
end function compute_C