-
Notifications
You must be signed in to change notification settings - Fork 44
/
Copy pathconv_premix.f90
1630 lines (1077 loc) · 41.2 KB
/
conv_premix.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) 2017-2018 Rich Townsend & 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 conv_premix
! Uses
use const_def
use star_private_def
use chem_def
use chem_lib
use num_lib
! No implicit typing
implicit none
! Parameter definitions
integer, parameter :: BURN_TYPE_NONE = 1
integer, parameter :: BURN_TYPE_H = 2
integer, parameter :: BURN_TYPE_HE = 3
integer, parameter :: BURN_TYPE_Z = 4
integer, parameter :: FIXED_PT_MODE = 5
integer, parameter :: FIXED_DT_MODE = 6
logical, parameter :: TRACE_ALL = .FALSE.
! Derived-type definitions
! The zone_info type stores information about the extent &
! properties of each convective zone
type zone_info
integer :: kc_t = 0 ! Cell index of top boundary
integer :: kc_b = 0 ! Cell index of bot boundary
real(dp) :: vc_t = 0._dp ! Convective velocity near top boundary
real(dp) :: vc_b = 0._dp ! Convective velocity near bot boundary
real(dp) :: vp_t = 0._dp ! Penetration velocity at top boundary
real(dp) :: vp_b = 0._dp ! Penetration velocity at bot boundary
real(dp) :: dt_t = 0._dp ! Clock for top boundary
real(dp) :: dt_b = 0._dp ! Clock for bottom boundary
integer :: burn_type = BURN_TYPE_NONE ! Type of burning in zone
logical :: sel_t = .FALSE. ! Flag to select top boundary for modification
logical :: sel_b = .FALSE. ! Flag to select bot boundary for modification
logical :: split_retreat = .FALSE. ! Flag to indicate the zone has split, or the trailing bdy retreated
real(dp), allocatable :: avg_xa(:) ! Initial average abundances
real(dp), allocatable :: davg_xa_dt(:) ! Initial rate of change of average abundances (due to burning)
integer :: avoid_inc_iso = 0 ! Isotope id for which increases should be avoided
end type zone_info
! The saved_data type saves cell and face data from star_info, to
! enable us to restore the latter to an earlier state. Only data in
! cells kc_t:kc_b are used, but the arrays are full-size
type saved_data
real(dp), allocatable :: xa(:,:)
real(dp), allocatable :: lnd(:)
real(dp), allocatable :: rho(:)
real(dp), allocatable :: lnPgas(:)
real(dp), allocatable :: Pgas(:)
real(dp), allocatable :: gradL_composition_term(:)
integer, allocatable :: update_mode(:)
integer :: kc_t = 0
integer :: kc_b = 0
end type saved_data
! Access specifers
private
public :: do_conv_premix
! Procedures
contains
subroutine do_conv_premix (s, ierr)
use star_utils, only: start_time, update_time
type(star_info), pointer :: s
integer, intent(out) :: ierr
logical, parameter :: TRACE_CONV_PREMIX = TRACE_ALL
integer :: update_mode(s%nz)
type(saved_data) :: sd
type(zone_info), allocatable :: zi(:)
integer :: j_it
real(dp), allocatable :: dr(:)
real(dp), allocatable :: dt_t(:)
real(dp), allocatable :: dt_b(:)
integer :: i_bdy
logical :: t_bdy
integer(8) :: time0
real(dp) :: total
ierr = 0
if (s% doing_timing) call start_time(s, time0, total)
! Initialize the update_mode values. These control how cells will
! be updated (i.e., how the thermodynamic state will be
! re-evaluated after abundances have changed). For untouched
! cells, update_mode is determined by the lnPgas_flag setting
!if (s%lnPgas_flag) then
! update_mode = FIXED_PT_MODE
!else
update_mode = FIXED_DT_MODE
!endif
! Allocate the saved data arrays
call alloc_saved_data_(s, sd)
! Initialize the zone info table
call init_zone_info_(s, zi)
! Loop until there are no boundaries left to advance
j_it = 0
iter_loop : do
j_it = j_it + 1
! Validate the current zone info table
call validate_zone_info_(s, zi)
! Check for completion
if (.NOT. (ANY(zi%sel_t .OR. zi%sel_b))) exit iter_loop
! Calculate mixing timescales for each boundary (the explicit
! allocations of dt_t and dt_b are required to workaround a
! gfortran bug, where allocate-on-assign will raise bogus array
! bounds violations)
dr = s%r(zi%kc_t+1) - s%r(zi%kc_b)
allocate(dt_t(SIZE(zi)))
allocate(dt_b(SIZE(zi)))
dt_t = MERGE(dr/zi%vc_t, HUGE(0._dp), MASK=zi%sel_t)
dt_b = MERGE(dr/zi%vc_b, HUGE(0._dp), MASK=zi%sel_b)
! Select the boundary with the shortest timescale
i_bdy = MINLOC(MIN(dt_t, dt_b), DIM=1)
t_bdy = dt_t(i_bdy) < dt_b(i_bdy)
if (TRACE_CONV_PREMIX) then
write(*,*) 'Selected i_bdy, t_bdy:', i_bdy, t_bdy
end if
! Advance the selected boundary, modifying the model and the
! zone info table
call advance_bdy_(s, update_mode, zi, i_bdy, t_bdy, sd, j_it)
! Loop around
deallocate(dt_t)
deallocate(dt_b)
end do iter_loop
! Finish
s% need_to_setvars = .true.
if (s% doing_timing) call update_time(s, time0, total, s% time_conv_premix)
return
end subroutine do_conv_premix
!****
subroutine advance_bdy_ (s, update_mode, zi, i_bdy, t_bdy, sd, j_it)
type(star_info), pointer :: s
integer, intent(inout) :: update_mode(:)
type(zone_info), allocatable, intent(inout) :: zi(:)
integer, intent(inout) :: i_bdy
logical, intent(in) :: t_bdy
type(saved_data), intent(inout) :: sd
integer, intent(in) :: j_it
logical, parameter :: TRACE_ADVANCE_BDY = TRACE_ALL
integer :: j_ad
character(256) :: filename
! Advance the top (bottom) boundary of the zone zi(i_bdy). t_bdy
! is .true. for the top boundary, .false. for the bottom. Because
! the zone info table can potentially change during the process
! (as zones split/merge), so can i_bdy
if (TRACE_ADVANCE_BDY) then
write(*,*) 'Advancing zone boundary: i_bdy, t_bdy, j_it=', i_bdy, t_bdy, j_it
end if
if (s% conv_premix_dump_snapshots) then
write(filename, 100) 'cpm-snap.', s%model_number, '.it_', j_it, '.ad_', 0, '.dat'
100 format(A,I5.5,A,I5.5,A,I5.5,A)
call dump_snapshot_(s, filename)
end if
j_ad = 0
advance_loop : do
j_ad = j_ad + 1
! Check whether we've reached the edge of the grid
if (t_bdy .AND. zi(i_bdy)%kc_t == 1) then
zi(i_bdy)%sel_t = .FALSE.
elseif (.NOT. t_bdy .AND. zi(i_bdy)%kc_b == s%nz) then
zi(i_bdy)%sel_b = .FALSE.
endif
! Check whether the advancing face is still selected
if ((t_bdy .AND. .NOT. zi(i_bdy)%sel_t) .OR. &
(.NOT. t_bdy .AND. .NOT. zi(i_bdy)%sel_b)) exit advance_loop
! Advance the boundary by one cell
call advance_bdy_by_one_cell_(s, update_mode, zi, i_bdy, t_bdy, sd)
! If necessary, dump a snapshot
if (s% conv_premix_dump_snapshots .AND. MOD(j_ad, 50) == 0) then
write(filename, 100) 'cpm-snap.', s%model_number, '.it_', j_it, '.ad_', j_ad, '.dat'
call dump_snapshot_(s, filename)
end if
! Loop around
end do advance_loop
if (s% conv_premix_dump_snapshots) then
write(filename, 100) 'cpm-snap.', s%model_number, '.it_', j_it, '.ad_', j_ad-1, '.dat'
call dump_snapshot_(s, filename)
end if
if (TRACE_ADVANCE_BDY) then
write(*,*) 'Done advancing zone boundary'
end if
! Finish
return
end subroutine advance_bdy_
!****
subroutine advance_bdy_by_one_cell_ (s, update_mode, zi, i_bdy, t_bdy, sd)
type(star_info), pointer :: s
integer, intent(inout) :: update_mode(:)
type(zone_info), allocatable, intent(inout) :: zi(:)
integer, intent(inout) :: i_bdy
logical, intent(in) :: t_bdy
type(saved_data), intent(inout) :: sd
integer, parameter :: N_SUB_INI = 1
integer, parameter :: N_SUB_MAX = 1024
logical, parameter :: TRACE_MIX_CELL = TRACE_ALL
logical, parameter :: TRACE_MIX_SUBCELL = TRACE_ALL
real(dp) :: dr
real(dp) :: vp
real(dp) :: delta_dt
real(dp) :: dt_limit
integer :: kc_new
integer :: n_sub
integer :: kc_t
integer :: kc_b
integer :: kf_t
integer :: kf_b
real(dp), allocatable :: xa_sub(:,:)
real(dp) :: m_sub
logical :: has_split
integer :: kc_sub
real(dp) :: m_zone
integer :: l
real(dp) :: avg_xa(s%species)
integer :: n_nc
integer :: kf
real(dp) :: davg_xa_dt(s%species)
! Advance the top (bottom) boundary of the zone zi(i_bdy) by one
! cell
! Calculate the clock increment for adding the new cell
if (t_bdy) then
kc_t = zi(i_bdy)%kc_t
dr = s%r(kc_t) - s%r(kc_t+1)
vp = zi(i_bdy)%vp_t
else
kc_b = zi(i_bdy)%kc_b
if (kc_b < s%nz) then
dr = s%r(kc_b) - s%r(kc_b+1)
else
dr = s%r(kc_b)
endif
vp = zi(i_bdy)%vp_b
endif
delta_dt = dr/vp
! Check whether we have enough time; if not, return
if (s%conv_premix_time_factor > 0._dp) then
dt_limit = s%conv_premix_time_factor*s%dt
else
dt_limit = HUGE(0._dp)
endif
if (t_bdy) then
if (zi(i_bdy)%dt_t + delta_dt > dt_limit) then
zi(i_bdy)%sel_t = .FALSE.
if (TRACE_MIX_CELL) then
write(*,*) ' Outcome: out of time', zi(i_bdy)%dt_t, delta_dt, dt_limit
end if
return
else
zi(i_bdy)%dt_t = zi(i_bdy)%dt_t + delta_dt
endif
else
if (zi(i_bdy)%dt_b + delta_dt > dt_limit) then
zi(i_bdy)%sel_b = .FALSE.
if (TRACE_MIX_CELL) then
write(*,*) ' Outcome: out of time', zi(i_bdy)%dt_b, delta_dt, dt_limit
end if
return
else
zi(i_bdy)%dt_b = zi(i_bdy)%dt_b + delta_dt
endif
endif
! Determine the index of the new cell we're going to add
if (t_bdy) then
kc_new = zi(i_bdy)%kc_t - 1
else
kc_new = zi(i_bdy)%kc_b + 1
endif
! Save the current model state over all cells we *might* modify
if (t_bdy) then
call save_model_(s, update_mode, kc_new, zi(i_bdy)%kc_b, sd)
else
call save_model_(s, update_mode, zi(i_bdy)%kc_t, kc_new, sd)
endif
! Loop over subcell refinement levels
n_sub = N_SUB_INI
refine_loop : do
! Set initial cell and face indices
kc_t = zi(i_bdy)%kc_t
kc_b = zi(i_bdy)%kc_b
kf_t = kc_t
kf_b = kc_b + 1
if (TRACE_MIX_CELL) then
write(*,*) 'Attempting to add cell via n_sub subcells:', n_sub
write(*,*) ' kc_t :', kc_t
write(*,*) ' kc_b :', kc_b
write(*,*) ' kc_new :', kc_new
end if
! Set update_mode over the cells being mixed
if (s%conv_premix_fix_pgas) then
update_mode(kc_t:kc_b) = FIXED_PT_MODE
else
update_mode(kc_t:kc_b) = FIXED_DT_MODE
endif
! Set gradL_composition_term to zero on interior faces plus the
! advancing face
if (t_bdy) then
s%gradL_composition_term(kf_t:kf_b-1) = 0._dp
else
s%gradL_composition_term(kf_t+1:kf_b) = 0._dp
endif
! Divide the new cell into n_sub virtual subcells, and store the
! (initially uniform) abundances of these subcells
if (ALLOCATED(xa_sub)) deallocate(xa_sub)
allocate(xa_sub(s%species,n_sub))
do l = 1, s%species
xa_sub(l,:) = s%xa(l,kc_new)
end do
m_sub = s%dm(kc_new)/n_sub
! Mix subcells into the zone, one by one
has_split = .FALSE.
subcell_loop : do kc_sub = 1, n_sub
if (TRACE_MIX_SUBCELL) then
write(*,*) 'In subcell loop:', kc_sub, n_sub
end if
! Determine average abundances across the zone and the 1:k_s
! subcells
m_zone = SUM(s%dm(kc_t:kc_b))
do l = 1, s%species
avg_xa(l) = (SUM(s%dm(kc_t:kc_b)*s%xa(l,kc_t:kc_b)) + &
SUM(m_sub*xa_sub(l,1:kc_sub))) / (m_zone + m_sub*kc_sub)
end do
avg_xa = MAX(MIN(avg_xa, 1._dp), 0._dp)
avg_xa = avg_xa/SUM(avg_xa)
! Update abundances using the average
do l = 1, s%species
s%xa(l,kc_t:kc_b) = avg_xa(l)
xa_sub(l,1:kc_sub) = avg_xa(l)
end do
call update_model_(s, update_mode, kc_t, kc_b)
! Determine how many interior faces are now non-convective
n_nc = COUNT(s%mlt_mixing_type(kf_t+1:kf_b-1) /= convective_mixing)
if (n_nc > 1 .AND. n_sub < N_SUB_MAX) then
! More than a single one: double the number of subcells
! and restart (as long as we're below the subcell limit)
n_sub = 2*n_sub
call restore_model_(s, update_mode, sd)
cycle refine_loop
elseif (n_nc > 0) then
! A single one (or more, if we're over the subcell
! limit): move the trailing boundary to its new position
! (defined as the closest non-convective face to the
! advancing boundary)
has_split = .TRUE.
if (t_bdy) then
search_down_loop : do kf = kf_t+1, kf_b-1
if (s%mlt_mixing_type(kf) /= convective_mixing) exit search_down_loop
end do search_down_loop
kc_b = kf - 1
kf_b = kc_b + 1
if (TRACE_MIX_SUBCELL) then
write(*,*) ' Moved lower boundary to', kc_b
end if
else
search_up_loop : do kf = kf_b-1, kf_t+1, -1
if (s%mlt_mixing_type(kf) /= convective_mixing) exit search_up_loop
end do search_up_loop
kc_t = kf
kf_t = kc_t
if (TRACE_MIX_SUBCELL) then
write(*,*) ' Moved upper boundary to', kc_t
end if
endif
end if
end do subcell_loop
! If we reach this point, all went well; exit
if (TRACE_MIX_SUBCELL) then
write(*,*) ' Completed mixing subcell'
end if
exit refine_loop
end do refine_loop
! Update the abundances in the new cell, and adjust the cell/face
! indices
s%xa(:,kc_new) = avg_xa
if (s%conv_premix_fix_pgas) then
update_mode(kc_new) = FIXED_PT_MODE
else
update_mode(kc_new) = FIXED_DT_MODE
endif
if (t_bdy) then
kc_t = kc_new
kf_t = kc_t
else
kc_b = kc_new
kf_b = kc_b + 1
endif
call update_model_(s, update_mode, kc_t, kc_b)
! Check whether an abundance increase has occurred; if so, revert
! back to the starting model and return
if (zi(i_bdy)%avoid_inc_iso /= 0) then
! Evaluate the rate of change of the zone-average abundances
davg_xa_dt = (avg_xa - zi(i_bdy)%avg_xa)/s%dt
! Check whether an increase will occur (davg_xa_dt is positive,
! and greater in magnitude than the rate of decrease due to
! burning)
associate (iso => zi(i_bdy)%avoid_inc_iso)
if (davg_xa_dt(iso) > MAX(-zi(i_bdy)%davg_xa_dt(iso), 0._dp)) then
call restore_model_(s, update_mode, sd)
if (t_bdy) then
zi(i_bdy)%sel_t = .FALSE.
else
zi(i_bdy)%sel_b = .FALSE.
endif
if (TRACE_MIX_CELL) then
write(*,*) ' Outcome: abundance reversal for isotope', iso
end if
return
end if
end associate
end if
! Determine whether the face just inside the advancing boundary
! has become/remained convective; if not, revert back to the
! starting model and return
if (t_bdy) then
if (s%mlt_mixing_type(kf_t+1) /= convective_mixing) then
call restore_model_(s, update_mode, sd)
zi(i_bdy)%sel_t = .FALSE.
if (TRACE_MIX_CELL) then
write(*,*) ' Outcome: non-convective at top'
end if
return
endif
else
if (s%mlt_mixing_type(kf_b-1) /= convective_mixing) then
call restore_model_(s, update_mode, sd)
zi(i_bdy)%sel_b = .FALSE.
if (TRACE_MIX_CELL) then
write(*,*) ' Outcome: non-convective at bottom'
end if
return
endif
end if
! Update the zone info table to reflect all the changes we made
call update_zone_info_(s, i_bdy, t_bdy, has_split, zi)
! Finish
return
end subroutine advance_bdy_by_one_cell_
!****
subroutine init_zone_info_ (s, zi)
type(star_info), pointer :: s
type(zone_info), allocatable, intent(out) :: zi(:)
integer :: i
! Initialize the zone info table for the whole star
call create_zone_info_(s, 1, s%nz, zi)
! Perform additional initializations
zone_loop : do i = 1, SIZE(zi)
! Set penetration velocities
call set_penetration_velocities_(s, zi(i))
! Set burn types and initial average abundances
call set_burn_data_(s, zi(i))
call set_abund_data_(s, zi(i))
! Initialize clocks
zi(i)%dt_t = 0._dp
zi(i)%dt_b = 0._dp
! Set flags
zi(i)%sel_t = zi(i)%kc_t /= 1
zi(i)%sel_b = zi(i)%kc_b /= s%nz
zi(i)%split_retreat = .FALSE.
end do zone_loop
! Finish
return
end subroutine init_zone_info_
!****
subroutine create_zone_info_ (s, kc_t, kc_b, zi)
type(star_info), pointer :: s
integer, intent(in) :: kc_t
integer, intent(in) :: kc_b
type(zone_info), allocatable, intent(out) :: zi(:)
logical :: in_conv
type(zone_info) :: zi_new
integer :: kf
! Create a zone info table spanning the cells kc_t:kc_b. For the
! purposes of convective premixing, a zone is defined as a regions
! of two or more cells with (i) convective faces inside the
! region, (ii) non-convective faces at the top and bottom
! boundaries of the region. The possible exceptions to these
! rules apply to zones which begin (end) at kc_t (kc_b); their
! upper (lower) boundary faces are not checked for being
! non-convective.
if (kc_t >= kc_b) then
write(*,*) 'conv_premix: Invalid cell range in create_zone_info_'
stop
endif
! Create the empty zone info table
allocate(zi(0))
! Loop through faces, looking for zone boundaries and adding zones
! (cf. locate_convection_boundaries in mix_info.f90)
kf = kc_b
in_conv = s%mlt_mixing_type(kf) == convective_mixing
if (in_conv) then
zi_new%kc_b = kf
zi_new%vc_b = s%mlt_vc(kf)
endif
face_loop : do kf = kc_b-1, kc_t+1, -1
if (in_conv .AND. s%mlt_mixing_type(kf) /= convective_mixing) then
! Transitioning out of a convective zone; complete
! definition of new zone info and add to to the table
zi_new%kc_t = kf
zi_new%vc_t = s%mlt_vc(kf+1)
zi = [zi,zi_new]
in_conv = .FALSE.
elseif (.NOT. in_conv .AND. s%mlt_mixing_type(kf) == convective_mixing) then
! Transitioning into a convective zone; begin definition of
! new zone info
zi_new%kc_b = kf
zi_new%vc_b = s%mlt_vc(kf)
in_conv = .TRUE.
endif
end do face_loop
kf = kc_t
if (in_conv) then
zi_new%kc_t = kf
zi_new%vc_t = s%mlt_vc(kf+1)
zi = [zi,zi_new]
end if
! Finish
return
end subroutine create_zone_info_
!****
subroutine update_zone_info_ (s, i_bdy, t_bdy, has_split, zi)
type(star_info), pointer :: s
integer, intent(inout) :: i_bdy
logical, intent(in) :: t_bdy
logical, intent(in) :: has_split
type(zone_info), allocatable, intent(inout) :: zi(:)
logical, parameter :: TRACE_UPDATE_ZONE = TRACE_ALL
type(zone_info), allocatable :: zi_new(:)
integer :: i
! Update the zone info table to account for a single-cell
! advancement in zone i_bdy, direction t_bdy
! Advance the zone by one cell
if (t_bdy) then
zi(i_bdy)%kc_t = zi(i_bdy)%kc_t - 1
zi(i_bdy)%vc_t = s%mlt_vc(zi(i_bdy)%kc_t+1)
else
zi(i_bdy)%kc_b = zi(i_bdy)%kc_b + 1
zi(i_bdy)%vc_b = s%mlt_vc(zi(i_bdy)%kc_b)
endif
call set_penetration_velocities_(s, zi(i_bdy))
! If necessary, truncate (or even delete) the adjacent zone we've
! encroached into
if (t_bdy .AND. i_bdy < SIZE(zi)) then
if (zi(i_bdy)%kc_t == zi(i_bdy+1)%kc_b) then
! Encroached into zone above
if (zi(i_bdy+1)%kc_b-zi(i_bdy+1)%kc_t > 1) then
! Truncate the zone
zi(i_bdy+1)%kc_b = zi(i_bdy+1)%kc_b - 1
zi(i_bdy+1)%vc_b = s%mlt_vc(zi(i_bdy+1)%kc_b)
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Truncated zone above to', zi(i_bdy+1)%kc_b, zi(i_bdy+1)%vc_b, &
s%mlt_mixing_type(zi(i_bdy+1)%kc_b), zi(i_bdy+1)%kc_b-zi(i_bdy+1)%kc_t+1
endif
else
! Delete the zone
zi = [zi(:i_bdy),zi(i_bdy+2:)]
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Deleted zone above'
endif
end if
end if
elseif (.NOT. t_bdy .AND. i_bdy > 1) then
if (zi(i_bdy)%kc_b == zi(i_bdy-1)%kc_t) then
! Encroached into zone below
if (zi(i_bdy-1)%kc_b-zi(i_bdy-1)%kc_t > 1) then
! Truncate the zone
zi(i_bdy-1)%kc_t = zi(i_bdy-1)%kc_t + 1
zi(i_bdy-1)%vc_t = s%mlt_vc(zi(i_bdy-1)%kc_t+1)
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Truncated zone below to', zi(i_bdy-1)%kc_t
endif
else
! Delete the zone
zi = [zi(:i_bdy-2),zi(i_bdy:)]
i_bdy = i_bdy - 1
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Deleted zone below'
endif
end if
end if
endif
! If the zone has split (or its trailing boundary has retreated),
! then replace the zone with as many new zones as necessary
if (has_split) then
! Create new zones to replace the zone
call create_zone_info_(s, zi(i_bdy)%kc_t, zi(i_bdy)%kc_b, zi_new)
! Set up parameters for the new zones
new_zone_loop : do i = 1, SIZE(zi_new)
! Clocks & selection flags
if (zi_new(i)%kc_t == zi(i_bdy)%kc_t) then
zi_new(i)%dt_t = zi(i_bdy)%dt_t
zi_new(i)%sel_t = zi(i_bdy)%sel_t
else
if (t_bdy) then
zi_new(i)%dt_t = zi(i_bdy)%dt_t
zi_new(i)%sel_t = .FALSE.
else
zi_new(i)%dt_t = zi(i_bdy)%dt_b
zi_new(i)%sel_t = .FALSE.
endif
endif
if (zi_new(i)%kc_b == zi(i_bdy)%kc_b) then
zi_new(i)%dt_b = zi(i_bdy)%dt_b
zi_new(i)%sel_b = zi(i_bdy)%sel_b
else
if (t_bdy) then
zi_new(i)%dt_b = zi(i_bdy)%dt_t
zi_new(i)%sel_b = .FALSE.
else
zi_new(i)%dt_b = zi(i_bdy)%dt_b
zi_new(i)%sel_b = .FALSE.
endif
endif
! Penetration velocities
call set_penetration_velocities_(s, zi_new(i))
! Burn data
call set_burn_data_(s, zi_new(i))
! Initial abundances
zi_new(i)%avg_xa = zi(i_bdy)%avg_xa
zi_new(i)%davg_xa_dt = zi(i_bdy)%davg_xa_dt
end do new_zone_loop
! Replace the zone with the new zones
zi = [zi(:i_bdy-1),zi_new,zi(i_bdy+1:)]
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Replaced zone with new zones', SIZE(zi), SIZE(zi_new)
end if
else
if (TRACE_UPDATE_ZONE) then
write(*,*) 'Zone did not split, updating indices',&
COUNT(s%mlt_mixing_type(zi(i_bdy)%kc_t+1:zi(i_bdy)%kc_b) /= convective_mixing)
end if
end if
! Finish
return
end subroutine update_zone_info_
!****
subroutine set_penetration_velocities_ (s, zi)
type(star_info), pointer :: s
type(zone_info), intent(inout) :: zi
real(dp) :: mu_i
real(dp) :: mu_e
integer :: kf
real(dp) :: z_s
! Set the penetration velocities for the zone info. These are
! calculated by evaluating the molecular weight contrast between
! the boundary cell and the adjacent (radiative) cell, and
! calculating a characteristic penetration velocity using the
! approach by Castellani (1971)