-
Notifications
You must be signed in to change notification settings - Fork 41
/
Copy pathsolver_support.f90
1560 lines (1373 loc) · 58.6 KB
/
solver_support.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) 2012-2019 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 solver_support
use star_private_def
use utils_lib, only: is_bad
use const_def
use num_def
implicit none
logical, parameter :: dbg = .false.
contains
subroutine set_xscale_info(s, nvar, ierr)
type (star_info), pointer :: s
integer, intent(in) :: nvar
integer, intent(out) :: ierr
integer :: i, j, k, nz, nvar_hydro
real(dp), parameter :: xscale_min = 1
real(dp) :: var_scale, lum_scale, vel_scale, omega_scale
include 'formats'
ierr = 0
if (dbg) write(*, *) 'set_xscale'
nvar_hydro = s% nvar_hydro
nz = s% nz
do k=1,nz
do i=1,nvar
if (i <= nvar_hydro) then ! structure variable
if (i == s% i_j_rot) then
s% x_scale(i,k) = 10d0*sqrt(s% cgrav(k)*s% m(k)*s% r_start(k))
else
s% x_scale(i,k) = max(xscale_min, abs(s% xh_start(i,k)))
end if
else ! abundance variable
s% x_scale(i,k) = max(s% xa_scale, s% xa_start(i-nvar_hydro,k))
end if
end do
end do
contains
subroutine dump_xscale
integer :: k, j, k0, k1
include 'formats'
!write(*,1) 's% xa_scale', s% xa_scale
do k=1,s% nz
do j=1,nvar
write(*,2) 'xscale ' // trim(s% nameofvar(j)), k, s% x_scale(j,k)
end do
write(*,'(A)')
end do
call mesa_error(__FILE__,__LINE__,'set_xscale')
end subroutine dump_xscale
end subroutine set_xscale_info
subroutine eval_equations(s, nvar, ierr)
use hydro_eqns, only: eval_equ
use mix_info, only: set_dxdt_mix
use star_utils, only: update_time, total_times
type (star_info), pointer :: s
integer, intent(in) :: nvar
integer, intent(out) :: ierr
integer :: cnt, i, j, k, nz
integer :: id
real(dp) :: dt, theta_dt
include 'formats'
ierr = 0
nz = s% nz
if (dbg) write(*, *) 'eval_equations'
dt = s% dt
if (s% solver_iter == 0) then
if (dbg) write(*, *) 'skip set_solver_vars on call before 1st iter'
else
if (dbg) write(*, *) 'call set_solver_vars'
call set_solver_vars(s, nvar, dt, ierr)
if (ierr /= 0) then
if (s% report_ierr) &
write(*,2) 'eval_equations: set_solver_vars returned ierr', ierr
return
end if
end if
call set_dxdt_mix(s)
if (ierr == 0) then
do k=1,nz
do j=1,nvar
s% equ(j,k) = 0d0
s% residual_weight(j,k) = 1d0
s% correction_weight(j,k) = 1d0
end do
end do
if (dbg) write(*, *) 'call eval_equ'
call eval_equ(s, nvar, ierr)
if (ierr /= 0) then
if (s% report_ierr) &
write(*, *) 'eval_equations: eval_equ returned ierr', ierr
return
end if
end if
if (ierr /= 0) return
if (.not. s% stop_for_bad_nums) return
cnt = 0
do i=1,nz
do j=1, nvar
if (is_bad_num(s% equ(j, i))) then
cnt = cnt + 1
s% retry_message = 'eval_equations: equ has a bad num'
if (s% report_ierr) then
write(*,4) 'eval_equations: equ has a bad num ' // trim(s% nameofequ(j)), &
j, i, nvar, s% equ(j, i)
write(*,2) 'cell', i
write(*,2) 'nz', s% nz
end if
if (s% stop_for_bad_nums) call mesa_error(__FILE__,__LINE__,'eval_equations')
end if
end do
end do
if (cnt > 0) then
ierr = -1
return
end if
contains
subroutine dump_eval_equ
integer :: k, j, k0, k1
include 'formats'
do k=1,s% nz
do j=1,nvar
write(*,2) 'dx ' // trim(s% nameofvar(j)), k, s% solver_dx(j, k)
end do
write(*,'(A)')
end do
call mesa_error(__FILE__,__LINE__,'dump_eval_equ')
end subroutine dump_eval_equ
end subroutine eval_equations
subroutine sizequ(s, nvar, equ_norm, equ_max, k_max, j_max, ierr) ! equ = residuals
type (star_info), pointer :: s
integer, intent(in) :: nvar
real(dp), intent(out) :: equ_norm, equ_max
integer, intent(out) :: k_max, j_max, ierr
integer :: j, k, num_terms, n, nz, nvar_hydro, nvar_chem, &
max_loc, skip_eqn1, skip_eqn2, skip_eqn3
real(dp) :: sumequ, absq, max_energy_resid, avg_energy_resid
logical :: dbg
include 'formats'
ierr = 0
equ_norm = 0d0
equ_max = 0d0
k_max = 0
j_max = 0
dbg = s% solver_check_everything
nvar_hydro = min(nvar, s% nvar_hydro)
nvar_chem = s% nvar_chem
nz = s% nz
n = nz
num_terms = 0
sumequ = 0
skip_eqn1 = 0
skip_eqn2 = 0
skip_eqn3 = 0
if (s% convergence_ignore_equL_residuals) skip_eqn1 = s% i_equL
if (s% convergence_ignore_alpha_RTI_residuals) skip_eqn2 = s% i_dalpha_RTI_dt
if (s% do_burn .or. s% do_mix) then
num_terms = num_terms + nvar*nz
if (skip_eqn1 > 0) num_terms = num_terms - nz
if (skip_eqn2 > 0) num_terms = num_terms - nz
if (skip_eqn3 > 0) num_terms = num_terms - nz
do k = 1, nz
do j = 1, nvar
if (j == skip_eqn1 .or. j == skip_eqn2 .or. j == skip_eqn3) cycle
if (is_bad(s% equ(j,k)) .or. is_bad(s% residual_weight(j,k))) then
ierr = 1
return
end if
absq = abs(s% equ(j,k)*s% residual_weight(j,k))
sumequ = sumequ + absq
if (absq > equ_max) then
equ_max = absq
j_max = j
k_max = k
end if
end do
end do
else
if (skip_eqn1 == 0 .and. skip_eqn2 == 0) then
num_terms = num_terms + nvar_hydro*nz
else if (skip_eqn1 > 0 .and. skip_eqn2 > 0) then
num_terms = num_terms + (nvar_hydro-2)*nz
else
num_terms = num_terms + (nvar_hydro-1)*nz
end if
do k = 1, nz
do j = 1, nvar_hydro
if (j == skip_eqn1 .or. j == skip_eqn2) cycle
absq = abs(s% equ(j,k)*s% residual_weight(j,k))
sumequ = sumequ + absq
if (is_bad(sumequ)) then
if (dbg) then
write(*,3) trim(s% nameofequ(j)) // ' sumequ', j, k, sumequ
call mesa_error(__FILE__,__LINE__,'sizeq 1')
end if
ierr = -1
if (s% report_ierr) &
write(*,3) 'bad equ(j,k)*s% residual_weight(j,k) ' // trim(s% nameofequ(j)), &
j, k, s% equ(j,k)*s% residual_weight(j,k)
if (s% stop_for_bad_nums) call mesa_error(__FILE__,__LINE__,'sizeq 2')
return
end if
if (absq > equ_max) then
equ_max = absq
j_max = j
k_max = k
end if
end do
end do
end if
if (s% do_burn .or. s% do_mix) then
num_terms = num_terms + nvar_chem*nz
do k = 1, nz
do j = nvar_hydro+1, nvar
absq = abs(s% equ(j,k)*s% residual_weight(j,k))
sumequ = sumequ + absq
if (absq > equ_max) then
equ_max = absq
j_max = j
k_max = k
end if
end do
end do
end if
equ_norm = sumequ/num_terms
if (dbg) write(*,4) trim(s% nameofequ(j_max)) // ' sizequ equ_max norm', &
k_max, s% solver_iter, s% model_number, equ_max, equ_norm
if (dbg) call dump_equ
return
call dump_equ
call mesa_error(__FILE__,__LINE__,'sizequ')
contains
subroutine dump_equ
integer :: k, j, k0, k1
include 'formats'
do k=1,s% nz
do j=1,nvar
write(*,3) 'equ ' // trim(s% nameofequ(j)), &
k, s% solver_iter, s% equ(j, k)
end do
write(*,'(A)')
!if (k == 6) exit
end do
end subroutine dump_equ
end subroutine sizequ
subroutine sizeB(s, nvar, B, max_correction, correction_norm, max_zone, max_var, ierr)
type (star_info), pointer :: s
integer, intent(in) :: nvar
real(dp), pointer, dimension(:,:) :: B ! (nvar, nz)
real(dp), intent(out) :: correction_norm ! a measure of the average correction
real(dp), intent(out) :: max_correction ! magnitude of the max correction
integer, intent(out) :: max_zone, max_var, ierr
integer :: k, i, nz, num_terms, j, n, nvar_hydro, jmax, num_xa_terms, &
skip1, skip2, skip3, skip4, skip5
real(dp) :: abs_corr, sum_corr, sum_xa_corr, x_limit, &
max_abs_correction, max_abs_correction_cv, max_abs_corr_for_k, max_abs_xa_corr_for_k
logical :: found_NaN, found_bad_num, report
real(dp), parameter :: frac = 0.1d0
logical, parameter :: dbg = .false.
logical, parameter :: check_for_bad_nums = .true.
logical, parameter :: save_max_abs_corr_for_k = .true.
include 'formats'
if (dbg) write(*, *) 'enter sizeB'
ierr = 0
nz = s% nz
n = nz
nvar_hydro = min(nvar, s% nvar_hydro)
if (s% include_L_in_correction_limits) then
skip1 = 0
do k=1,nz
s% correction_weight(s% i_lum,k) = 1d0/(frac*s% L_start(1) + abs(s% L(k)))
end do
else
skip1 = s% i_lum
end if
if (s% u_flag .and. s% include_u_in_correction_limits) then
skip2 = 0
do k=1,nz
s% correction_weight(s% i_u,k) = 1d0/(frac*s% csound_start(k) + abs(s% u(k)))
end do
else if (s% v_flag .and. s% include_v_in_correction_limits) then
skip2 = 0
do k=1,nz
s% correction_weight(s% i_v,k) = 1d0/(frac*s% csound_start(k) + abs(s% v(k)))
end do
else if (s% u_flag) then
skip2 = s% i_u
else
skip2 = s% i_v
end if
if (s% RSP2_flag .and. s% include_w_in_correction_limits) then
skip3 = 0
do k=1,nz
if (abs(s% w(k)) < 1d0) then
s% correction_weight(s% i_w,k) = 0d0
else
s% correction_weight(s% i_w,k) = 1d0/(1d3 + abs(s% w(k))) ! don't sweat the small stuff
end if
end do
else
skip3 = s% i_w
end if
skip4 = s% i_Hp
skip5 = 0
max_zone = 0
max_var = 0
num_terms = 0
num_xa_terms = 0
sum_corr = 0
sum_xa_corr = 0
max_correction = 0
max_abs_correction = 0
x_limit = s% correction_xa_limit
found_NaN = .false.
found_bad_num = .false.
report = s% report_ierr
cell_loop: do k = 1, nz
max_abs_corr_for_k = 0
max_abs_xa_corr_for_k = 0
if (s% do_burn .or. s% do_mix) then
jmax = nvar
else
jmax = nvar_hydro
end if
var_loop: do j = 1, jmax
if (j == skip1 .or. &
j == skip2 .or. &
j == skip3 .or. &
j == skip4 .or. &
j == skip5 .or. &
j == s% i_alpha_RTI) cycle
if (check_for_bad_nums) then
if (is_bad_num(B(j,k)*s% correction_weight(j,k))) then
found_bad_num = .true.
if (report) write(*,2) 'sizeB: bad num for correction ' // &
s% nameofvar(j), k, B(j,k)*s% correction_weight(j,k)
if (s% stop_for_bad_nums) then
found_NaN = .true.
write(*,3) s% nameofvar(j) // ' B(j,k)*s% correction_weight(j,k)', &
j, k, B(j,k)*s% correction_weight(j,k)
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
max_zone = k
max_var = j
exit cell_loop
cycle
end if
end if
if (j > nvar_hydro) then
if (s% xa_start(j-nvar_hydro,k) < x_limit) cycle
end if
abs_corr = abs(B(j,k)*s% correction_weight(j,k))
if (is_bad_num(abs_corr)) then
found_bad_num = .true.
if (report) write(*,3) 'B(j,k)*s% correction_weight(j,k)', &
j, k, B(j,k)*s% correction_weight(j,k)
if (s% stop_for_bad_nums) found_NaN = .true.
end if
if (abs_corr > max_abs_corr_for_k &
.and. .not. (j > nvar_hydro .and. s% ignore_species_in_max_correction)) &
max_abs_corr_for_k = abs_corr
if (abs_corr > max_abs_correction &
.and. .not. (j > nvar_hydro .and. s% ignore_species_in_max_correction)) then
max_correction = B(j,k)*s% correction_weight(j,k)
max_abs_correction = abs_corr
max_zone = k
max_var = j
end if
if (j > nvar_hydro) then
num_xa_terms = num_xa_terms + 1
sum_xa_corr = sum_xa_corr + abs_corr
if (abs_corr > max_abs_xa_corr_for_k) &
max_abs_xa_corr_for_k = abs_corr
else
num_terms = num_terms + 1
sum_corr = sum_corr + abs_corr
end if
end do var_loop
if (num_xa_terms > 0) then
num_terms = num_terms + 1
sum_corr = sum_corr + sum_xa_corr/num_xa_terms
end if
if (s% do_burn .or. s% do_mix) then
species_loop: do j = nvar_hydro+1, nvar
i = j - s% nvar_hydro
if (check_for_bad_nums) then
if (is_bad_num(B(j,k)*s% correction_weight(j,k))) then
found_bad_num = .true.
if (report) write(*,3) 'chem B(j,k)*s% correction_weight(j,k)', j, k, &
B(j,k)*s% correction_weight(j,k)
if (s% stop_for_bad_nums) then
found_NaN = .true.
write(*,3) 'chem B(j,k)*s% correction_weight(j,k)', &
j, k, B(j,k)*s% correction_weight(j,k)
call mesa_error(__FILE__,__LINE__,'sizeB')
max_zone = k
max_var = j
exit cell_loop
end if
end if
end if
! recall that correction dx = B*xscale, so B is a relative correction
if (s% xa_start(i,k) >= x_limit) then
abs_corr = abs(B(j,k)*s% correction_weight(j,k))
if (.not. s% ignore_species_in_max_correction) then
if (abs_corr > max_abs_corr_for_k) max_abs_corr_for_k = abs_corr
if (abs_corr > max_abs_correction) then
max_abs_correction = abs_corr
max_correction = B(j,k)*s% correction_weight(j,k)
max_zone = k
max_var = j
end if
end if
sum_corr = sum_corr + abs_corr
num_terms = num_terms + 1
end if
end do species_loop
end if
s% max_abs_xa_corr(k) = max_abs_xa_corr_for_k
end do cell_loop
if (found_bad_num) then
ierr = -1
if (found_NaN .and. s% stop_for_bad_nums) then
write(*,*) 'found bad num'
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
if (.not. dbg) return
end if
if (is_bad_num(sum_corr)) then
ierr = -1
if (s% stop_for_bad_nums) then
if (report) write(*,*) 'sum_corr', sum_corr
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
if (.not. dbg) return
write(*,*) 'sum_corr', sum_corr
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
correction_norm = sum_corr/num_terms !sqrt(sum_corr/num_terms)
if (dbg) then
write(*,2) 'sizeB: iter, correction_norm, max_correction', &
s% solver_iter, correction_norm, max_correction
if (max_correction > 1d50 .or. is_bad_num(correction_norm)) then
call show_stuff
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
end if
if (s% solver_show_correction_info) call show_stuff
abs_corr = max_abs_correction
s% abs_max_corr2 = s% abs_max_corr1; s% abs_max_corr1 = abs_corr
s% max_var2 = s% max_var1; s% max_var1 = max_var
s% max_zone2 = s% max_zone1; s% max_zone1 = max_zone
if (ierr /= 0) call mesa_error(__FILE__,__LINE__,'ierr in sizeB')
if (is_bad_num(max_correction)) then
ierr = -1
if (s% stop_for_bad_nums) then
if (report) write(*,*) 'max_correction', max_correction
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
if (.not. dbg) return
write(*,*) 'max_correction', max_correction
call mesa_error(__FILE__,__LINE__,'sizeB')
end if
if (s% solver_iter < 3) return
! check for flailing
if ( &
abs_corr > s% tol_max_correction .and. &
abs_corr > s% abs_max_corr1 .and. s% abs_max_corr1 > s% abs_max_corr2 .and. &
max_zone == s% max_zone1 .and. s% max_zone1 == s% max_zone2 .and. &
max_var == s% max_var1 .and. s% max_var1 == s% max_var2) then
if (s% solver_show_correction_info) then
write(*,*) 'give up because diverging'
end if
max_correction = 1d99
end if
contains
subroutine show_stuff
integer :: j, k
real(dp) :: dx, prev, new
include 'formats'
if (s% solver_iter == 1) then
write(*,'(A)')
write(*,'(4a7,12a16,99a13)') &
'model', 'iter', 'var', 'zone', &
'corr norm', 'max corr', 'xscale', &
'dx', 'new-prev', 'new', 'prev', &
'mass loc', 'log dt/yr', 'lgE', 'lgT', 'lgRho'
end if
k = max_zone
j = max_var
if (j > nvar_hydro) then
prev = s% xa_start(j - nvar_hydro,k)
else
prev = s% xh_start(j,k)
end if
dx = B(j,k)*s% correction_weight(j,k)*s% x_scale(j,k)
new = prev + dx
write(*,'(2i7,a7,i7,12e16.8,99f13.8)') &
s% model_number, s% solver_iter, trim(s% nameofvar(max_var)), k, &
correction_norm, B(j,k)*s% correction_weight(j,k), s% x_scale(j,k), &
dx, new - prev, new, prev, &
s% m(k)/Msun, log10(s% dt/secyer), &
s% lnE(k)/ln10, s% lnT(k)/ln10, &
s% lnd(k)/ln10
end subroutine show_stuff
subroutine dump_B
integer :: k, j, k0, k1
include 'formats'
do k=1,s% nz
do j=1,nvar
write(*,2) 'B ' // trim(s% nameofequ(j)), k, B(j, k)
end do
write(*,'(A)')
end do
call mesa_error(__FILE__,__LINE__,'dump_equ')
end subroutine dump_B
end subroutine sizeB
! the proposed change to dx is B*xscale*correction_factor
! edit correction_factor and/or B as necessary so that the new dx will be valid.
! set ierr nonzero if things are beyond repair.
subroutine Bdomain(s, nvar, B, correction_factor, ierr)
use const_def, only: dp
use chem_def, only: chem_isos
use star_utils, only: current_min_xa_hard_limit, rand
use rsp_def, only: EFL0
type (star_info), pointer :: s
integer, intent(in) :: nvar
real(dp), pointer, dimension(:,:) :: B ! (nvar, nz)
real(dp), intent(inout) :: correction_factor
integer, intent(out) :: ierr
integer :: id, i, j, k, nz, species, bad_j, bad_k, &
i_alpha_RTI, i_w_div_wc
real(dp) :: alpha, min_alpha, new_xa, old_xa, dxa, eps, min_xa_hard_limit, &
old_E, dE, new_E, old_lnd, dlnd, new_lnd, dw, new_w, &
dw_div_wc, old_w_div_wc, new_w_div_wc, dconv_vel, old_conv_vel, new_conv_vel, &
dalpha_RTI, new_alpha_RTI, old_alpha_RTI, log_conv_vel_v0, &
dlum_surf, old_lum_surf, new_lum_surf
include 'formats'
ierr = 0
min_alpha = 1d0
nz = s% nz
if (s% RSP2_flag) & ! clip change in w to maintain non-negativity.
call clip_so_non_negative(s% i_w, 0d0)
if (s% RTI_flag) & ! clip change in alpha_RTI to maintain non-negativity.
call clip_so_non_negative(s% i_alpha_RTI, 0d0)
if (s% i_lum>=0 .and. s% scale_max_correction_for_negative_surf_lum) then
!ensure surface luminosity does not become negative
dlum_surf = B(s% i_lum,1)*s% x_scale(s% i_lum,1)
old_lum_surf = s% xh_start(s% i_lum,1) + s% solver_dx(s% i_lum,1)
new_lum_surf = old_lum_surf + dlum_surf
if (new_lum_surf < 0d0 .and. old_lum_surf > 0d0) then
correction_factor = min(correction_factor, &
-s% max_frac_for_negative_surf_lum*old_lum_surf/dlum_surf)
end if
end if
!if (s% w_div_wc_flag) then
! do k=1,nz
! old_xa = s% xh_start(s% i_w_div_wc,k) + dx(s% i_w_div_wc,k)
! dxa = B(s% i_w_div_wc,k)*xscale(s% i_w_div_wc,k)*correction_factor
! if (dxa > 0.9d0*(1-old_xa)) then
! if(k==91)write(*,*) "problems", old_xa, dxa, old_xa+dxa
! correction_factor = min(correction_factor,(0.9d0*(1-old_xa))/B(s% i_w_div_wc,k)*xscale(s%i_w_div_wc,k))
! dxa = B(s% i_w_div_wc,k)*xscale(s% i_w_div_wc,k)*correction_factor
! if(k==91)write(*,*) "need to reduce correction", k, old_xa+dxa
! end if
! end do
!end if
if ((.not. s% do_solver_damping_for_neg_xa) .or. &
(.not. (s% do_mix .or. s% do_burn))) then
if (min_alpha < 1d0) then
min_alpha = max(min_alpha, s% corr_coeff_limit)
correction_factor = min_alpha*correction_factor
end if
return
end if
bad_j = 0
bad_k = 0
if (nvar > s% nvar_hydro) then
species = s% species
min_xa_hard_limit = current_min_xa_hard_limit(s)
eps = -0.5d0*min_xa_hard_limit ! allow xa to be this much below 0d0
do k=1,nz
do j=1,species
i = j + s% nvar_hydro
old_xa = s% xa_start(j,k) + s% solver_dx(i,k)
if (old_xa <= 1d-90) cycle
dxa = B(i,k)*s% x_scale(i,k)*correction_factor
new_xa = old_xa + dxa
if (new_xa >= 0d0) cycle
alpha = -(old_xa + eps)/dxa
! so dxa*alpha = -old_xa - eps,
! and therefore old_xa + alpha*dxa = -eps = 0.5*min_xa_hard_limit
if (alpha < min_alpha) then
min_alpha = alpha
bad_j = j
bad_k = k
end if
end do
end do
end if
min_alpha = max(min_alpha, s% corr_coeff_limit)
correction_factor = min_alpha*correction_factor
if (s% trace_solver_damping .and. min_alpha < 1d0 .and. bad_j > 0) then
write(*,4) 'solver damping to avoid negative mass fractions: ' // &
trim(chem_isos% name(s% chem_id(bad_j))), bad_k, &
s% model_number, s% solver_iter, min_alpha
end if
contains
subroutine clip_so_non_negative(i,minval)
integer, intent(in) :: i
real(dp), intent(in) :: minval
real(dp) :: dval, old_val, new_val
do k = 1, s% nz
dval = B(i,k)*s% x_scale(i,k)*correction_factor
old_val = s% xh_start(i,k) + s% solver_dx(i,k)
new_val = old_val + dval
if (dval >= 0) cycle
if (new_val >= 0d0) cycle
dval = minval - old_val
B(i,k) = dval/(s% x_scale(i,k)*correction_factor)
end do
end subroutine clip_so_non_negative
end subroutine Bdomain
subroutine inspectB(s, nvar, B, ierr)
type (star_info), pointer :: s
integer, intent(in) :: nvar
real(dp), pointer, dimension(:,:) :: B ! (nvar, nz)
integer, intent(out) :: ierr
integer :: id
integer, parameter :: inspectB_iter_stop = -1
include 'formats'
if (dbg) write(*, *) 'inspectB', s% solver_iter
ierr = 0
if (s% solver_iter == inspectB_iter_stop) then
call dumpB
call mesa_error(__FILE__,__LINE__,'debug: inspectB')
end if
contains
subroutine dumpB
integer :: k, j, k0, k1
include 'formats'
do k=1,s% nz
do j=1,nvar
write(*,2) 'B ' // trim(s% nameofvar(j)), k, B(j, k)
write(*,2) 'xscale ' // trim(s% nameofvar(j)), k, s% x_scale(j, k)
write(*,2) 'dx ' // trim(s% nameofvar(j)), k, s% solver_dx(j, k)
end do
write(*,'(A)')
end do
call mesa_error(__FILE__,__LINE__,'dumpB')
end subroutine dumpB
end subroutine inspectB
! about to declare victory... but may want to do another iteration
! 1 means force another iteration
! 0 means don't need to force another
! -1 means failure. solver returns with non-convergence.
integer function force_another_iteration(s, iter, itermin)
type (star_info), pointer :: s
integer, intent(in) :: iter ! have finished this many iterations and have converged
integer, intent(in) :: itermin ! this is the requested minimum. iter may be < itermin.
integer :: k, res
include 'formats'
if (iter < itermin) then
force_another_iteration = 1
return
end if
force_another_iteration = 0
if (s% k_below_just_added > 1 .and. &
s% num_surf_revisions < s% max_num_surf_revisions .and. &
abs(s% lnS(1) - s% surf_lnS) > &
s% max_abs_rel_change_surf_lnS*max(s% lnS(1),s% surf_lnS)) then
s% surf_lnS = s% lnS(1)
s% num_surf_revisions = s% num_surf_revisions + 1
force_another_iteration = 1
s% used_extra_iter_in_solver_for_accretion = .true.
return
end if
end function force_another_iteration
subroutine set_solver_vars(s, nvar, dt, ierr)
type (star_info), pointer :: s
integer, intent(in) :: nvar
real(dp), intent(in) :: dt
integer, intent(out) :: ierr
s% num_solver_setvars = s% num_solver_setvars + 1
call set_vars_for_solver(s, nvar, 1, s% nz, dt, ierr)
end subroutine set_solver_vars
subroutine set_vars_for_solver(s, nvar, nzlo, nzhi, dt, ierr)
use const_def, only: secyer, Msun, Lsun, Rsun
use star_utils, only: set_rmid, set_dm_bar, set_m_and_dm, set_rv_info
use star_utils, only: current_min_xa_hard_limit, current_sum_xa_hard_limit, &
lookup_nameofvar
use hydro_vars, only: set_hydro_vars
use hydro_rotation, only: set_i_rot, set_omega
use mix_info, only: get_convection_sigmas
use chem_def
type (star_info), pointer :: s
integer, intent(in) :: nvar, nzlo, nzhi
real(dp), intent(in) :: dt
integer, intent(out) :: ierr
logical, parameter :: &
skip_basic_vars = .true., &
skip_micro_vars = .false., &
skip_kap = .false., &
skip_neu = .false., &
skip_net = .false., &
skip_eos = .false., &
skip_mlt = .false., &
skip_grads = .true., &
skip_rotation = .true., &
skip_m_grav_and_grav = .true., &
skip_brunt = .true., &
skip_mixing_info = .true., &
skip_set_cz_bdy_mass = .true., &
skip_other_cgrav = .true.
logical :: do_chem, try_again, do_edit_lnR, report_dx
integer :: i, j, k, kk, klo, khi, i_var, &
i_lnd, i_lnT, i_lnR, i_lum, i_w, i_Hp, i_v, &
i_u, i_alpha_RTI, i_w_div_wc, i_j_rot, &
fe56, nvar_chem, species, nz, nvar_hydro
real(dp), dimension(:, :), pointer :: xh_start, xa_start
integer :: op_err, kbad, &
cnt, max_fixes, loc(2), k_lo, k_hi, k_const_mass
real(dp) :: r2, xavg, du, u00, um1, dx_for_i_var, x_for_i_var, &
dq_sum, xa_err_norm, d_dxdt_dx, min_xa_hard_limit, sum_xa_hard_limit
logical :: do_lnd, do_lnT, do_lnR, do_lum, do_w, &
do_u, do_v, do_alpha_RTI, do_w_div_wc, do_j_rot
include 'formats'
ierr = 0
nz = s% nz
nvar_chem = s% nvar_chem
species = s% species
nvar_hydro = s% nvar_hydro
d_dxdt_dx = 1d0/s% dt
xh_start => s% xh_start
xa_start => s% xa_start
report_dx = &
s% solver_test_partials_dx_0 > 0d0 .and. &
s% solver_test_partials_k > 0 .and. &
s% solver_call_number == s% solver_test_partials_call_number .and. &
s% solver_test_partials_iter_number == s% solver_iter .and. &
len_trim(s% solver_test_partials_show_dx_var_name) > 0
if (report_dx) then
k = s% solver_test_partials_k
i_var = lookup_nameofvar(s, s% solver_test_partials_show_dx_var_name)
if (i_var > 0) then
if (i_var > nvar_hydro) then
dx_for_i_var = s% solver_dx(i_var,k)
x_for_i_var = xa_start(i_var-nvar_hydro,k) + s% solver_dx(i_var,k)
else
dx_for_i_var = s% solver_dx(i_var,k)
x_for_i_var = xh_start(i_var,k) + s% solver_dx(i_var,k)
end if
write(*,3) 'dx, x for var name ' // &
trim(s% solver_test_partials_show_dx_var_name), &
k, s% solver_iter, dx_for_i_var, x_for_i_var
end if
end if
do_chem = (s% do_burn .or. s% do_mix)
if (dbg .and. .not. skip_mixing_info) write(*,2) 'redo mix info iter', s% solver_iter
min_xa_hard_limit = current_min_xa_hard_limit(s)
sum_xa_hard_limit = current_sum_xa_hard_limit(s)
i_lnd = s% i_lnd
i_lnT = s% i_lnT
i_lnR = s% i_lnR
i_lum = s% i_lum
i_w = s% i_w
i_Hp = s% i_Hp
i_v = s% i_v
i_u = s% i_u
i_alpha_RTI = s% i_alpha_RTI
i_w_div_wc = s% i_w_div_wc
i_j_rot = s% i_j_rot
do_lnd = i_lnd > 0 .and. i_lnd <= nvar
do_lnT = i_lnT > 0 .and. i_lnT <= nvar
do_lnR = i_lnR > 0 .and. i_lnR <= nvar
do_lum = i_lum > 0 .and. i_lum <= nvar
do_w = i_w > 0 .and. i_w <= nvar
do_v = i_v > 0 .and. i_v <= nvar
do_u = i_u > 0 .and. i_u <= nvar
do_alpha_RTI = i_alpha_RTI > 0 .and. i_alpha_RTI <= nvar
do_w_div_wc = i_w_div_wc > 0 .and. i_w_div_wc <= nvar
do_j_rot = i_j_rot > 0 .and. i_j_rot <= nvar
fe56 = s% net_iso(ife56)
if (fe56 /= 0) fe56 = nvar_hydro+fe56
if (nvar > nvar_hydro) then
do k=1,nz
do j=1,species
s% xa_sub_xa_start(j,k) = s% solver_dx(j+nvar_hydro,k)
s% xa(j,k) = xa_start(j,k) + s% solver_dx(j+nvar_hydro,k)
end do
end do
max_fixes = 0 !5
do cnt=1,max_fixes
loc = minloc(s% xa(1:species,1:nz))
j = loc(1)
k = loc(2)
if (s% xa(j,k) >= 1d-3*min_xa_hard_limit) exit ! too good to fix
if (s% xa(j,k) < min_xa_hard_limit) then
if (s% report_ierr) then
khi = nz
do kk=k+1,nz
if (s% xa(j,kk) < min_xa_hard_limit) cycle
khi = kk-1; exit
end do
klo = 1
do kk=k-1,1,-1
if (s% xa(j,kk) < min_xa_hard_limit) cycle
klo = kk+1; exit
end do
do k=klo,khi
write(*,2) &
'negative ' // trim(chem_isos% name(s% chem_id(j))), &
k, s% xa(j,k), xa_start(j,k), s% solver_dx(nvar_hydro+j,k), s% m(k)/Msun
end do
end if
s% retry_message = 'some abundance < min_xa_hard_limit'
ierr = -1
return
exit ! too bad to fix
end if
if (k == 1) then
k_lo = 1; k_hi = 2
else if (k == nz) then
k_lo = nz-1; k_hi = nz
else if (s% sig(k) > s% sig(k+1)) then
k_lo = k-1; k_hi = k
else if (s% sig(k+1) > 0) then
k_lo = k; k_hi = k+1
else
exit
end if
try_again = .true.
do while (try_again .and. sum(s% xa(j,k_lo:k_hi)*s% dq(k_lo:k_hi)) < 0)
try_again = .false.
if (k_lo > 1) then
if (s% sig(k_lo) > 0) then
k_lo = k_lo-1
try_again = .true.
end if
end if
if (k_hi < nz) then
if (s% sig(k_hi+1) > 0) then
k_hi = k_hi+1
try_again = .true.
end if
end if
end do
!write(*,3) 'no extend', k_lo, k_hi
if (.not. try_again) exit
dq_sum = sum(s% dq(k_lo:k_hi))
if (s% report_ierr) then
write(*,5) 'fix xa(j,k_lo:k_hi)', j, k_lo, k, k_hi, s% xa(j,k), &
sum(s% xa(j,k_lo:k_hi)*s% dq(k_lo:k_hi))/dq_sum, dq_sum
!stop