chombo-discharge
Loading...
Searching...
No Matches
CD_EBHelmholtzOp.H
Go to the documentation of this file.
1/* chombo-discharge
2 * Copyright © 2021 SINTEF Energy Research.
3 * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4 */
5
12#ifndef CD_EBHelmholtzOp_H
13#define CD_EBHelmholtzOp_H
14
15// Std includes
16#include <map>
17
18// Chombo includes
19#include <BaseEBBC.H>
20#include <LevelTGA.H>
21#include <BaseIVFAB.H>
22#include <VCAggStencil.H>
23
24// Our includes
25#include <CD_Timer.H>
26#include <CD_Location.H>
28#include <CD_EBCoarAve.H>
29#include <CD_EBReflux.H>
30#include <CD_EBMGRestrict.H>
31#include <CD_EBMGProlong.H>
32#include <CD_EBHelmholtzEBBC.H>
33#include <CD_EBHelmholtzDomainBC.H>
34#include <CD_NamespaceHeader.H>
35
40class EBHelmholtzOp : public LevelTGAHelmOp<LevelData<EBCellFAB>, EBFluxFAB>
41{
42public:
46 enum class Smoother
47 {
48 NoRelax,
49 PointJacobi,
50 GauSaiRedBlack,
51 GauSaiMultiColor,
52 };
53
57 EBHelmholtzOp() = delete;
58
64
70
103 const EBLevelGrid& a_eblgFine,
104 const EBLevelGrid& a_eblg,
105 const EBLevelGrid& a_eblgCoFi,
106 const EBLevelGrid& a_eblgCoar,
114 const RealVect& a_probLo,
115 const Real& a_dx,
116 const int& a_refToFine,
117 const int& a_refToCoar,
118 const bool& a_hasFine,
119 const bool& a_hasCoar,
120 const bool& a_hasMGObjects,
121 const Real& a_alpha,
122 const Real& a_beta,
128 const Smoother& a_smoother,
129 const Real& a_relaxFactor);
130
134 virtual ~EBHelmholtzOp();
135
142
148 operator=(const EBHelmholtzOp&& a_oper) = delete;
149
153 void
155
159 void
161
166 void
168
173 void
175
180 void
182
187 void
189
197 void
201
207 getAcoef();
208
214 getBcoef();
215
222
229 void
231
239 void
241
253 void
256 const EBCellFAB& a_resid,
257 const EBCellFAB& a_Acoef,
258 const EBFluxFAB& a_Bcoef,
260 const Box& a_cellBox,
261 const DataIndex& a_dit) const noexcept;
262
275 void
278 const EBCellFAB& a_resid,
279 const EBCellFAB& a_Acoef,
280 const EBFluxFAB& a_Bcoef,
282 const Box& a_cellBox,
283 const DataIndex& a_dit,
284 const int& a_redBlack) const noexcept;
285
298 void
301 const EBCellFAB& a_resid,
302 const EBCellFAB& a_Acoef,
303 const EBFluxFAB& a_Bcoef,
305 const Box& a_cellBox,
306 const DataIndex& a_dit,
307 const IntVect& a_color) const noexcept;
308
316 void
320 const bool a_homogeneousPhysBc) override;
321
328 void
330
337 void
339
344 void
346
352 void
354
362 void
366 const IntVect& a_color) const;
367
377 void
380 const LevelData<EBCellFAB>* const a_phiCoar,
381 const bool a_homogeneousPhysBC,
382 const bool a_homogeneousCFBC);
383
391 void
393
402 void
405 const EBCellFAB& a_Acoef,
406 const EBFluxFAB& a_Bcoef,
408 const Box& a_cellBox,
409 const DataIndex& a_dit,
410 const bool a_homogeneousPhysBC) const noexcept;
411
420 void
423 const EBCellFAB& a_Acoef,
424 const EBFluxFAB& a_Bcoef,
426 const Box& a_cellBox,
427 const DataIndex& a_dit,
428 const bool a_homogeneousPhysBC) const noexcept;
429
440 void
442 const EBFluxFAB& a_Bcoef,
443 const Box& a_cellBox,
444 const DataIndex& a_dit,
445 const bool a_homogeneousPhysBc) const noexcept;
446
455 void
457
466 void
468 const EBCellFAB& a_phi,
469 const EBCellFAB& a_Acoef,
470 const EBFluxFAB& a_Bcoef,
473 const Box& a_cellBox,
474 const DataIndex& a_dit,
475 const bool a_homogeneousPhysBC) const noexcept;
476
482 void
484
489 void
491
500 void
502
508 void
510
517 void
519
523 Real
524 dotProduct(const LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override final;
525
532 void
533 incr(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs, const Real a_scale) override final;
534
543 void
547 const Real a_a,
548 const Real a_b) override final;
549
555 void
556 scale(LevelData<EBCellFAB>& a_lhs, const Real& a_scale) override final;
557
563 void
565
571 Real
572 norm(const LevelData<EBCellFAB>& a_rhs, const int a_order) override final;
573
578 void
579 setToZero(LevelData<EBCellFAB>& a_lhs) override final;
580
587 void
589
596 void
598
605 void
608 const LevelData<EBCellFAB>& a_rhs) override final;
609
615 void
617
621 int
623
634 void
641
653 void
660
671 void
679
688 void
694
704 void
711
720 void
725
734 void
740
749 void
755
761 void
763
770 void
774
780 void
782
789 void
791
797 void
799
804 void
806
810 void
812
817 void
819
828 void
834
838 void
840
844 void
846
851 getFlux() const;
852
858
864
869
874
879
884
889
894
899
904
909
914
919
924
929
935
941
947
952
957
962
967
972
977
982
987
992
997
1002
1007
1012
1017
1022
1027
1032
1037
1042
1047
1052
1057
1062
1067
1072
1077
1082
1087
1092
1097
1102
1107
1112
1119
1125
1130
1135
1140
1145
1150
1155
1160
1165
1170
1177 void
1179
1186 void
1188
1195 void
1197
1201 void
1203
1207 void
1209
1213 void
1215
1219 void
1221
1231
1241
1250 void
1255 const int a_dir);
1256
1265 void
1270 const int a_dir);
1271
1280 void
1285 const int a_dir);
1286
1292 void
1294
1302 void
1307};
1308
1309#include <CD_NamespaceFooter.H>
1310
1311#endif
Declaration of conservative coarsening utility.
Declaration of a prolongation operator for EB geometric multigrid.
Declaration of a restriction operator for EB geometric multigrid.
Declaration of a class that can interpolate more ghost cells near the coarse-fine boundary near the E...
Declaration of a class which can reflux over the coarse-fine interface.
Declaration of cell positions.
Implementation of CD_Timer.H.
Class which replaces data at coarse level of refinement with average at fine level of refinement.
Definition CD_EBCoarAve.H:31
Base class for passing domain boundary conditions into EBHelmholtzOp.
Definition CD_EBHelmholtzDomainBC.H:29
Base class for passing EB boundary conditions into EBHelmholtzOp.
Definition CD_EBHelmholtzEBBC.H:29
Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho.
Definition CD_EBHelmholtzOp.H:41
VoFStencil getFaceCentroidFluxStencil(const FaceIndex &a_face, const DataIndex &a_dit) const
Get the face-centroid flux stencil.
Definition CD_EBHelmholtzOp.cpp:2127
void relaxGSMultiColor(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Multi-colored gauss-seidel relaxation.
Definition CD_EBHelmholtzOp.cpp:1829
void setAcoAndBco(const RefCountedPtr< LevelData< EBCellFAB > > &a_Acoef, const RefCountedPtr< LevelData< EBFluxFAB > > &a_Bcoef, const RefCountedPtr< LevelData< BaseIVFAB< Real > > > &a_BcoefIrreg)
Update with new A and B coefficients.
Definition CD_EBHelmholtzOp.cpp:165
void reflux(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > &a_finerOp)
Reflux algorithm.
Definition CD_EBHelmholtzOp.cpp:2304
LevelData< EBFluxFAB > & getFlux() const
Returns m_flux. This is used in refluxing routines.
Definition CD_EBHelmholtzOp.cpp:268
static constexpr int m_nComp
Number of components that we solve for (always one..)
Definition CD_EBHelmholtzOp.H:863
void AMRResidualNC(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp) override final
Compute AMR residual on coarsest.
Definition CD_EBHelmholtzOp.cpp:1069
RefCountedPtr< EBCoarAve > m_coarAve
Conservative coarsener.
Definition CD_EBHelmholtzOp.H:1071
void axby(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, const Real a_a, const Real a_b) override final
Set a_lhs = a*x + b*y.
Definition CD_EBHelmholtzOp.cpp:643
void create(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Create data which clones the layout of the other.
Definition CD_EBHelmholtzOp.cpp:536
void gauSaiMultiColor(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color) const
Multi-colored Gauss-Seidel kernel. Public because MFHelmholtzOp may want to use use.
static constexpr int m_comp
Component that we solve for.
Definition CD_EBHelmholtzOp.H:868
void turnOnCFInterp()
Turn on BCs.
Definition CD_EBHelmholtzOp.cpp:212
EBLevelGrid m_eblgCoFi
Coarsened of m_eblg.
Definition CD_EBHelmholtzOp.H:1011
const RefCountedPtr< LevelData< BaseIVFAB< Real > > > & getBcoefIrreg()
Get the Helmholtz B-coefficient on the EB.
Definition CD_EBHelmholtzOp.cpp:193
void AMRProlong(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection) override final
Prolongation onto AMR level.
Definition CD_EBHelmholtzOp.cpp:1114
Smoother m_smoother
Relaxation method.
Definition CD_EBHelmholtzOp.H:883
void applyOpRegular(EBCellFAB &a_Lphi, EBCellFAB &a_phi, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBC) const noexcept
Apply operator in regular cells.
Definition CD_EBHelmholtzOp.cpp:1224
LayoutData< BaseIVFAB< VoFStencil > > m_relaxStencils
Operator stencils in irregular cells (and ones that border irregular cells if using a centroid discre...
Definition CD_EBHelmholtzOp.H:1118
bool m_hasCoar
True if there's a coarser level.
Definition CD_EBHelmholtzOp.H:913
void AMRResidualNF(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC) override final
Compute AMR residual on finest AMR level.
Definition CD_EBHelmholtzOp.cpp:1054
void computeFaceCentroidFlux(EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
Compute face-centered fluxes.
Definition CD_EBHelmholtzOp.cpp:2250
void applyDomainFlux(EBCellFAB &a_phi, const EBFluxFAB &a_Bcoef, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBc) const noexcept
Apply domain flux.
Definition CD_EBHelmholtzOp.cpp:1278
void applyOp(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoar, const bool a_homogeneousPhysBC, const bool a_homogeneousCFBC)
Apply operator on this level. This is a more general version which can turn on/off homogeneous and CF...
Definition CD_EBHelmholtzOp.cpp:1150
void allocateFlux() const noexcept
Allocate m_flux.
Definition CD_EBHelmholtzOp.cpp:252
void computeDiagWeight()
Calculate the weight of the diagonal term.
Definition CD_EBHelmholtzOp.cpp:1942
void scale(LevelData< EBCellFAB > &a_lhs, const Real &a_scale) override final
Scale data. Returns a_lhs = a_lhs*a_scale.
Definition CD_EBHelmholtzOp.cpp:655
EBHelmholtzOp()=delete
Disallowed default constructor.
virtual ~EBHelmholtzOp()
Dtor.
Definition CD_EBHelmholtzOp.cpp:198
Copier m_exchangeCopierFine
Pre-built exchange copier.
Definition CD_EBHelmholtzOp.H:996
void relax(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, int a_iterations) override final
Relaxation method. This does smoothing for the system L(correction) = residual.
Definition CD_EBHelmholtzOp.cpp:1620
void turnOffCoarsening()
Turn off coarsening operation.
Definition CD_EBHelmholtzOp.cpp:236
void relaxPointJacobi(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition CD_EBHelmholtzOp.cpp:1654
RefCountedPtr< EBHelmholtzEBBC > m_ebBc
Domain bc object.
Definition CD_EBHelmholtzOp.H:1056
void makeAggStencil()
Compute aggregated stencils.
Definition CD_EBHelmholtzOp.cpp:2057
void homogeneousCFInterp(LevelData< EBCellFAB > &a_phi)
Do homogeneous coarse-fine interpolation.
Definition CD_EBHelmholtzOp.cpp:1569
EBLevelGrid m_eblg
Grid.
Definition CD_EBHelmholtzOp.H:1006
LayoutData< RefCountedPtr< VCAggStencil > > m_aggRelaxStencil
For making irregular stencil applications go faster.
Definition CD_EBHelmholtzOp.H:1124
EBMGRestrict m_restrictOp
Restriction operator for AMR levels.
Definition CD_EBHelmholtzOp.H:1026
void pointJacobiKernel(EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit) const noexcept
Point Jacobi kernel.
Definition CD_EBHelmholtzOp.cpp:1694
bool m_refluxFree
Use reflux-free formulation or not.
Definition CD_EBHelmholtzOp.H:893
RefCountedPtr< EBHelmholtzDomainBC > m_domainBc
Domain bc object.
Definition CD_EBHelmholtzOp.H:1051
Real dotProduct(const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Compute the dot product??
Definition CD_EBHelmholtzOp.cpp:568
RefCountedPtr< EBReflux > m_fluxReg
Flux register.
Definition CD_EBHelmholtzOp.H:1066
void preCond(LevelData< EBCellFAB > &a_corr, const LevelData< EBCellFAB > &a_residual) override final
Precondition system before bottom solve.
Definition CD_EBHelmholtzOp.cpp:525
int m_numSmoothPreCond
Number of smoothings in the preconditioner.
Definition CD_EBHelmholtzOp.H:928
bool m_doInterpCF
Do coarse-fine interpolation or not.
Definition CD_EBHelmholtzOp.H:934
int m_refToFine
Refinement factor to fine level.
Definition CD_EBHelmholtzOp.H:923
Smoother
Relaxation method for the operators.
Definition CD_EBHelmholtzOp.H:47
LevelData< EBFluxFAB > * m_flux
For holding fluxes.
Definition CD_EBHelmholtzOp.H:1096
bool m_hasMGObjects
True if there is a multigrid level below this operator.
Definition CD_EBHelmholtzOp.H:903
Real m_dx
Grid resolution;.
Definition CD_EBHelmholtzOp.H:971
bool m_profile
Profile the operator.
Definition CD_EBHelmholtzOp.H:898
void turnOffExchange()
Turn off exchange operation.
Definition CD_EBHelmholtzOp.cpp:220
void relaxGSRedBlack(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition CD_EBHelmholtzOp.cpp:1720
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
Weights of diagonal alpha terms.
Definition CD_EBHelmholtzOp.H:1129
void fillGrad(const LevelData< EBCellFAB > &a_phi) override final
Not called, I think.
Definition CD_EBHelmholtzOp.cpp:1549
VoFStencil getFaceCenterFluxStencil(const FaceIndex &a_face, const DataIndex &a_dit) const
Get the face-centered flux stencil.
Definition CD_EBHelmholtzOp.cpp:2108
void restrictResidual(LevelData< EBCellFAB > &a_resCoar, LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs) override final
Restrict residual onto coarse level.
Definition CD_EBHelmholtzOp.cpp:768
void divideByIdentityCoef(LevelData< EBCellFAB > &a_rhs) override final
Divide by the a-coefficient.
Definition CD_EBHelmholtzOp.cpp:1523
void assign(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Assign data.
Definition CD_EBHelmholtzOp.cpp:544
RealVect m_probLo
Lower-left corner of domain.
Definition CD_EBHelmholtzOp.H:986
RefCountedPtr< EBMultigridInterpolator > m_interpolator
Interpolator object.
Definition CD_EBHelmholtzOp.H:1061
void computeRelaxationCoefficient()
Calculate relaxation coefficient.
Definition CD_EBHelmholtzOp.cpp:1997
void computeFaceCenteredFlux(EBFaceFAB &a_fluxCenter, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
Compute face-centered fluxes.
Definition CD_EBHelmholtzOp.cpp:2222
void AMROperatorNF(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, bool a_homogeneousPhysBC) override final
Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no finer levels.
Definition CD_EBHelmholtzOp.cpp:1012
EBMGProlong m_prolongOp
Prolongation operator for AMR levels.
Definition CD_EBHelmholtzOp.H:1036
void refluxFreeAMROperator(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
Apply the AMR operator, i.e. compute L(phi) in an AMR context.
Definition CD_EBHelmholtzOp.cpp:840
void AMROperatorNC(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp) override final
Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels.
Definition CD_EBHelmholtzOp.cpp:1023
IntVect m_ghostRhs
Ghost cells for rhs (note, the operator rhs)
Definition CD_EBHelmholtzOp.H:956
int refToCoarser() override final
Return coarsening factor to coarser level (1 if there is no coarser level);.
Definition CD_EBHelmholtzOp.cpp:793
LayoutData< BaseIFFAB< VoFStencil > > m_centroidFluxStencil[SpaceDim]
Face centroid flux stencil. Defined on all faces connecting one or more irregular vofs.
Definition CD_EBHelmholtzOp.H:1111
EBHelmholtzOp & operator=(const EBHelmholtzOp &a_oper)=delete
No copy assigment allowed.
Interval m_interval
Interval.
Definition CD_EBHelmholtzOp.H:878
void applyOpNoBoundary(LevelData< EBCellFAB > &a_ans, const LevelData< EBCellFAB > &a_phi) override final
Apply operator but turn off all BCs.
Definition CD_EBHelmholtzOp.cpp:1539
void assignLocal(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Local assignment function.
Definition CD_EBHelmholtzOp.cpp:560
bool m_doExchange
Turn on/off exchange operation.
Definition CD_EBHelmholtzOp.H:940
EBMGProlong m_prolongOpMG
Prolongation operator for MG levels.
Definition CD_EBHelmholtzOp.H:1041
Real norm(const LevelData< EBCellFAB > &a_rhs, const int a_order) override final
Compute norm of data.
Definition CD_EBHelmholtzOp.cpp:684
void AMRResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp) override final
Compute residual on this level. AMR version.
Definition CD_EBHelmholtzOp.cpp:1036
LayoutData< BaseIFFAB< FaceStencil > > m_interpStencil[SpaceDim]
Face centroid interpolation stencil.
Definition CD_EBHelmholtzOp.H:1106
void residual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, const bool a_homogeneousPhysBc) override
Compute residual on this level.
Definition CD_EBHelmholtzOp.cpp:512
void fillDomainFlux(EBFluxFAB &a_flux, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit)
Fill domain flux. This fills the flux on the domain face using centered differencing ala applyDomainF...
Definition CD_EBHelmholtzOp.cpp:1365
EBHelmholtzOp(const EBHelmholtzOp &&a_other)=delete
Disallowed move constructor.
RefCountedPtr< LevelData< EBFluxFAB > > m_Bcoef
B-coefficient in Helmholtz equation.
Definition CD_EBHelmholtzOp.H:1081
void deallocateFlux() const noexcept
Deallocate m_flux.
Definition CD_EBHelmholtzOp.cpp:260
Timer m_timer
Timer so user can profile.
Definition CD_EBHelmholtzOp.H:873
Location::Cell m_dataLocation
Data centering.
Definition CD_EBHelmholtzOp.H:888
void createCoarser(LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted) override final
Create coarsened data.
Definition CD_EBHelmholtzOp.cpp:746
EBLevelGrid m_eblgCoar
Coarse level grid (if the operator has a coarse level)
Definition CD_EBHelmholtzOp.H:1016
IntVect m_ghostPhi
Ghost cells for phi.
Definition CD_EBHelmholtzOp.H:951
void computeFlux(EBFaceFAB &a_fluxCentroid, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
Compute centroid fluxes in a direction.
Definition CD_EBHelmholtzOp.cpp:2203
EBLevelGrid m_eblgFine
Fine level grid (if the operator has a fine level)
Definition CD_EBHelmholtzOp.H:1001
void turnOffCFInterp()
Turn off BCs.
Definition CD_EBHelmholtzOp.cpp:204
Real m_relaxFactor
Successive over-relaxation factor (Must be >= 1).
Definition CD_EBHelmholtzOp.H:976
int m_refToCoar
Refinement factor to coarse level.
Definition CD_EBHelmholtzOp.H:918
Real m_alpha
Alpha-coefficient.
Definition CD_EBHelmholtzOp.H:961
LevelData< EBCellFAB > m_relCoef
Relaxation coefficient.
Definition CD_EBHelmholtzOp.H:1091
void diagonalScale(LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted) override final
Do diagonal scaling.
Definition CD_EBHelmholtzOp.cpp:1501
void setToZero(LevelData< EBCellFAB > &a_lhs) override final
Set data to zero.
Definition CD_EBHelmholtzOp.cpp:738
void setAlphaAndBeta(const Real &a_alpha, const Real &a_beta) override final
Set alpha coefficient and beta coefficient (can change as diffusion solvers progress)
Definition CD_EBHelmholtzOp.cpp:498
RefCountedPtr< LevelData< BaseFab< bool > > > m_validCells
Valid grid cells (might be nullptr on MG levels)
Definition CD_EBHelmholtzOp.H:1046
Copier m_exchangeCopier
Pre-built exchange copier.
Definition CD_EBHelmholtzOp.H:991
void interpolateCF(LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > *a_phiCoar, const bool a_homogeneousCFBC)
Apply coarse-fine boundary conditions.
Definition CD_EBHelmholtzOp.cpp:1593
void coarsenCell(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiFine)
Coarsen data from fine to coar level.
Definition CD_EBHelmholtzOp.cpp:2342
void gauSaiRedBlackKernel(EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const int &a_redBlack) const noexcept
Red-black Gauss-Seidel kernel.
Definition CD_EBHelmholtzOp.cpp:1767
LayoutData< BaseIFFAB< Real > > m_interpolant[SpaceDim]
Interpolant for when we want centroid fluxes.
Definition CD_EBHelmholtzOp.H:1101
void turnOnExchange()
Turn on exchange operation.
Definition CD_EBHelmholtzOp.cpp:228
void inhomogeneousCFInterp(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar)
Inhomogeneous coarse-fine interpolation.
Definition CD_EBHelmholtzOp.cpp:1581
void coarsenFlux(LevelData< EBFluxFAB > &a_flux, const LevelData< EBFluxFAB > &a_fineFlux)
Coarsen fluxes on the fine level onto this level.
Definition CD_EBHelmholtzOp.cpp:2350
bool m_hasFine
True if there's a finer level.
Definition CD_EBHelmholtzOp.H:908
bool m_doCoarsen
Turn on/off exchange operation.
Definition CD_EBHelmholtzOp.H:946
void assignCopier(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Copier &a_copier) override final
Assign lhs.
Definition CD_EBHelmholtzOp.cpp:552
void gauSaiMultiColorKernel(EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const IntVect &a_color) const noexcept
Multi-color Gauss-Seidel kernel.
Definition CD_EBHelmholtzOp.cpp:1875
std::map< std::pair< int, Side::LoHiSide >, Box > m_sideBox
Domain boxes on each side.
Definition CD_EBHelmholtzOp.H:1164
LayoutData< VoFIterator > m_vofIterStenc
VoFIterator which iterates over all cells that are 1) a cut-cell or 2) borders a cut-cell.
Definition CD_EBHelmholtzOp.H:1149
void scaleLocal(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) const noexcept
Local scaling function. Multiplies the left-hand side by the right-hand side.
Definition CD_EBHelmholtzOp.cpp:663
const RefCountedPtr< LevelData< EBCellFAB > > & getAcoef()
Get the Helmholtz A-coefficient on cell centers.
Definition CD_EBHelmholtzOp.cpp:181
RealVect m_vecDx
Vector resolution.
Definition CD_EBHelmholtzOp.H:981
void incr(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Real a_scale) override final
Increment operator.
Definition CD_EBHelmholtzOp.cpp:635
LayoutData< VoFIterator > m_vofIterMulti
VoFIterator for "multi-cells".
Definition CD_EBHelmholtzOp.H:1144
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
Weights of diagonal beta terms.
Definition CD_EBHelmholtzOp.H:1134
EBHelmholtzOp(const EBHelmholtzOp &a_other)=delete
Disallowed copy constructor.
Vector< IntVect > m_colors
"Colors" for the relaxation methods
Definition CD_EBHelmholtzOp.H:1169
const RefCountedPtr< LevelData< EBFluxFAB > > & getBcoef()
Get the Helmholtz B-coefficient on faces.
Definition CD_EBHelmholtzOp.cpp:187
RefCountedPtr< LevelData< EBCellFAB > > m_Acoef
A-coefficient in Helmholtz equation.
Definition CD_EBHelmholtzOp.H:1076
void turnOnCoarsening()
Turn on coarsening operation.
Definition CD_EBHelmholtzOp.cpp:244
LayoutData< VoFIterator > m_vofIterDomLo[SpaceDim]
VoF iterators for lo domain side.
Definition CD_EBHelmholtzOp.H:1154
LayoutData< VoFIterator > m_vofIterIrreg
VoFIterator for irregular cells.
Definition CD_EBHelmholtzOp.H:1139
void AMROperator(LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp) override final
Apply the AMR operator, i.e. compute L(phi) in an AMR context.
Definition CD_EBHelmholtzOp.cpp:799
void AMRRestrict(LevelData< EBCellFAB > &a_residualCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection, bool a_skip_res) override final
Restrict residual.
Definition CD_EBHelmholtzOp.cpp:1086
EBMGRestrict m_restrictOpMG
Restriction operator if this is an MG level.
Definition CD_EBHelmholtzOp.H:1031
void createCoarsened(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat) override final
Create coarsening of data holder.
Definition CD_EBHelmholtzOp.cpp:756
Real m_beta
Beta-coefficient.
Definition CD_EBHelmholtzOp.H:966
EBHelmholtzOp & operator=(const EBHelmholtzOp &&a_oper)=delete
No move assigment allowed.
void applyOpIrregular(EBCellFAB &a_Lphi, const EBCellFAB &a_phi, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const BaseIVFAB< Real > &a_alphaDiagWeight, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBC) const noexcept
Apply operator in irregular cells.
Definition CD_EBHelmholtzOp.cpp:1419
void AMRUpdateResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection) override final
Update AMR residual.
Definition CD_EBHelmholtzOp.cpp:1122
LayoutData< VoFIterator > m_vofIterDomHi[SpaceDim]
VoF iterators for hi domain side.
Definition CD_EBHelmholtzOp.H:1159
void buildCopier(Copier &a_copier, const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override
Build copier.
Definition CD_EBHelmholtzOp.cpp:2358
void defineStencils()
Define stencils.
Definition CD_EBHelmholtzOp.cpp:284
void prolongIncrement(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_correctCoarse) override final
Prolongation method.
Definition CD_EBHelmholtzOp.cpp:785
const LevelData< EBCellFAB > & getRelaxationCoeff() const noexcept
Returns m_relCoef. Used by MFHelmholtzOp.
Definition CD_EBHelmholtzOp.cpp:276
EBLevelGrid m_eblgCoarMG
Coarser grids (multigrid level)
Definition CD_EBHelmholtzOp.H:1021
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_BcoefIrreg
B-coefficient in Helmholtz equation, but on EB faces.
Definition CD_EBHelmholtzOp.H:1086
Class for prolongation of multigrid residual onto a finer grid.
Definition CD_EBMGProlong.H:30
Class for restricting multigrid residual onto a coarser grid.
Definition CD_EBMGRestrict.H:30
Multigrid interpolator class.
Definition CD_EBMGLeastSquaresInterpolator.H:48
Class which can do refluxing across a coarse-fine interface.
Definition CD_EBReflux.H:37
Class which is used for run-time monitoring of events.
Definition CD_Timer.H:31
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
Namespace for encapsulating various data centerings.
Definition CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition CD_Location.H:30