-
Notifications
You must be signed in to change notification settings - Fork 345
/
AMReX_MLCellLinOp.H
2106 lines (1889 loc) · 80.3 KB
/
AMReX_MLCellLinOp.H
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
#ifndef AMREX_ML_CELL_LINOP_H_
#define AMREX_ML_CELL_LINOP_H_
#include <AMReX_Config.H>
#include <AMReX_MLLinOp.H>
#include <AMReX_iMultiFab.H>
#include <AMReX_MultiFabUtil.H>
#include <AMReX_YAFluxRegister.H>
#include <AMReX_MLLinOp_K.H>
#include <AMReX_MLMG_K.H>
#ifndef BL_NO_FORT
#include <AMReX_MLLinOp_F.H>
#endif
namespace amrex {
template <typename MF>
class MLCellLinOpT // NOLINT(cppcoreguidelines-virtual-class-destructor)
: public MLLinOpT<MF>
{
public:
using FAB = typename MF::fab_type;
using RT = typename MF::value_type;
using BCType = LinOpBCType;
using BCMode = typename MLLinOpT<MF>::BCMode;
using StateMode = typename MLLinOpT<MF>::StateMode;
using Location = typename MLLinOpT<MF>::Location;
MLCellLinOpT ();
~MLCellLinOpT () override = default;
MLCellLinOpT (const MLCellLinOpT<MF>&) = delete;
MLCellLinOpT (MLCellLinOpT<MF>&&) = delete;
MLCellLinOpT<MF>& operator= (const MLCellLinOpT<MF>&) = delete;
MLCellLinOpT<MF>& operator= (MLCellLinOpT<MF>&&) = delete;
void define (const Vector<Geometry>& a_geom,
const Vector<BoxArray>& a_grids,
const Vector<DistributionMapping>& a_dmap,
const LPInfo& a_info = LPInfo(),
const Vector<FabFactory<FAB> const*>& a_factory = {});
void setLevelBC (int amrlev, const MF* levelbcdata,
const MF* robinbc_a = nullptr,
const MF* robinbc_b = nullptr,
const MF* robinbc_f = nullptr) final;
using MLLinOpT<MF>::setLevelBC;
bool needsUpdate () const override {
return MLLinOpT<MF>::needsUpdate();
}
void update () override;
virtual bool isCrossStencil () const { return true; }
virtual bool isTensorOp () const { return false; }
void updateSolBC (int amrlev, const MF& crse_bcdata) const;
void updateCorBC (int amrlev, const MF& crse_bcdata) const;
virtual void applyBC (int amrlev, int mglev, MF& in, BCMode bc_mode, StateMode s_mode,
const MLMGBndryT<MF>* bndry=nullptr, bool skip_fillboundary=false) const;
BoxArray makeNGrids (int grid_size) const;
void restriction (int, int, MF& crse, MF& fine) const override;
void interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const override;
void interpAssign (int amrlev, int fmglev, MF& fine, MF& crse) const override;
void interpolationAmr (int famrlev, MF& fine, const MF& crse,
IntVect const& nghost) const override;
void averageDownSolutionRHS (int camrlev, MF& crse_sol, MF& crse_rhs,
const MF& fine_sol, const MF& fine_rhs) override;
void apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
StateMode s_mode, const MLMGBndryT<MF>* bndry=nullptr) const override;
void smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
bool skip_fillboundary=false) const final;
void solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
const MF* crse_bcdata=nullptr) override;
void correctionResidual (int amrlev, int mglev, MF& resid, MF& x, const MF& b,
BCMode bc_mode, const MF* crse_bcdata=nullptr) final;
// The assumption is crse_sol's boundary has been filled, but not fine_sol.
void reflux (int crse_amrlev,
MF& res, const MF& crse_sol, const MF&,
MF&, MF& fine_sol, const MF&) const final;
void compFlux (int amrlev, const Array<MF*,AMREX_SPACEDIM>& fluxes,
MF& sol, Location loc) const override;
void compGrad (int amrlev, const Array<MF*,AMREX_SPACEDIM>& grad,
MF& sol, Location loc) const override;
void applyMetricTerm (int amrlev, int mglev, MF& rhs) const final;
void unapplyMetricTerm (int amrlev, int mglev, MF& rhs) const final;
Vector<RT> getSolvabilityOffset (int amrlev, int mglev,
MF const& rhs) const override;
void fixSolvabilityByOffset (int amrlev, int mglev, MF& rhs,
Vector<RT> const& offset) const override;
void prepareForSolve () override;
RT xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool local) const final;
virtual void Fapply (int amrlev, int mglev, MF& out, const MF& in) const = 0;
virtual void Fsmooth (int amrlev, int mglev, MF& sol, const MF& rhs, int redblack) const = 0;
virtual void FFlux (int amrlev, const MFIter& mfi,
const Array<FAB*,AMREX_SPACEDIM>& flux,
const FAB& sol, Location loc, int face_only=0) const = 0;
virtual void addInhomogNeumannFlux (int /*amrlev*/,
const Array<MF*,AMREX_SPACEDIM>& /*grad*/,
MF const& /*sol*/,
bool /*mult_bcoef*/) const {}
RT normInf (int amrlev, MF const& mf, bool local) const override;
void averageDownAndSync (Vector<MF>& sol) const override;
void avgDownResAmr (int clev, MF& cres, MF const& fres) const override;
struct BCTL {
BoundCond type;
RT location;
};
protected:
bool m_has_metric_term = false;
Vector<std::unique_ptr<MLMGBndryT<MF>> > m_bndry_sol;
Vector<std::unique_ptr<BndryRegisterT<MF>> > m_crse_sol_br;
Vector<std::unique_ptr<MLMGBndryT<MF>> > m_bndry_cor;
Vector<std::unique_ptr<BndryRegisterT<MF>> > m_crse_cor_br;
// In case of agglomeration, coarse MG grids on amr level 0 are
// not simply coarsened from fine MG grids. So we need to build
// bcond and bcloc for each MG level.
using RealTuple = Array<RT,2*BL_SPACEDIM>;
using BCTuple = Array<BoundCond,2*BL_SPACEDIM>;
class BndryCondLoc
{
public:
BndryCondLoc (const BoxArray& ba, const DistributionMapping& dm, int ncomp);
void setLOBndryConds (const Geometry& geom, const Real* dx,
const Vector<Array<BCType,AMREX_SPACEDIM> >& lobc,
const Vector<Array<BCType,AMREX_SPACEDIM> >& hibc,
int ratio, const RealVect& interior_bloc,
const Array<Real,AMREX_SPACEDIM>& domain_bloc_lo,
const Array<Real,AMREX_SPACEDIM>& domain_bloc_hi);
const Vector<BCTuple>& bndryConds (const MFIter& mfi) const noexcept {
return bcond[mfi];
}
const Vector<RealTuple>& bndryLocs (const MFIter& mfi) const noexcept {
return bcloc[mfi];
}
const BCTuple& bndryConds (const MFIter& mfi, int icomp) const noexcept {
return bcond[mfi][icomp];
}
const RealTuple& bndryLocs (const MFIter& mfi, int icomp) const noexcept {
return bcloc[mfi][icomp];
}
GpuArray<BCTL,2*AMREX_SPACEDIM> const* getBCTLPtr (const MFIter& mfi) const noexcept {
return bctl[mfi];
}
private:
LayoutData<Vector<BCTuple> > bcond;
LayoutData<Vector<RealTuple> > bcloc;
LayoutData<GpuArray<BCTL,2*AMREX_SPACEDIM>*> bctl;
Gpu::DeviceVector<GpuArray<BCTL,2*AMREX_SPACEDIM> > bctl_dv;
int m_ncomp;
};
Vector<Vector<std::unique_ptr<BndryCondLoc> > > m_bcondloc;
Vector<std::unique_ptr<MF> > m_robin_bcval;
// used to save interpolation coefficients of the first interior cells
mutable Vector<Vector<BndryRegisterT<MF>> > m_undrrelxr;
// boundary cell flags for covered, not_covered, outside_domain
Vector<Vector<Array<MultiMask,2*AMREX_SPACEDIM> > > m_maskvals;
Vector<std::unique_ptr<iMultiFab> > m_norm_fine_mask;
mutable Vector<YAFluxRegisterT<MF>> m_fluxreg;
private:
void defineAuxData ();
void defineBC ();
void computeVolInv () const;
mutable Vector<Vector<RT> > m_volinv; // used by solvability fix
};
template <typename T>
struct MLMGABCTag {
Array4<T> fab;
Array4<T const> bcval_lo;
Array4<T const> bcval_hi;
Array4<int const> mask_lo;
Array4<int const> mask_hi;
T bcloc_lo;
T bcloc_hi;
Box bx;
BoundCond bctype_lo;
BoundCond bctype_hi;
int blen;
int comp;
int dir;
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
Box const& box() const noexcept { return bx; }
};
template <typename T>
struct MLMGPSTag {
Array4<T> flo;
Array4<T> fhi;
Array4<int const> mlo;
Array4<int const> mhi;
T bcllo;
T bclhi;
Box bx;
BoundCond bctlo;
BoundCond bcthi;
int blen;
int comp;
int dir;
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
Box const& box() const noexcept { return bx; }
};
#ifdef AMREX_USE_EB
template <typename T>
struct MLMGPSEBTag {
Array4<T> flo;
Array4<T> fhi;
Array4<T const> ap;
Array4<int const> mlo;
Array4<int const> mhi;
T bcllo;
T bclhi;
Box bx;
BoundCond bctlo;
BoundCond bcthi;
int blen;
int comp;
int dir;
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
Box const& box() const noexcept { return bx; }
};
#endif
template <typename MF>
MLCellLinOpT<MF>::BndryCondLoc::BndryCondLoc (const BoxArray& ba,
const DistributionMapping& dm,
int ncomp)
: bcond(ba, dm),
bcloc(ba, dm),
bctl(ba, dm),
bctl_dv(bctl.local_size()*ncomp),
m_ncomp(ncomp)
{
auto* dp = bctl_dv.data();
for (MFIter mfi(bcloc); mfi.isValid(); ++mfi) {
bcond[mfi].resize(ncomp);
bcloc[mfi].resize(ncomp);
bctl[mfi] = dp;
dp += ncomp;
}
}
template <typename MF>
void
MLCellLinOpT<MF>::BndryCondLoc::
setLOBndryConds (const Geometry& geom, const Real* dx,
const Vector<Array<BCType,AMREX_SPACEDIM> >& lobc,
const Vector<Array<BCType,AMREX_SPACEDIM> >& hibc,
int ratio, const RealVect& interior_bloc,
const Array<Real,AMREX_SPACEDIM>& domain_bloc_lo,
const Array<Real,AMREX_SPACEDIM>& domain_bloc_hi)
{
const Box& domain = geom.Domain();
#ifdef AMREX_USE_OMP
#pragma omp parallel
#endif
for (MFIter mfi(bcloc); mfi.isValid(); ++mfi)
{
const Box& bx = mfi.validbox();
for (int icomp = 0; icomp < m_ncomp; ++icomp) {
RealTuple & bloc = bcloc[mfi][icomp];
BCTuple & bctag = bcond[mfi][icomp];
MLMGBndryT<MF>::setBoxBC(bloc, bctag, bx, domain,
lobc[icomp], hibc[icomp],
dx, ratio, interior_bloc,
domain_bloc_lo, domain_bloc_hi,
geom.isPeriodicArray());
}
}
Gpu::PinnedVector<GpuArray<BCTL,2*AMREX_SPACEDIM> > hv;
hv.reserve(bctl_dv.size());
for (MFIter mfi(bctl); mfi.isValid(); ++mfi)
{
for (int icomp = 0; icomp < m_ncomp; ++icomp) {
GpuArray<BCTL,2*AMREX_SPACEDIM> tmp;
for (int m = 0; m < 2*AMREX_SPACEDIM; ++m) {
tmp[m].type = bcond[mfi][icomp][m];
tmp[m].location = bcloc[mfi][icomp][m];
}
hv.push_back(std::move(tmp));
}
}
Gpu::copyAsync(Gpu::hostToDevice, hv.begin(), hv.end(), bctl_dv.begin());
Gpu::streamSynchronize();
}
template <typename MF>
MLCellLinOpT<MF>::MLCellLinOpT ()
{
this->m_ixtype = IntVect::TheCellVector();
}
template <typename MF>
void
MLCellLinOpT<MF>::define (const Vector<Geometry>& a_geom,
const Vector<BoxArray>& a_grids,
const Vector<DistributionMapping>& a_dmap,
const LPInfo& a_info,
const Vector<FabFactory<FAB> const*>& a_factory)
{
MLLinOpT<MF>::define(a_geom, a_grids, a_dmap, a_info, a_factory);
defineAuxData();
defineBC();
}
template <typename MF>
void
MLCellLinOpT<MF>::defineAuxData ()
{
BL_PROFILE("MLCellLinOp::defineAuxData()");
m_undrrelxr.resize(this->m_num_amr_levels);
m_maskvals.resize(this->m_num_amr_levels);
m_fluxreg.resize(this->m_num_amr_levels-1);
m_norm_fine_mask.resize(this->m_num_amr_levels-1);
const int ncomp = this->getNComp();
for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
{
m_undrrelxr[amrlev].resize(this->m_num_mg_levels[amrlev]);
for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
{
m_undrrelxr[amrlev][mglev].define(this->m_grids[amrlev][mglev],
this->m_dmap[amrlev][mglev],
1, 0, 0, ncomp);
}
}
for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
{
m_maskvals[amrlev].resize(this->m_num_mg_levels[amrlev]);
for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
{
for (OrientationIter oitr; oitr; ++oitr)
{
const Orientation face = oitr();
const int ngrow = 1;
const int extent = this->isCrossStencil() ? 0 : 1; // extend to corners
m_maskvals[amrlev][mglev][face].define(this->m_grids[amrlev][mglev],
this->m_dmap[amrlev][mglev],
this->m_geom[amrlev][mglev],
face, 0, ngrow, extent, 1, true);
}
}
}
for (int amrlev = 0; amrlev < this->m_num_amr_levels-1; ++amrlev)
{
const IntVect ratio{this->m_amr_ref_ratio[amrlev]};
m_fluxreg[amrlev].define(this->m_grids[amrlev+1][0],
this->m_grids[amrlev][0],
this->m_dmap[amrlev+1][0],
this->m_dmap[amrlev][0],
this->m_geom[amrlev+1][0],
this->m_geom[amrlev][0],
ratio, amrlev+1, ncomp);
m_norm_fine_mask[amrlev] = std::make_unique<iMultiFab>
(makeFineMask(this->m_grids[amrlev][0], this->m_dmap[amrlev][0],
this->m_grids[amrlev+1][0],
ratio, 1, 0));
}
#if (AMREX_SPACEDIM != 3)
m_has_metric_term = !this->m_geom[0][0].IsCartesian() && this->info.has_metric_term;
#endif
}
template <typename MF>
void
MLCellLinOpT<MF>::defineBC ()
{
BL_PROFILE("MLCellLinOp::defineBC()");
const int ncomp = this->getNComp();
m_bndry_sol.resize(this->m_num_amr_levels);
m_crse_sol_br.resize(this->m_num_amr_levels);
m_bndry_cor.resize(this->m_num_amr_levels);
m_crse_cor_br.resize(this->m_num_amr_levels);
m_robin_bcval.resize(this->m_num_amr_levels);
for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
{
m_bndry_sol[amrlev] = std::make_unique<MLMGBndryT<MF>>(this->m_grids[amrlev][0],
this->m_dmap[amrlev][0],
ncomp,
this->m_geom[amrlev][0]);
}
for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
{
const int in_rad = 0;
const int out_rad = 1;
const int extent_rad = 2;
const int crse_ratio = this->m_amr_ref_ratio[amrlev-1];
BoxArray cba = this->m_grids[amrlev][0];
cba.coarsen(crse_ratio);
m_crse_sol_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
(cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
}
for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
{
const int in_rad = 0;
const int out_rad = 1;
const int extent_rad = 2;
const int crse_ratio = this->m_amr_ref_ratio[amrlev-1];
BoxArray cba = this->m_grids[amrlev][0];
cba.coarsen(crse_ratio);
m_crse_cor_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
(cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
m_crse_cor_br[amrlev]->setVal(RT(0.0));
}
// This has be to done after m_crse_cor_br is defined.
for (int amrlev = 1; amrlev < this->m_num_amr_levels; ++amrlev)
{
m_bndry_cor[amrlev] = std::make_unique<MLMGBndryT<MF>>
(this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, this->m_geom[amrlev][0]);
MF bc_data(this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, 1);
bc_data.setVal(0.0);
m_bndry_cor[amrlev]->setBndryValues(*m_crse_cor_br[amrlev], 0, bc_data, 0, 0, ncomp,
IntVect(this->m_amr_ref_ratio[amrlev-1]));
Vector<Array<LinOpBCType,AMREX_SPACEDIM> > bclohi
(ncomp,Array<LinOpBCType,AMREX_SPACEDIM>{{AMREX_D_DECL(BCType::Dirichlet,
BCType::Dirichlet,
BCType::Dirichlet)}});
m_bndry_cor[amrlev]->setLOBndryConds(bclohi, bclohi, this->m_amr_ref_ratio[amrlev-1], RealVect{});
}
m_bcondloc.resize(this->m_num_amr_levels);
for (int amrlev = 0; amrlev < this->m_num_amr_levels; ++amrlev)
{
m_bcondloc[amrlev].resize(this->m_num_mg_levels[amrlev]);
for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
{
m_bcondloc[amrlev][mglev] = std::make_unique<BndryCondLoc>(this->m_grids[amrlev][mglev],
this->m_dmap[amrlev][mglev],
ncomp);
}
}
}
template <typename MF>
void
MLCellLinOpT<MF>::setLevelBC (int amrlev, const MF* a_levelbcdata, const MF* robinbc_a,
const MF* robinbc_b, const MF* robinbc_f)
{
BL_PROFILE("MLCellLinOp::setLevelBC()");
AMREX_ALWAYS_ASSERT(amrlev >= 0 && amrlev < this->m_num_amr_levels);
const int ncomp = this->getNComp();
MF zero;
IntVect ng(1);
if (this->hasHiddenDimension()) { ng[this->hiddenDirection()] = 0; }
if (a_levelbcdata == nullptr) {
zero.define(this->m_grids[amrlev][0], this->m_dmap[amrlev][0], ncomp, ng);
zero.setVal(RT(0.0));
} else {
AMREX_ALWAYS_ASSERT(a_levelbcdata->nGrowVect().allGE(ng));
}
const MF& bcdata = (a_levelbcdata == nullptr) ? zero : *a_levelbcdata;
int br_ref_ratio = -1;
if (amrlev == 0)
{
if (this->needsCoarseDataForBC())
{
AMREX_ALWAYS_ASSERT(!this->hasHiddenDimension());
br_ref_ratio = this->m_coarse_data_crse_ratio > 0 ? this->m_coarse_data_crse_ratio : 2;
if (m_crse_sol_br[amrlev] == nullptr && br_ref_ratio > 0)
{
const int in_rad = 0;
const int out_rad = 1;
const int extent_rad = 2;
const int crse_ratio = br_ref_ratio;
BoxArray cba = this->m_grids[amrlev][0];
cba.coarsen(crse_ratio);
m_crse_sol_br[amrlev] = std::make_unique<BndryRegisterT<MF>>
(cba, this->m_dmap[amrlev][0], in_rad, out_rad, extent_rad, ncomp);
}
if (this->m_coarse_data_for_bc != nullptr) {
AMREX_ALWAYS_ASSERT(this->m_coarse_data_crse_ratio > 0);
const Box& cbx = amrex::coarsen(this->m_geom[0][0].Domain(), this->m_coarse_data_crse_ratio);
m_crse_sol_br[amrlev]->copyFrom(*this->m_coarse_data_for_bc, 0, 0, 0, ncomp,
this->m_geom[0][0].periodicity(cbx));
} else {
m_crse_sol_br[amrlev]->setVal(RT(0.0));
}
m_bndry_sol[amrlev]->setBndryValues(*m_crse_sol_br[amrlev], 0,
bcdata, 0, 0, ncomp, IntVect(br_ref_ratio));
br_ref_ratio = this->m_coarse_data_crse_ratio;
}
else
{
m_bndry_sol[amrlev]->setPhysBndryValues(bcdata,0,0,ncomp);
br_ref_ratio = 1;
}
}
else
{
m_bndry_sol[amrlev]->setPhysBndryValues(bcdata,0,0,ncomp);
br_ref_ratio = this->m_amr_ref_ratio[amrlev-1];
}
m_bndry_sol[amrlev]->setLOBndryConds(this->m_lobc, this->m_hibc, br_ref_ratio, this->m_coarse_bc_loc);
const Real* dx = this->m_geom[amrlev][0].CellSize();
for (int mglev = 0; mglev < this->m_num_mg_levels[amrlev]; ++mglev)
{
m_bcondloc[amrlev][mglev]->setLOBndryConds(this->m_geom[amrlev][mglev], dx,
this->m_lobc, this->m_hibc,
br_ref_ratio, this->m_coarse_bc_loc,
this->m_domain_bloc_lo, this->m_domain_bloc_hi);
}
if (this->hasRobinBC()) {
AMREX_ASSERT(robinbc_a != nullptr && robinbc_b != nullptr && robinbc_f != nullptr);
m_robin_bcval[amrlev] = std::make_unique<MF>(this->m_grids[amrlev][0],
this->m_dmap[amrlev][0],
ncomp*3, 1);
const Box& domain = this->m_geom[amrlev][0].Domain();
MFItInfo mfi_info;
if (Gpu::notInLaunchRegion()) { mfi_info.SetDynamic(true); }
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(*m_robin_bcval[amrlev], mfi_info); mfi.isValid(); ++mfi) {
Box const& vbx = mfi.validbox();
Array4<RT const> const& ra = robinbc_a->const_array(mfi);
Array4<RT const> const& rb = robinbc_b->const_array(mfi);
Array4<RT const> const& rf = robinbc_f->const_array(mfi);
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
const Box& blo = amrex::adjCellLo(vbx, idim);
const Box& bhi = amrex::adjCellHi(vbx, idim);
bool outside_domain_lo = !(domain.contains(blo));
bool outside_domain_hi = !(domain.contains(bhi));
if ((!outside_domain_lo) && (!outside_domain_hi)) { continue; }
for (int icomp = 0; icomp < ncomp; ++icomp) {
Array4<RT> const& rbc = (*m_robin_bcval[amrlev])[mfi].array(icomp*3);
if (this->m_lobc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_lo)
{
AMREX_HOST_DEVICE_PARALLEL_FOR_3D(blo, i, j, k,
{
rbc(i,j,k,0) = ra(i,j,k,icomp);
rbc(i,j,k,1) = rb(i,j,k,icomp);
rbc(i,j,k,2) = rf(i,j,k,icomp);
});
}
if (this->m_hibc_orig[icomp][idim] == LinOpBCType::Robin && outside_domain_hi)
{
AMREX_HOST_DEVICE_PARALLEL_FOR_3D(bhi, i, j, k,
{
rbc(i,j,k,0) = ra(i,j,k,icomp);
rbc(i,j,k,1) = rb(i,j,k,icomp);
rbc(i,j,k,2) = rf(i,j,k,icomp);
});
}
}
}
}
}
}
template <typename MF>
void
MLCellLinOpT<MF>::update ()
{
if (MLLinOpT<MF>::needsUpdate()) { MLLinOpT<MF>::update(); }
}
template <typename MF>
void
MLCellLinOpT<MF>::updateSolBC (int amrlev, const MF& crse_bcdata) const
{
BL_PROFILE("MLCellLinOp::updateSolBC()");
AMREX_ALWAYS_ASSERT(amrlev > 0);
const int ncomp = this->getNComp();
m_crse_sol_br[amrlev]->copyFrom(crse_bcdata, 0, 0, 0, ncomp,
this->m_geom[amrlev-1][0].periodicity());
m_bndry_sol[amrlev]->updateBndryValues(*m_crse_sol_br[amrlev], 0, 0, ncomp,
IntVect(this->m_amr_ref_ratio[amrlev-1]));
}
template <typename MF>
void
MLCellLinOpT<MF>::updateCorBC (int amrlev, const MF& crse_bcdata) const
{
BL_PROFILE("MLCellLinOp::updateCorBC()");
AMREX_ALWAYS_ASSERT(amrlev > 0);
const int ncomp = this->getNComp();
m_crse_cor_br[amrlev]->copyFrom(crse_bcdata, 0, 0, 0, ncomp,
this->m_geom[amrlev-1][0].periodicity());
m_bndry_cor[amrlev]->updateBndryValues(*m_crse_cor_br[amrlev], 0, 0, ncomp,
IntVect(this->m_amr_ref_ratio[amrlev-1]));
}
template <typename MF>
void
MLCellLinOpT<MF>::applyBC (int amrlev, int mglev, MF& in, BCMode bc_mode, StateMode,
const MLMGBndryT<MF>* bndry, bool skip_fillboundary) const
{
BL_PROFILE("MLCellLinOp::applyBC()");
// No coarsened boundary values, cannot apply inhomog at mglev>0.
BL_ASSERT(mglev == 0 || bc_mode == BCMode::Homogeneous);
BL_ASSERT(bndry != nullptr || bc_mode == BCMode::Homogeneous);
const int ncomp = this->getNComp();
const int cross = isCrossStencil();
const int tensorop = isTensorOp();
if (!skip_fillboundary) {
in.FillBoundary(0, ncomp, this->m_geom[amrlev][mglev].periodicity(), cross);
}
int flagbc = bc_mode == BCMode::Inhomogeneous;
const int imaxorder = this->maxorder;
const Real* dxinv = this->m_geom[amrlev][mglev].InvCellSize();
const RT dxi = static_cast<RT>(dxinv[0]);
const RT dyi = (AMREX_SPACEDIM >= 2) ? static_cast<RT>(dxinv[1]) : RT(1.0);
const RT dzi = (AMREX_SPACEDIM == 3) ? static_cast<RT>(dxinv[2]) : RT(1.0);
const auto& maskvals = m_maskvals[amrlev][mglev];
const auto& bcondloc = *m_bcondloc[amrlev][mglev];
FAB foofab(Box::TheUnitBox(),ncomp);
const auto& foo = foofab.const_array();
MFItInfo mfi_info;
if (Gpu::notInLaunchRegion()) { mfi_info.SetDynamic(true); }
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(cross || tensorop || Gpu::notInLaunchRegion(),
"non-cross stencil not support for gpu");
const int hidden_direction = this->hiddenDirection();
#ifdef AMREX_USE_GPU
if ((cross || tensorop) && Gpu::inLaunchRegion())
{
Vector<MLMGABCTag<RT>> tags;
tags.reserve(in.local_size()*AMREX_SPACEDIM*ncomp);
for (MFIter mfi(in); mfi.isValid(); ++mfi) {
const Box& vbx = mfi.validbox();
const auto& iofab = in.array(mfi);
const auto & bdlv = bcondloc.bndryLocs(mfi);
const auto & bdcv = bcondloc.bndryConds(mfi);
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
if (idim != hidden_direction) {
const Orientation olo(idim,Orientation::low);
const Orientation ohi(idim,Orientation::high);
const auto& bvlo = (bndry != nullptr) ?
bndry->bndryValues(olo).const_array(mfi) : foo;
const auto& bvhi = (bndry != nullptr) ?
bndry->bndryValues(ohi).const_array(mfi) : foo;
for (int icomp = 0; icomp < ncomp; ++icomp) {
tags.emplace_back(MLMGABCTag<RT>{iofab, bvlo, bvhi,
maskvals[olo].const_array(mfi),
maskvals[ohi].const_array(mfi),
bdlv[icomp][olo], bdlv[icomp][ohi],
amrex::adjCell(vbx,olo),
bdcv[icomp][olo], bdcv[icomp][ohi],
vbx.length(idim), icomp, idim});
}
}
}
}
ParallelFor(tags,
[=] AMREX_GPU_DEVICE (int i, int j, int k, MLMGABCTag<RT> const& tag) noexcept
{
if (tag.dir == 0)
{
mllinop_apply_bc_x(0, i, j, k, tag.blen, tag.fab,
tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
imaxorder, dxi, flagbc, tag.comp);
mllinop_apply_bc_x(1, i+tag.blen+1, j, k, tag.blen, tag.fab,
tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
imaxorder, dxi, flagbc, tag.comp);
}
#if (AMREX_SPACEDIM > 1)
else
#if (AMREX_SPACEDIM > 2)
if (tag.dir == 1)
#endif
{
mllinop_apply_bc_y(0, i, j, k, tag.blen, tag.fab,
tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
imaxorder, dyi, flagbc, tag.comp);
mllinop_apply_bc_y(1, i, j+tag.blen+1, k, tag.blen, tag.fab,
tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
imaxorder, dyi, flagbc, tag.comp);
}
#if (AMREX_SPACEDIM > 2)
else {
mllinop_apply_bc_z(0, i, j, k, tag.blen, tag.fab,
tag.mask_lo, tag.bctype_lo, tag.bcloc_lo, tag.bcval_lo,
imaxorder, dzi, flagbc, tag.comp);
mllinop_apply_bc_z(1, i, j, k+tag.blen+1, tag.blen, tag.fab,
tag.mask_hi, tag.bctype_hi, tag.bcloc_hi, tag.bcval_hi,
imaxorder, dzi, flagbc, tag.comp);
}
#endif
#endif
});
} else
#endif
if (cross || tensorop)
{
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(in, mfi_info); mfi.isValid(); ++mfi)
{
const Box& vbx = mfi.validbox();
const auto& iofab = in.array(mfi);
const auto & bdlv = bcondloc.bndryLocs(mfi);
const auto & bdcv = bcondloc.bndryConds(mfi);
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
{
if (hidden_direction == idim) { continue; }
const Orientation olo(idim,Orientation::low);
const Orientation ohi(idim,Orientation::high);
const Box blo = amrex::adjCellLo(vbx, idim);
const Box bhi = amrex::adjCellHi(vbx, idim);
const int blen = vbx.length(idim);
const auto& mlo = maskvals[olo].array(mfi);
const auto& mhi = maskvals[ohi].array(mfi);
const auto& bvlo = (bndry != nullptr) ? bndry->bndryValues(olo).const_array(mfi) : foo;
const auto& bvhi = (bndry != nullptr) ? bndry->bndryValues(ohi).const_array(mfi) : foo;
for (int icomp = 0; icomp < ncomp; ++icomp) {
const BoundCond bctlo = bdcv[icomp][olo];
const BoundCond bcthi = bdcv[icomp][ohi];
const RT bcllo = bdlv[icomp][olo];
const RT bclhi = bdlv[icomp][ohi];
if (idim == 0) {
mllinop_apply_bc_x(0, blo, blen, iofab, mlo,
bctlo, bcllo, bvlo,
imaxorder, dxi, flagbc, icomp);
mllinop_apply_bc_x(1, bhi, blen, iofab, mhi,
bcthi, bclhi, bvhi,
imaxorder, dxi, flagbc, icomp);
} else if (idim == 1) {
mllinop_apply_bc_y(0, blo, blen, iofab, mlo,
bctlo, bcllo, bvlo,
imaxorder, dyi, flagbc, icomp);
mllinop_apply_bc_y(1, bhi, blen, iofab, mhi,
bcthi, bclhi, bvhi,
imaxorder, dyi, flagbc, icomp);
} else {
mllinop_apply_bc_z(0, blo, blen, iofab, mlo,
bctlo, bcllo, bvlo,
imaxorder, dzi, flagbc, icomp);
mllinop_apply_bc_z(1, bhi, blen, iofab, mhi,
bcthi, bclhi, bvhi,
imaxorder, dzi, flagbc, icomp);
}
}
}
}
}
else
{
#ifdef BL_NO_FORT
amrex::Abort("amrex_mllinop_apply_bc not available when BL_NO_FORT=TRUE");
#else
if constexpr (std::is_same_v<Real,RT>) {
#ifdef AMREX_USE_OMP
#pragma omp parallel
#endif
for (MFIter mfi(in, mfi_info); mfi.isValid(); ++mfi)
{
const Box& vbx = mfi.validbox();
const auto & bdlv = bcondloc.bndryLocs(mfi);
const auto & bdcv = bcondloc.bndryConds(mfi);
const RealTuple & bdl = bdlv[0];
const BCTuple & bdc = bdcv[0];
for (OrientationIter oitr; oitr; ++oitr)
{
const Orientation ori = oitr();
int cdr = ori;
RT bcl = bdl[ori];
int bct = bdc[ori];
const auto& fsfab = (bndry != nullptr) ? bndry->bndryValues(ori)[mfi] : foofab;
const Mask& m = maskvals[ori][mfi];
amrex_mllinop_apply_bc(BL_TO_FORTRAN_BOX(vbx),
BL_TO_FORTRAN_ANYD(in[mfi]),
BL_TO_FORTRAN_ANYD(m),
cdr, bct, bcl,
BL_TO_FORTRAN_ANYD(fsfab),
imaxorder, dxinv, flagbc, ncomp, cross);
}
}
} else {
amrex::Abort("Not supported");
}
#endif
}
}
template <typename MF>
BoxArray
MLCellLinOpT<MF>::makeNGrids (int grid_size) const
{
const Box& dombx = this->m_geom[0].back().Domain();
const BoxArray& old_ba = this->m_grids[0].back();
const int N = old_ba.size();
Vector<Box> bv;
bv.reserve(N);
for (int i = 0; i < N; ++i)
{
Box b = old_ba[i];
b.coarsen(grid_size);
b.refine(grid_size);
IntVect sz = b.size();
const IntVect nblks {AMREX_D_DECL(sz[0]/grid_size, sz[1]/grid_size, sz[2]/grid_size)};
IntVect big = b.smallEnd() + grid_size - 1;
b.setBig(big);
#if (AMREX_SPACEDIM == 3)
for (int kk = 0; kk < nblks[2]; ++kk) {
#endif
#if (AMREX_SPACEDIM >= 2)
for (int jj = 0; jj < nblks[1]; ++jj) {
#endif
for (int ii = 0; ii < nblks[0]; ++ii)
{
IntVect shft{AMREX_D_DECL(ii*grid_size,jj*grid_size,kk*grid_size)};
Box bb = amrex::shift(b,shft);
bb &= dombx;
bv.push_back(bb);
}
#if (AMREX_SPACEDIM >= 2)
}
#endif
#if (AMREX_SPACEDIM == 3)
}
#endif
}
std::sort(bv.begin(), bv.end());
bv.erase(std::unique(bv.begin(), bv.end()), bv.end());
BoxList bl(std::move(bv));
return BoxArray{std::move(bl)};
}
template <typename MF>
void
MLCellLinOpT<MF>::restriction (int amrlev, int cmglev, MF& crse, MF& fine) const
{
const int ncomp = this->getNComp();
IntVect ratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[cmglev-1];
amrex::average_down(fine, crse, 0, ncomp, ratio);
}
template <typename MF>
void
MLCellLinOpT<MF>::interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const
{
const int ncomp = this->getNComp();
Dim3 ratio3 = {2,2,2};
IntVect ratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[fmglev];
AMREX_D_TERM(ratio3.x = ratio[0];,
ratio3.y = ratio[1];,
ratio3.z = ratio[2];);
#ifdef AMREX_USE_GPU
if (Gpu::inLaunchRegion() && fine.isFusingCandidate()) {
auto const& finema = fine.arrays();
auto const& crsema = crse.const_arrays();
ParallelFor(fine, IntVect(0), ncomp,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
{
int ic = amrex::coarsen(i,ratio3.x);
int jc = amrex::coarsen(j,ratio3.y);
int kc = amrex::coarsen(k,ratio3.z);
finema[box_no](i,j,k,n) += crsema[box_no](ic,jc,kc,n);
});
Gpu::streamSynchronize();
} else
#endif
{
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(fine,TilingIfNotGPU()); mfi.isValid(); ++mfi)
{
const Box& bx = mfi.tilebox();
Array4<RT const> const& cfab = crse.const_array(mfi);
Array4<RT> const& ffab = fine.array(mfi);
AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
{
int ic = amrex::coarsen(i,ratio3.x);
int jc = amrex::coarsen(j,ratio3.y);
int kc = amrex::coarsen(k,ratio3.z);
ffab(i,j,k,n) += cfab(ic,jc,kc,n);
});
}
}
}
template <typename MF>
void
MLCellLinOpT<MF>::interpAssign (int amrlev, int fmglev, MF& fine, MF& crse) const
{
const int ncomp = this->getNComp();
const Geometry& crse_geom = this->Geom(amrlev,fmglev+1);
const IntVect refratio = (amrlev > 0) ? IntVect(2) : this->mg_coarsen_ratio_vec[fmglev];
const IntVect ng = crse.nGrowVect();
MF cfine;
const MF* cmf;
if (amrex::isMFIterSafe(crse, fine))
{
crse.FillBoundary(crse_geom.periodicity());
cmf = &crse;
}
else
{
BoxArray cba = fine.boxArray();
cba.coarsen(refratio);
cfine.define(cba, fine.DistributionMap(), ncomp, ng);
cfine.setVal(RT(0.0));
cfine.ParallelCopy(crse, 0, 0, ncomp, IntVect(0), ng, crse_geom.periodicity());
cmf = & cfine;
}
bool isEB = fine.hasEBFabFactory();
ignore_unused(isEB);