-
Notifications
You must be signed in to change notification settings - Fork 9
/
finite.c
6979 lines (6037 loc) · 216 KB
/
finite.c
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
/*! \file finite.c
\brief Routines related to finite difference and grid
*/
#include "ko.h"
//**********************************************************************
/*! \fn int reduce_order_check(ldouble *pm2,ldouble *pm1,ldouble *p0,ldouble *pp1,ldouble *pp2,int ix,int iy,int iz)
\brief reduces reconstruction order to first order for small time steps or within certain radius
\param[in] pm2, pm1, p0, pp1, pp2 Pointers to primitives in cell and neighbors
\param[in] ix,iy,iz Cell indices
*/
//**********************************************************************
int
reduce_order_check(ldouble *pm2,ldouble *pm1,ldouble *p0,ldouble *pp1,ldouble *pp2,int ix,int iy,int iz)
{
int reconstrpar=0;
#ifdef REDUCEORDERTEMP
ldouble t0,tp1,tm1,tmin;
t0 =calc_PEQ_Tfromurho(p0[UU], p0[RHO]);
tp1=calc_PEQ_Tfromurho(pp1[UU],pp1[RHO]);
tm1=calc_PEQ_Tfromurho(pm1[UU],pm1[RHO]);
tmin=my_min(t0,my_min(tp1,tm1));
if(tmin<REDUCEORDERTEMP)
reconstrpar=1;
#endif
#ifdef REDUCEORDERRADIUS
ldouble xxBL[4];
get_xx_arb(ix,iy,iz,xxBL,BLCOORDS);
if(xxBL[1]<REDUCEORDERRADIUS)
reconstrpar=1;
#endif
#ifdef REDUCEORDERATBH
ldouble xxBL[4];
get_xx_arb(ix,iy,iz,xxBL,BLCOORDS);
if(xxBL[1]<rhorizonBL)
reconstrpar=1;
#endif
return reconstrpar;
}
//**********************************************************************
/*! \fn ldouble reduce_minmod_theta(ldouble *pm2,ldouble *pm1,ldouble *p0,ldouble *pp1,ldouble *pp2,int ix,int iy,int iz)
\brief reduces minmod theta value to 1 near axis and inner boundary
\param[in] pm2, pm1, p0, pp1, pp2 Pointers to primitives in cell and neighbors
\param[in] ix,iy,iz Cell indices
*/
//**********************************************************************
ldouble
reduce_minmod_theta(ldouble *pm2,ldouble *pm1,ldouble *p0,ldouble *pp1,ldouble *pp2,int ix,int iy,int iz)
{
ldouble theta=MINMOD_THETA;
#if (REDUCEMINMODTHETA==1) //reduce near the axis
int gix,giy,giz;
giy=iy+TOJ;
int limit = NCCORRECTPOLAR + 1;
if(giy<limit || giy> (TNY-limit-1))
theta=1.0;
#endif
#if (REDUCEMINMODTHETA==2) //reduce near the inner boundary
int gix,giy,giz;
gix=ix+TOI;
int limit = NCELLSREDUCEMINMODTHETA;
if(gix<limit)
theta=1.0;
#endif
return theta;
}
//**********************************************************************
/*! \fn int avg2point(ldouble *um2,ldouble *um1,ldouble *u0,ldouble *up1,ldouble *up2,ldouble *ul,ldouble *ur,ldouble dxm2,ldouble dxm1,ldouble dx0,ldouble dxp1,ldouble dxp2,int param,ldouble theta)
\brief Interpolates primitives to the left and right walls of current cell i
@param[in] um2, um1, u0, up1, up2 values of primitive in the i-2, i-1, i, i+1, i+2 cells
@param[out] ul, ur interpolated primitives at the left and right walls of cell i
@param[in] dxm2, dxm1, dx0, dxp1, dxp2 sizes of the five cells
@param[in] param reconstrpar
@param[in] minmod_theta MINMOD_THETA
Several interpolation schemes are available.
INT_ORDER=0: basic donor cell\n
INT_ORDER=1: Minmod (FLUXLIMITER=0), Monotonized Central (FLUXLIMITER=1), Superbee (FLUXLIMITER=4)\n
INT_ORDER=2: Piecewise Parabolic Method (PPM)\n
*/
//**********************************************************************
int
avg2point(ldouble *um2,ldouble *um1,ldouble *u0,ldouble *up1,ldouble *up2,
ldouble *ul,ldouble *ur,
ldouble dxm2,ldouble dxm1,ldouble dx0,ldouble dxp1,ldouble dxp2,
int param,ldouble theta)
{
ldouble r0[NV],rm1[NV],rp1[NV];
int int_order_local = INT_ORDER;
if(param!=0) //reduce integration order
{
int_order_local = INT_ORDER-param;
} // if(param!=0)
// default to donor cell
if(int_order_local<0)
int_order_local = 0;
if(int_order_local==0) //DONOR CELL
{
int iv;
for(iv=0;iv<NV;iv++)
{
ur[iv]=u0[iv];
ul[iv]=u0[iv];
}
} // else if(INT_ORDER==0)
else if(int_order_local==1) //linear
{
ldouble diffpar=theta;
int iv;
for(iv=0;iv<NV;iv++)
{
// Slope limiter code rewritten by Ramesh. No function-calls, no divisions.
ldouble slope;
ldouble deltam = u0[iv]-um1[iv];
ldouble deltap = up1[iv]-u0[iv];
if (deltam * deltap <= 0.)
{
// We are at a local maximum or minimum. Use zero slope (i.e., donor cell)
ur[iv] = u0[iv];
ul[iv] = u0[iv];
}
else
{
if (deltam > 0.)
{
// Slopes are positive. Choose limiter appropriately
if (FLUXLIMITER == 0) // MinMod
{
slope = my_min(my_min(theta*deltam, 0.5*(deltam+deltap)), theta*deltap); // theta=1 is MinMod, theta=2 is MC
}
else if (FLUXLIMITER == 1) // MC
{
slope = my_min(my_min(2*deltam, 0.5*(deltam+deltap)), 2*deltap);
}
else if (FLUXLIMITER == 2) // Osher -- discouraged since it is not symmetric
{
printf("Error: Osher slope limiter is discouraged since it is not symmetric\n");
exit(-2);
}
else if (FLUXLIMITER == 3) // Koren -- discouraged since it is not symmetric
{
printf("Error: Koren slope limiter is discouraged since it is not symmetric\n");
exit(-3);
}
else // Superbee
{
slope = my_max(my_min(2*deltam, deltap), my_min(deltam, 2*deltap));
}
}
else
{
// Slopes are negative. Choose limiter appropriately
if (FLUXLIMITER == 0) // MinMod
{
slope = my_max(my_max(theta*deltam, 0.5*(deltam+deltap)), theta*deltap); // theta=1 is MinMod, theta=2 is MC
}
else if (FLUXLIMITER == 1) // MC
{
slope = my_max(my_max(2*deltam, 0.5*(deltam+deltap)), 2*deltap);
}
else if (FLUXLIMITER == 2) // Osher -- discouraged since it is not symmetric
{
printf("Error: Osher slope limiter is discouraged since it is not symmetric\n");
exit(-2);
}
else if (FLUXLIMITER == 3) // Koren -- discouraged since it is not symmetric
{
printf("Error: Koren slope limiter is discouraged since it is not symmetric\n");
exit(-3);
}
else // Superbee
{
slope = my_min(my_max(2*deltam, deltap), my_max(deltam, 2*deltap));
}
}
ur[iv] = u0[iv] + 0.5*slope;
ul[iv] = u0[iv] - 0.5*slope;
}
if(isnan(ur[iv]) || isnan (ul[iv])) printf("interperr %d %e %e %e %e %e\n",iv,um2[iv],um1[iv],u0[iv],up1[iv],up2[iv]);
//u0 remains intact - in linear reconstruction cell averaged equals cell centered
} // for(iv=0;iv<NV;iv++)
} // else if(INT_ORDER==1)
else if(int_order_local==2) //parabolic PPM
{
//The following is based on Colella & Woodward (J. Comp. Phys. 54, 174, 1984).
//It uses five points: m2, m1, 0, p1, p2.
//The code has been checked and verified by Ramesh: July 14, 2017
// Define various quantities that apear in the formula
ldouble dxp2_plus_dxp1 = dxp2 + dxp1;
ldouble dxp2_plus_dxp1_inv = 1. / dxp2_plus_dxp1;
ldouble dxp1_plus_dx0 = dxp1 + dx0;
ldouble dxp1_plus_dx0_inv = 1. / dxp1_plus_dx0;
ldouble dx0_plus_dxm1 = dx0 + dxm1;
ldouble dx0_plus_dxm1_inv = 1. / dx0_plus_dxm1;
ldouble dxm1_plus_dxm2 = dxm1 + dxm2;
ldouble dxm1_plus_dxm2_inv = 1. / dxm1_plus_dxm2;
ldouble dxm1_plus_dx0_plus_dxp1_inv = 1. / (dxm1+dx0+dxp1);
ldouble dx0_plus_dxp1_plus_dxp2_inv = 1. / (dx0+dxp1+dxp2);
ldouble dxm2_plus_dxm1_plus_dx0_inv = 1. / (dxm2+dxm1+dx0);
ldouble dx0_plus_twodxm1 = dx0 + 2. * dxm1;
ldouble dxp1_plus_twodx0 = dxp1 + 2. * dx0;
ldouble dxm1_plus_twodxm2 = dxm1 + 2. * dxm2;
ldouble twodxp1_plus_dx0 = 2. * dxp1 + dx0;
ldouble twodxp2_plus_dxp1 = 2. * dxp2 + dxp1;
ldouble twodx0_plus_dxm1 = 2. * dx0 + dxm1;
//ldouble l,r;
int iv;
ldouble dri[NV],drim1[NV],drip1[NV];
for(iv=0;iv<NV;iv++)
{
// dri, drip1, drim1 are the slopes delta a_j, delta a_{j+1}, delta a_{j-1} in eq (1.7) of C&W
dri[iv] = dx0 * dxm1_plus_dx0_plus_dxp1_inv *
(dx0_plus_twodxm1 * dxp1_plus_dx0_inv * (up1[iv]-u0[iv]) +
twodxp1_plus_dx0 * dx0_plus_dxm1_inv * (u0[iv]-um1[iv]));
drip1[iv] = dxp1 * dx0_plus_dxp1_plus_dxp2_inv *
(dxp1_plus_twodx0 * dxp2_plus_dxp1_inv * (up2[iv]-up1[iv]) +
twodxp2_plus_dxp1 * dxp1_plus_dx0_inv * (up1[iv]-u0[iv]));
drim1[iv] = dxm1 * dxm2_plus_dxm1_plus_dx0_inv *
(dxm1_plus_twodxm2 * dx0_plus_dxm1_inv * (u0[iv]-um1[iv]) +
twodx0_plus_dxm1 * dxm1_plus_dxm2_inv * (um1[iv]-um2[iv]));
// Limit the slopes to be monotonic. This is eq (1.8) in C&W. (Note a minor typo in C&W: one of their _{j-1} should be _{j+1})
if( (up1[iv]-u0[iv]) * (u0[iv]-um1[iv]) > 0.)
{
dri[iv] = my_min(fabs(dri[iv]), my_min(2. * fabs(u0[iv] - um1[iv]), 2. * fabs(u0[iv] - up1[iv]))) * my_sign(dri[iv]);
}
else
{
dri[iv]=0.;
}
if( (up2[iv]-up1[iv]) * (up1[iv]-u0[iv]) > 0.)
{
drip1[iv] = my_min(fabs(drip1[iv]), my_min(2. * fabs(up1[iv] - u0[iv]), 2. * fabs(up1[iv] - up2[iv]))) * my_sign(drip1[iv]);
}
else
{
drip1[iv]=0.;
}
if( (u0[iv]-um1[iv]) * (um1[iv]-um2[iv]) > 0.)
{
drim1[iv] = my_min(fabs(drim1[iv]), my_min(2. * fabs(um1[iv] - um2[iv]), 2. * fabs(um1[iv] - u0[iv]))) * my_sign(drim1[iv]);
}
else
{
drim1[iv]=0.;
}
}
// Work on the right face of cell j
ldouble Z1, Z2, DX_inv;
Z1 = dx0_plus_dxm1 / dxp1_plus_twodx0;
Z2 = dxp2_plus_dxp1 / twodxp1_plus_dx0;
DX_inv = 1. / (dxm1+dx0+dxp1+dxp2);
for(iv=0;iv<NV;iv++)
{
// This is a_{j+1/2) in eq (1.6) of Colella & Woodward
ur[iv] = u0[iv] + dx0 * dxp1_plus_dx0_inv * (up1[iv]-u0[iv]) +
DX_inv * ((2.*dxp1*dx0) * dxp1_plus_dx0_inv * (Z1-Z2) * (up1[iv]-u0[iv]) -
dx0 * Z1 * drip1[iv] + dxp1 * Z2 * dri[iv]);
}
// Next work on the left face of cell j
Z1 = dxm1_plus_dxm2 / dx0_plus_twodxm1;
Z2 = dxp1_plus_dx0 / twodx0_plus_dxm1;
DX_inv = 1. / (dxm2+dxm1+dx0+dxp1);
for(iv=0;iv<NV;iv++)
{
// This is a_{j-1/2} in eq (1.6) of Colella & Woodward
ul[iv] = um1[iv] + dxm1 * dx0_plus_dxm1_inv * (u0[iv]-um1[iv]) +
DX_inv * ((2.*dx0*dxm1) * dx0_plus_dxm1_inv * (Z1-Z2) * (u0[iv]-um1[iv]) -
dxm1 * Z1 * dri[iv] + dx0 * Z2 * drim1[iv]);
}
// Make sure that the parabola remains monotonic. The following is equivalent to eq (1.10) in C&W, though it looks different
for(iv=0;iv<NV;iv++)
{
if((ur[iv]-u0[iv])*(u0[iv]-ul[iv])<=0.)
{
ul[iv] = u0[iv];
ur[iv] = u0[iv];
}
if((ur[iv] - ul[iv]) * (ul[iv] - (3. * u0[iv] - 2. * ur[iv])) < 0.)
{
ul[iv] = 3. * u0[iv] - 2.*ur[iv];
}
if((ur[iv] - ul[iv]) * ((3. * u0[iv] - 2. * ul[iv]) - ur[iv]) < 0.)
{
ur[iv] = 3. * u0[iv] - 2. * ul[iv];
}
//pass up reconstructed value at center - only if reconstructing average -> center
//To check consistency, make sure that u0[iv] in the LHS below is equal to u0[iv]
//u0[iv]=ul[iv]+.5*(ur[iv]-ul[iv] + 6.*(u0[iv]-.5*(ul[iv]+ur[iv]))*(1.-.5));
}
} // else if(INT_ORDER==2)
return 0;
}
//**********************************************************************
/*! \fn int calc_wavespeeds()
\brief Calculates and saves wavespeeds for cells within the domain plus ghost cells
Wavespeeds are in code units for the basic uniformly spaced coordinates x used in KORAL\n
For GRMHD problems, assumes that the wave dispersion relation in the fluid frame is isotropic,
with wave speed given by \f$v^2 = c_s^2 + v_A^2 - c_s^2 v_A^2\f$.
This is then transformed to the code coordinate frame.
Calculates left, right and maximum wavespeeds for hydro and radiation for each cell (ix, iy, iz).
These velocities are saved in array aaa.\n
The values in aaa are then transferred to global arrays:\n
-Hydro left/right: ahdxl, ahdxr, ahdyl, ahdyr, ahdzl, ahdzr\n
-Hydro maximum absolute: ahdx, ahdy, ahdz\n
-Radiation left/right: aradxl, aradxr, aradyl, aradyr, aradzl, aradzr\n
-Radiation maximum absolute: aradx, arady, aradz\n
*/
//**********************************************************************
int
calc_wavespeeds()
{
int ii;
#pragma omp parallel for
for(ii=0;ii<Nloop_1;ii++) //domain plus lim (=1 usually) ghost cells
{
int ix,iy,iz;
ix=loop_1[ii][0];
iy=loop_1[ii][1];
iz=loop_1[ii][2];
#ifndef MPI4CORNERS
if(if_outsidegc(ix,iy,iz)==1) continue; //avoid corners
#endif
ldouble aaa[24];
if(!is_cell_active(ix,iy,iz))
continue;
calc_wavespeeds_lr(ix,iy,iz,aaa);
save_wavespeeds(ix,iy,iz,aaa);
}
return 0;
}
//**********************************************************************
/*! \fn int save_wavespeeds(int ix, int iy, int iz, ldouble *aaa)
\brief saves characteristic wavespeeds from aaa[] to the global arrays. Also compute timesteps.
\param[in] ix,iy,iz cell indices
\param[in] aaa array with cell wavespeeds
*/
//**********************************************************************
int
save_wavespeeds(int ix,int iy,int iz, ldouble *aaa)
{
ldouble aaaxhd,aaaxrad,aaayhd,aaayrad,aaazhd,aaazrad;
ldouble aaaxrad2,aaayrad2,aaazrad2;
//hydro wavespeeds
set_u_scalar(ahdxl,ix,iy,iz,aaa[0]);
set_u_scalar(ahdxr,ix,iy,iz,aaa[1]);
set_u_scalar(ahdyl,ix,iy,iz,aaa[2]);
set_u_scalar(ahdyr,ix,iy,iz,aaa[3]);
set_u_scalar(ahdzl,ix,iy,iz,aaa[4]);
set_u_scalar(ahdzr,ix,iy,iz,aaa[5]);
aaaxhd=my_max(fabs(aaa[0]),fabs(aaa[1]));
aaayhd=my_max(fabs(aaa[2]),fabs(aaa[3]));
aaazhd=my_max(fabs(aaa[4]),fabs(aaa[5]));
set_u_scalar(ahdx,ix,iy,iz,aaaxhd);
set_u_scalar(ahdy,ix,iy,iz,aaayhd);
set_u_scalar(ahdz,ix,iy,iz,aaazhd);
//original radiative wavespeeds in slots aaa[6::11] are used later
//speeds in slots[12::] are limited by the optical depth
//used to calculate the fluxes
set_u_scalar(aradxl,ix,iy,iz,aaa[6+6]);
set_u_scalar(aradxr,ix,iy,iz,aaa[6+7]);
set_u_scalar(aradyl,ix,iy,iz,aaa[6+8]);
set_u_scalar(aradyr,ix,iy,iz,aaa[6+9]);
set_u_scalar(aradzl,ix,iy,iz,aaa[6+10]);
set_u_scalar(aradzr,ix,iy,iz,aaa[6+11]);
aaaxrad=my_max(fabs(aaa[6+6]),fabs(aaa[6+7]));
aaayrad=my_max(fabs(aaa[6+8]),fabs(aaa[6+9]));
aaazrad=my_max(fabs(aaa[6+10]),fabs(aaa[6+11]));
set_u_scalar(aradx,ix,iy,iz,aaaxrad);
set_u_scalar(arady,ix,iy,iz,aaayrad);
set_u_scalar(aradz,ix,iy,iz,aaazrad);
//searching for the maximal unlimited wavespeed for setting the timestep
// here radiative wavespeeds not limited by the optical depth
aaaxrad=my_max(fabs(aaa[6]),fabs(aaa[7]));
aaayrad=my_max(fabs(aaa[8]),fabs(aaa[9]));
aaazrad=my_max(fabs(aaa[10]),fabs(aaa[11]));
ldouble dx=get_size_x(ix,0);
ldouble dy=get_size_x(iy,1);
ldouble dz=get_size_x(iz,2);
ldouble wsx,wsy,wsz;
wsx=aaaxhd;
wsy=aaayhd;
wsz=aaazhd;
#ifdef RADIATION
#ifndef SKIPRADEVOLUTION
wsx=my_max(aaaxhd,aaaxrad);
wsy=my_max(aaayhd,aaayrad);
wsz=my_max(aaazhd,aaazrad);
#endif
#endif
//determine the time step
//only domain cells
if(if_indomain(ix,iy,iz)==1)
{
ldouble tstepden,ws_ph;
if(NZ>1 && NY>1)
tstepden=wsx/dx + wsy/dy + wsz/dz;
else if(NZ==1 && NY>1)
tstepden=wsx/dx + wsy/dy;
else if(NY==1 && NZ>1)
tstepden=wsx/dx + wsz/dz;
else
tstepden=wsx/dx;
tstepden/=TSTEPLIM;
set_u_scalar(cell_tstepden,ix,iy,iz,tstepden);
//global variables for maximum/minimum (inverse) cell timestep
////#pragma omp critical
if(tstepden>tstepdenmax) tstepdenmax=tstepden;
if(tstepden<tstepdenmin) tstepdenmin=tstepden;
}
return 0;
}
//**********************************************************************
/*! \fn int save_timesteps
\brief saves individual timesteps to be used inside time stepping
Does not work with MPI
*/
//**********************************************************************
int
save_timesteps()
{
int ii;
ldouble dtminloc = BIG;
//#pragma omp parallel for reduction(min:dtminloc) // "min" option in reduction is not available in older OpenMP implementations
for(ii = 0; ii < Nloop_0; ii++)
{
int ix, iy, iz;
ix = loop_0[ii][0];
iy = loop_0[ii][1];
iz = loop_0[ii][2];
ldouble cell_dt_local = 1. / get_u_scalar(cell_tstepden, ix, iy, iz);
// timestep is shortest of current cell and neighbors.
// ANDREW: is this doing the right thing on the x boundaries?
#ifdef SHORTERTIMESTEP
ldouble dtm1, dtp1, dt;
if(ix > 0)
dtm1 = 1. / get_u_scalar(cell_tstepden, ix-1, iy, iz);
else
dtm1 = BIG;
if(ix < NX)
dtp1 = 1. / get_u_scalar(cell_tstepden, ix+1, iy, iz);
else
dtp1 = BIG;
dt = 1. / get_u_scalar(cell_tstepden, ix, iy, iz);
cell_dt_local = my_min(my_min(dtm1, dtp1), dt);
#endif
set_u_scalar(cell_dt,ix,iy,iz,cell_dt_local);
//find the shortest
//global variables for shortest cell timestep
if(cell_dt_local < dtminloc)
dtminloc = cell_dt_local;
}
return 0;
}
//**********************************************************************
/*! \fn int calc_u2p(int type, int setflags)
\brief Calculates all primitives from global u
\param[in] type, not currently used
\param[in] setflags, should always=1 to set flags for cell fixups
*/
//**********************************************************************
int
calc_u2p(int type, int setflags)
{
int ii;
//timer start
struct timespec temp_clock;
my_clock_gettime(&temp_clock);
start_u2ptime=(ldouble)temp_clock.tv_sec+(ldouble)temp_clock.tv_nsec/1.e9;
//calculates the primitives
#pragma omp parallel for schedule (static)
for(ii=0;ii<Nloop_0;ii++) //domain only
{
int ix,iy,iz;
ix=loop_0[ii][0];
iy=loop_0[ii][1];
iz=loop_0[ii][2];
//skip if cell is passive
if(!is_cell_active(ix,iy,iz))
continue;
calc_primitives(ix,iy,iz,type,setflags);
}
//fixup here hd and rad after inversions
cell_fixup(FIXUP_U2PMHD);
#ifdef RADIATION
cell_fixup(FIXUP_U2PRAD);
#endif
//re-set boundary conditions
set_bc(global_time,0);
//timer stop
my_clock_gettime(&temp_clock);
end_u2ptime=(ldouble)temp_clock.tv_sec+(ldouble)temp_clock.tv_nsec/1.e9;
return 0;
}
//**********************************************************************
/*! \fn int do_correct()
\brief Corrects conserved quantities on polar axis and neutron star surface
*/
//**********************************************************************
int
do_correct()
{
int ix,iy,iz,ii;
//correct the polar axis
#ifdef CORRECT_POLARAXIS
correct_polaraxis();
#endif
#ifdef CORRECT_POLARAXIS_3D
///correct_polaraxis(); //first fill at fixed phi including magnetic field
correct_polaraxis_3d(); //then overwrite velocities with phi averages
#endif
#ifdef SMOOTH_POLARAXIS
smooth_polaraxis();
#endif
//correct the neutron star surface
#ifdef CORRECT_NSSURFACE
correct_nssurface();
#endif
return 0;
}
//**********************************************************************
/*! \fn int op_explicit(ldouble t, ldouble dtin)
\brief Explicit evolution from time t to t+dtin
\param[in] t, dtin initial time and time step
Saves initial primitives and conserveds in ppreexplict, upreexplicit\n
Interpolates primitives to cell faces using the required interpolation scheme,
and calculates left- and right-biased conserveds and fluxes\n
Calculates wave speeds\n
Computes combined fluxes at the faces using the selected approximate Riemann solver\n
*/
//**********************************************************************
int
op_explicit(ldouble t, ldouble dtin)
{
int ix,iy,iz,iv,ii;
ldouble dt;
int perform_sweep = 1; //Brandon - Added for special conditions where a sweep should not be performed under SPECIAL_BC_CHECK
#ifdef SPECIAL_BC_CHECK
int giix,giiy,giiz;
#endif
// Save conserveds and primitives over domain + ghost (no corners)
copyi_u(1.,u,upreexplicit); //conserved quantities before explicit update
copyi_u(1.,p,ppreexplicit); //primitive quantities before explicit update
// calculates H/R and velocity averages
calc_avgs_throughout();
// calculates wavespeeds over the domain and ghost cells
calc_wavespeeds();
#ifndef SKIPEVOLUTION
//**********************************************************************
// Next interpolate to the cell walls and calculate left and right-biased fluxes
#pragma omp parallel for private(ii,ix,iy,iz,iv) schedule (static)
for(ii=0;ii<Nloop_1;ii++) //domain plus lim (=1 usually) ghost cells
{
ix=loop_1[ii][0];
iy=loop_1[ii][1];
iz=loop_1[ii][2];
#ifdef SPECIAL_BC_CHECK
giix = ix + TOI;
giiy = iy + TOJ;
giiz = iz + TOK;
#endif
#ifndef MPI4CORNERS
if(if_outsidegc(ix,iy,iz)==1) continue; //avoid corners
#endif
#ifdef SPECIAL_BC_CHECK
#include PR_BC_SPECIAL_LOOP
if(ret_val == 4)
{
continue; //Exclude 'corners' in stream region
}
#endif
//create arrays for interpolating conserved quantities
struct geometry geom;
//ldouble x0[3],x0l[3],x0r[3],xm1[3],xp1[3];
ldouble fd_r0[NV],fd_rm1[NV],fd_rp1[NV];
ldouble fd_u0[NV],fd_up1[NV],fd_up2[NV],fd_um1[NV],fd_um2[NV];
ldouble fd_p0[NV],fd_pp1[NV],fd_pp2[NV],fd_pm1[NV],fd_pm2[NV],fd_pm3[NV],fd_pp3[NV];
ldouble fd_s0[NV],fd_sp1[NV],fd_sp2[NV],fd_sm1[NV],fd_sm2[NV];
ldouble fd_uLr[NV],fd_uLl[NV],fd_uRl[NV],fd_uRr[NV];
ldouble fd_pLr[NV],fd_pLl[NV],fd_pRl[NV],fd_pRr[NV];
ldouble fd_pr[NV],fd_pl[NV];
ldouble fd_sLr[NV],fd_sLl[NV],fd_sRl[NV],fd_sRr[NV];
ldouble fd_fstarl[NV],fd_fstarr[NV],fd_dul[3*NV],fd_dur[3*NV],fd_pdiffl[NV],fd_pdiffr[NV];
ldouble a0[2],am1[2],ap1[2],al,ar,amax,dx;
ldouble ffRl[NV],ffRr[NV],ffLl[NV],ffLr[NV];
ldouble ffl[NV],ffr[NV];
ldouble dx0, dxm2, dxm1, dxp1, dxp2;
ldouble minmod_theta=MINMOD_THETA;
int reconstrpar;
int dol,dor;
int iv;
#ifdef TRANSMITTING_YBC
ldouble fd_prT[NV],fd_plT[NV]; //Brandon for getting fluxes from GC instead of D+GC as usual
#endif
//**********************************************************************
// x 'sweep'
//**********************************************************************
perform_sweep = 1;
#ifdef MPI4CORNERS
if(NX>1 && iy>=-1 && iy<NY+1 && iz>=-1 && iz<NZ+1 && perform_sweep == 1) //needed to calculate face fluxes for flux-CT divB enforcement
#else
if(NX>1 && iy>=0 && iy<NY && iz>=0 && iz<NZ && perform_sweep == 1)
#endif
{
dol=dor=1;
if(ix<0) dol=0;
if(ix>=NX) dor=0;
#ifdef SPECIAL_BC_CHECK //Don't do l/r fluxes when at GC - Brandon
#ifndef SEARCH_STREAM_BOUNDARY
if(TNY>1 && TNZ==1)
{
if(giiy >= STREAM_IYT && giiy <= STREAM_IYB)
{
if(giix == STREAM_IX || giix == (STREAM_IX+1) || giix == (STREAM_IX+2) || giix == (STREAM_IX+3)) dor=0;
}
}
else if(TNY==1 && TNZ>1)
{
if(giiz >= STREAM_IZT && giiz <= STREAM_IZB)
{
if(giix == STREAM_IX || giix == (STREAM_IX+1) || giix == (STREAM_IX+2) || giix == (STREAM_IX+3)) dor=0;
}
}
else if(TNY>1 && TNZ>1)
{
#ifndef STREAM_RING
if(giiy >= STREAM_IYT && giiy <= STREAM_IYB && giiz >= STREAM_IZT && giiz <= STREAM_IZB)
#else
if(giiy >= STREAM_IYT && giiy <= STREAM_IYB)
#endif
{
if(giix == STREAM_IX || giix == (STREAM_IX+1) || giix == (STREAM_IX+2) || giix == (STREAM_IX+3)) dor=0;
}
}
#endif
#endif
// is_cell_active is currently always 1
// skip flux calculation if not needed
if((ix==0 && is_cell_active(ix,iy,iz)==0) || (ix>0 && is_cell_active(ix,iy,iz)==0 && is_cell_active(ix-1,iy,iz)==0))
dol=0;
if((ix==NX-1 && is_cell_active(ix,iy,iz)==0) || (ix<NX-1 && is_cell_active(ix,iy,iz)==0 && is_cell_active(ix+1,iy,iz)==0))
dor=0;
// dx0, dxm1, dxp1 are x-sizes (wall to wall) of cells ix, ix-1, ix+1, dxm2m, dxp2 are sizes of cells ix-2, ix+2
dx0=get_size_x(ix,0);
dxm1=get_size_x(ix-1,0);
dxp1=get_size_x(ix+1,0);
if(INT_ORDER>1)
{
dxm2=get_size_x(ix-2,0);
dxp2=get_size_x(ix+2,0);
}
for(iv=0;iv<NV;iv++)
{
//fd_p0, fd_pp1, fd_pm1 are primitives at current, left and right cells, fd_pm2, fd_pp2 are for next two cells
fd_p0[iv] =get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix+1,iy,iz);
fd_pm1[iv]=get_u(p,iv,ix-1,iy,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix-2,iy,iz);
fd_pp2[iv]=get_u(p,iv,ix+2,iy,iz);
}
}
reconstrpar=0;
#ifdef REDUCEORDERWHENNEEDED
reconstrpar = reduce_order_check(fd_pm2,fd_pm1,fd_p0,fd_pp1,fd_pp2,ix,iy,iz);
#endif
minmod_theta=MINMOD_THETA;
#ifdef REDUCEMINMODTHETA // reduce minmod_theta near axis or inner boundary
minmod_theta = reduce_minmod_theta(fd_pm2,fd_pm1,fd_p0,fd_pp1,fd_pp2,ix,iy,iz);
#endif
// Interpolate primitives to the left and right walls of current cell: fd_pl, fd_pr
avg2point(fd_pm2,fd_pm1,fd_p0,fd_pp1,fd_pp2,fd_pl,fd_pr,dxm2,dxm1,dx0,dxp1,dxp2,reconstrpar,minmod_theta);
// if(ix>0)
if(dol) //no need to calculate at left face of first GC if dol=0
{
// Left wall of current cell: compute fluxes and save in array ffl[NV]
fill_geometry_face(ix,iy,iz,0,&geom);
check_floors_mhd(fd_pl,VELPRIM,&geom);
f_flux_prime(fd_pl,0,ix,iy,iz,ffl,1);
}
// if(ix<NX)
if(dor) //no need to calculate at right face of first GC if dor=0
{
// Right wall of current cell: compute fluxes and save in array ffr[NV]
fill_geometry_face(ix+1,iy,iz,0,&geom);
check_floors_mhd(fd_pr,VELPRIM,&geom);
f_flux_prime(fd_pr,0,ix+1,iy,iz,ffr,0);
}
//save interpolated values to memory
//Note that l and r of a given cell ix are the left and right wall of that cell
//whereas L and R of given ix are quantities to the left and right of wall ix
for(iv=0;iv<NV;iv++)
{
// Save fd_pl in array pbRx (Primitive_R) of wall ix
// Save fd_pr in array pbLx (Primitive_L) of wall ix+1
set_ubx(pbRx,iv,ix,iy,iz,fd_pl[iv]);
set_ubx(pbLx,iv,ix+1,iy,iz,fd_pr[iv]);
if(dol)
// Save ffl in array flRx (F_R) of wall ix
set_ubx(flRx,iv,ix,iy,iz,ffl[iv]);
if(dor)
// Save ffr in array flLx (F_L) of wall ix+1
set_ubx(flLx,iv,ix+1,iy,iz,ffr[iv]);
}
} // if(NX>1 && iy>=0 && iy<NY && iz>=0 && iz<NZ...)
//**********************************************************************
//y 'sweep'
//**********************************************************************
perform_sweep = 1;
#ifdef SPECIAL_BC_CHECK
//if(giix > STREAM_IX)
//if(giix > STREAM_IX && giix <= (STREAM_IX+1))
if(giix > STREAM_IX && giix <= (STREAM_IX+3))
{
#ifdef MPI4CORNERS
if(TNY>1 && TNZ==1)
{
if(giiy > (STREAM_IYT-1) && giiy < (STREAM_IYB+1)) perform_sweep = 0;
}
else if(TNY>1 && TNZ>1)
{
#ifndef STREAM_RING
if(giiz > (STREAM_IZT-1) && giiz < (STREAM_IZB-1))
{
if(giiy > (STREAM_IYT-1) && giiy < (STREAM_IYB+1)) perform_sweep = 0;
}
#else
if(giiy > (STREAM_IYT-1) && giiy < (STREAM_IYB+1)) perform_sweep = 0;
#endif
}
#else
if(TNY>1 && TNZ==1)
{
if(giiy >= (STREAM_IYT-1) && giiy <= (STREAM_IYB+1)) perform_sweep = 0;
}
else if(TNY>1 && TNZ>1)
{
#ifndef STREAM_RING
if(giiz >= (STREAM_IZT-1) && giiz <= (STREAM_IZB-1))
{
if(giiy >= (STREAM_IYT-1) && giiy <= (STREAM_IYB+1)) perform_sweep = 0;
}
#else
if(giiy >= (STREAM_IYT-1) && giiy <= (STREAM_IYB+1)) perform_sweep = 0;
#endif
}
#endif
}
#endif
#ifdef MPI4CORNERS
if(NY>1 && ix>=-1 && ix<NX+1 && iz>=-1 && iz<NZ+1 && perform_sweep == 1)
#else
if(NY>1 && ix>=0 && ix<NX && iz>=0 && iz<NZ && perform_sweep == 1)
#endif
{
dol=dor=1;
if(iy<0) dol=0;
if(iy>=NY) dor=0;
// is_cell_active is currently always 1
// skip flux calculation if not needed
if((iy==0 && is_cell_active(ix,iy,iz)==0) || (iy>0 && is_cell_active(ix,iy,iz)==0 && is_cell_active(ix,iy-1,iz)==0))
dol=0;
if((iy==NY-1 && is_cell_active(ix,iy,iz)==0) || (iy<NY-1 && is_cell_active(ix,iy,iz)==0 && is_cell_active(ix,iy+1,iz)==0))
dor=0;
// dx0, dxm1, dxp1 are y-sizes (wall to wall) of cells iy, iy-1, iy+1, dxm2m, dxp2 are sizes of cells iy-2, iy+2
dx0=get_size_x(iy,1);
dxm1=get_size_x(iy-1,1);
dxp1=get_size_x(iy+1,1);
if(INT_ORDER>1)
{
dxm2=get_size_x(iy-2,1);
dxp2=get_size_x(iy+2,1);
}
for(iv=0;iv<NV;iv++)
{
//fd_p0, fd_pp1, fd_pm1 are primitives at current, left and right cells, fd_pm2, fd_pp2 are for next two cells
fd_p0[iv]=get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
#ifdef TRANSMITTING_YBC
//copy iy=-1 to iy=0 for calculating fluxes
if(TJ==0 && iy==0)
{
fd_p0[iv]=get_u(p,iv,ix,iy-1,iz); //iy=0 -> iy=-1
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
}
if(TJ==0 && iy==1)
{
fd_p0[iv]=get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-2,iz); //iy=0 -> iy=-1
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
}
if(TJ==0 && iy==2)
{
fd_p0[iv]=get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-3,iz); //iy=0 -> iy=-1
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
}
//copy iy=NY to iy=NY-1 for calculating fluxes
if(TJ==NTY-1 && iy==(NY-1) )
{
fd_p0[iv]=get_u(p,iv,ix,iy+1,iz); //iy=NY-1 -> iy=NY
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
}
if(TJ==NTY-1 && iy==(NY-2) )
{
fd_p0[iv]=get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix,iy+2,iz); //iy=NY-1 -> iy=NY
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+2,iz);
}
}
if(TJ==NTY-1 && iy==(NY-3) )
{
fd_p0[iv]=get_u(p,iv,ix,iy,iz);
fd_pp1[iv]=get_u(p,iv,ix,iy+1,iz);
fd_pm1[iv]=get_u(p,iv,ix,iy-1,iz);
if(INT_ORDER>1)
{
fd_pm2[iv]=get_u(p,iv,ix,iy-2,iz);
fd_pp2[iv]=get_u(p,iv,ix,iy+3,iz); //iy=NY-1 -> iy=NY
}
}
#endif //TRANSMITTING_YBC
}
reconstrpar=0;
#ifdef REDUCEORDERWHENNEEDED
reconstrpar = reduce_order_check(fd_pm2,fd_pm1,fd_p0,fd_pp1,fd_pp2,ix,iy,iz);
#endif