chombo-discharge
CD_MFHelmholtzOp.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_MFHelmholtzOp_H
13 #define CD_MFHelmholtzOp_H
14 
15 // Std includes
16 #include <map>
17 #include <chrono>
18 
19 // Chombo includes
20 #include <MFCellFAB.H>
21 #include <MFFluxFAB.H>
22 
23 // Our includes
24 #include <CD_Location.H>
26 #include <CD_MFCoarAve.H>
27 #include <CD_MFReflux.H>
28 #include <CD_MFLevelGrid.H>
29 #include <CD_MFBaseIVFAB.H>
30 #include <CD_MFHelmholtzJumpBC.H>
31 #include <CD_EBHelmholtzOp.H>
32 #include <CD_MFHelmholtzDomainBCFactory.H>
35 #include <CD_NamespaceHeader.H>
36 
42 class MFHelmholtzOp : public LevelTGAHelmOp<LevelData<MFCellFAB>, MFFluxFAB>
43 {
44 public:
48  enum class Smoother
49  {
50  PointJacobi,
51  GauSaiRedBlack,
52  GauSaiMultiColor,
53  };
54 
58  MFHelmholtzOp() = delete;
59 
63  MFHelmholtzOp(const MFHelmholtzOp& a_op) = delete;
64 
68  MFHelmholtzOp(const MFHelmholtzOp&& a_op) = delete;
69 
103  MFHelmholtzOp(const Location::Cell a_dataLocation,
104  const MFLevelGrid& a_mflgFine,
105  const MFLevelGrid& a_mflg,
106  const MFLevelGrid& a_mflgCoFi,
107  const MFLevelGrid& a_mflgCoar,
108  const MFLevelGrid& a_mflgCoarMG,
109  const MFMultigridInterpolator& a_interpolator,
110  const MFReflux& a_fluxReg,
111  const MFCoarAve& a_coarAve,
112  const RefCountedPtr<MFHelmholtzDomainBCFactory>& a_domainBC,
113  const RefCountedPtr<MFHelmholtzEBBCFactory>& a_ebBC,
114  const RefCountedPtr<MFHelmholtzJumpBCFactory>& a_jumpBcFactory,
115  const RealVect& a_probLo,
116  const Real& a_dx,
117  const int& a_refToFine,
118  const int& a_refToCoar,
119  const bool& a_hasFine,
120  const bool& a_hasCoar,
121  const bool& a_hasMGObjects,
122  const bool& a_isMGOperator,
123  const Real& a_alpha,
124  const Real& a_beta,
125  const RefCountedPtr<LevelData<MFCellFAB>>& a_Acoef,
126  const RefCountedPtr<LevelData<MFFluxFAB>>& a_Bcoef,
127  const RefCountedPtr<LevelData<MFBaseIVFAB>>& a_BcoefIrreg,
128  const IntVect& a_ghostPhi,
129  const IntVect& a_ghostRhs,
130  const int& a_jumpOrder,
131  const int& a_jumpWeight,
132  const Smoother& a_smoother);
133 
137  ~MFHelmholtzOp();
138 
143  operator=(const MFHelmholtzOp& a_oper) = delete;
144 
149  operator=(const MFHelmholtzOp&& a_oper) = delete;
150 
157  void
158  setAcoAndBco(const RefCountedPtr<LevelData<MFCellFAB>>& a_Acoef,
159  const RefCountedPtr<LevelData<MFFluxFAB>>& a_Bcoef,
160  const RefCountedPtr<LevelData<MFBaseIVFAB>>& a_BcoefIrreg);
161 
166  const RefCountedPtr<LevelData<MFCellFAB>>&
167  getAcoef();
168 
173  const RefCountedPtr<LevelData<MFFluxFAB>>&
174  getBcoef();
175 
180  const RefCountedPtr<LevelData<MFBaseIVFAB>>&
181  getBcoefIrreg();
182 
186  void
187  setJump(RefCountedPtr<LevelData<BaseIVFAB<Real>>>& a_jump);
188 
192  RefCountedPtr<EBHelmholtzOp>&
193  getHelmOp(const int a_phase);
194 
198  int
199  refToCoarser() override final;
200 
206  void
207  setAlphaAndBeta(const Real& a_alpha, const Real& a_beta) override final;
208 
213  void
214  divideByIdentityCoef(LevelData<MFCellFAB>& a_rhs) override final;
215 
219  void
220  applyOpNoBoundary(LevelData<MFCellFAB>& a_ans, const LevelData<MFCellFAB>& a_phi) override final;
221 
226  void
227  fillGrad(const LevelData<MFCellFAB>& a_phi) override final;
228 
237  void
238  getFlux(MFFluxFAB& a_flux,
239  const LevelData<MFCellFAB>& a_data,
240  const Box& a_grid,
241  const DataIndex& a_dit,
242  Real a_scale) override final;
243 
251  void
252  residual(LevelData<MFCellFAB>& a_residual,
253  const LevelData<MFCellFAB>& a_phi,
254  const LevelData<MFCellFAB>& a_rhs,
255  const bool a_homogeneousPhysBc) override final;
256 
263  void
264  preCond(LevelData<MFCellFAB>& a_corr, const LevelData<MFCellFAB>& a_residual) override final;
265 
273  void
274  applyOp(LevelData<MFCellFAB>& a_Lphi, const LevelData<MFCellFAB>& a_phi, bool a_homogeneousPhysBc) override final;
275 
282  template <typename Duration = std::chrono::microseconds>
283  Vector<long long>
284  computeOperatorLoads(LevelData<MFCellFAB>& a_phi, const int a_numApply) noexcept;
285 
292  void
293  relax(LevelData<MFCellFAB>& a_correction, const LevelData<MFCellFAB>& a_residual, int a_iterations) override final;
294 
301  void
302  relaxPointJacobi(LevelData<MFCellFAB>& a_correction, const LevelData<MFCellFAB>& a_residual, const int a_iterations);
303 
310  void
311  relaxGSRedBlack(LevelData<MFCellFAB>& a_correction, const LevelData<MFCellFAB>& a_residual, const int a_iterations);
312 
319  void
320  relaxGSMultiColor(LevelData<MFCellFAB>& a_correction, const LevelData<MFCellFAB>& a_residual, const int a_iterations);
321 
327  void
328  create(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs) override final;
329 
336  void
337  createCoarser(LevelData<MFCellFAB>& a_coarse, const LevelData<MFCellFAB>& a_fine, bool a_ghosted) override final;
338 
345  void
346  createCoarsened(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs, const int& a_refRat) override final;
347 
354  void
355  incr(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs, Real a_scale) override final;
356 
360  Real
361  dotProduct(const LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_2) override final;
362 
368  void
369  scale(LevelData<MFCellFAB>& a_lhs, const Real& a_scale) override final;
370 
375  void
376  setToZero(LevelData<MFCellFAB>& a_lhs) override final;
377 
383  void
384  assign(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs) override final;
385 
394  void
395  assignCopier(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs, const Copier& a_copier) override final;
396 
402  void
403  assignLocal(LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs) override final;
404 
411  void
412  buildCopier(Copier& a_copier, const LevelData<MFCellFAB>& a_lhs, const LevelData<MFCellFAB>& a_rhs) override;
413 
419  Real
420  norm(const LevelData<MFCellFAB>& a_lhs, int a_order) override final;
421 
430  void
431  axby(LevelData<MFCellFAB>& a_lhs,
432  const LevelData<MFCellFAB>& a_x,
433  const LevelData<MFCellFAB>& a_y,
434  const Real a_a,
435  const Real a_b) override final;
436 
443  void
444  restrictResidual(LevelData<MFCellFAB>& a_resCoar,
445  LevelData<MFCellFAB>& a_phi,
446  const LevelData<MFCellFAB>& a_rhs) override final;
447 
453  void
454  prolongIncrement(LevelData<MFCellFAB>& a_phi, const LevelData<MFCellFAB>& a_correctCoarse) override final;
455 
462  void
463  AMRUpdateResidual(LevelData<MFCellFAB>& a_residual,
464  const LevelData<MFCellFAB>& a_correction,
465  const LevelData<MFCellFAB>& a_coarseCorrection) override final;
466 
475  void
476  AMRRestrict(LevelData<MFCellFAB>& a_residualCoarse,
477  const LevelData<MFCellFAB>& a_residual,
478  const LevelData<MFCellFAB>& a_correction,
479  const LevelData<MFCellFAB>& a_coarseCorrection,
480  bool a_skip_res) override final;
481 
487  void
488  AMRProlong(LevelData<MFCellFAB>& a_correction, const LevelData<MFCellFAB>& a_coarseCorrection) override final;
489 
500  void
501  AMRResidual(LevelData<MFCellFAB>& a_residual,
502  const LevelData<MFCellFAB>& a_phiFine,
503  const LevelData<MFCellFAB>& a_phi,
504  const LevelData<MFCellFAB>& a_phiCoar,
505  const LevelData<MFCellFAB>& a_rhs,
506  bool a_homogeneousPhysBC,
507  AMRLevelOp<LevelData<MFCellFAB>>* a_finerOp) override final;
508 
517  void
518  AMRResidualNF(LevelData<MFCellFAB>& a_residual,
519  const LevelData<MFCellFAB>& a_phi,
520  const LevelData<MFCellFAB>& a_phiCoar,
521  const LevelData<MFCellFAB>& a_rhs,
522  bool a_homogeneousPhysBC) override final;
523 
533  void
534  AMRResidualNC(LevelData<MFCellFAB>& a_residual,
535  const LevelData<MFCellFAB>& a_phiFine,
536  const LevelData<MFCellFAB>& a_phi,
537  const LevelData<MFCellFAB>& a_rhs,
538  bool a_homogeneousPhysBC,
539  AMRLevelOp<LevelData<MFCellFAB>>* a_finerOp) override final;
540 
549  void
550  AMROperatorNF(LevelData<MFCellFAB>& a_Lphi,
551  const LevelData<MFCellFAB>& a_phi,
552  const LevelData<MFCellFAB>& a_phiCoar,
553  bool a_homogeneousPhysBC) override final;
554 
563  void
564  AMROperatorNC(LevelData<MFCellFAB>& a_Lphi,
565  const LevelData<MFCellFAB>& a_phi,
566  const LevelData<MFCellFAB>& a_phiCoar,
567  bool a_homogeneousPhysBC,
568  AMRLevelOp<LevelData<MFCellFAB>>* a_finerOp) override final;
569 
570 protected:
574  static constexpr int m_comp = 0;
575 
579  static constexpr int m_nComp = 1;
580 
585 
590 
594  Vector<IntVect> m_colors;
595 
599  std::map<int, RefCountedPtr<EBHelmholtzOp>> m_helmOps;
600 
604  RefCountedPtr<MFHelmholtzJumpBC> m_jumpBC;
605 
609  RefCountedPtr<LevelData<BaseIVFAB<Real>>> m_jump;
610 
614  std::map<int, RefCountedPtr<LevelData<BaseIVFAB<Real>>>> m_dirichletBcValues;
615 
620 
625 
630 
635 
640 
645 
650 
654  RefCountedPtr<LevelData<MFCellFAB>> m_Acoef;
655 
659  RefCountedPtr<LevelData<MFFluxFAB>> m_Bcoef;
660 
664  RefCountedPtr<LevelData<MFBaseIVFAB>> m_BcoefIrreg;
665 
670 
675 
680 
685 
690 
694  IntVect m_ghostPhi;
695 
699  IntVect m_ghostRhs;
700 
704  bool m_hasCoar;
705 
709  bool m_hasFine;
710 
715 
721  void
722  updateJumpBC(const LevelData<MFCellFAB>& a_phi, const bool a_homogeneousPhysBC);
723 
728  void
729  exchangeGhost(const LevelData<MFCellFAB>& a_phi) const;
730 
738  void
739  interpolateCF(const LevelData<MFCellFAB>& a_phi, const LevelData<MFCellFAB>* a_phiCoar, const bool a_homogeneousCF);
740 
750  void
751  applyOp(LevelData<MFCellFAB>& a_Lphi,
752  const LevelData<MFCellFAB>& a_phi,
753  const LevelData<MFCellFAB>* const a_phiCoar,
754  const bool a_homogeneousPhysBC,
755  const bool a_homogeneousCFBC);
756 
767  void
768  AMROperator(LevelData<MFCellFAB>& a_Lphi,
769  const LevelData<MFCellFAB>& a_phiFine,
770  const LevelData<MFCellFAB>& a_phi,
771  const LevelData<MFCellFAB>& a_phiCoar,
772  const bool a_homogeneousPhysBC,
773  AMRLevelOp<LevelData<MFCellFAB>>* a_finerOp) override;
774 };
775 
776 #include <CD_NamespaceFooter.H>
777 
778 #include <CD_MFHelmholtzOpImplem.H>
779 
780 #endif
Declaration of Helmholtz multigrid operators.
Declaration of cell positions.
Declaration of a multiphase BaseIVFAB<Real>
Wrapper class for holding multifluid EBCoarAves.
Declaration of an EB boundary condition factory class for MFHelmholtzOp.
Declaration of a factory class for making MFHelmholtzJumpBC objects for use in a multifluid MFHelmhol...
Declaration of class for computing "jump interface" boundary conditions for multifluid Helmholtz code...
Implementation of CD_MFHelmholtzOp.H.
Declaration of a wrapper for wrapping multifluid EBLevelGrids.
Declaration of a wrapper class for holding multifluid EBMultigridInterpolators.
Declaration of a class for refluxing in a multiphase context.
Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho.
Definition: CD_EBHelmholtzOp.H:41
Multiphase BaseIVFAB<Real>.
Definition: CD_MFBaseIVFAB.H:29
Class for coarsening data in a multifluid context.
Definition: CD_MFCoarAve.H:26
Class for computing "jump interface" boundary conditions for multifluid code.
Definition: CD_MFHelmholtzJumpBC.H:42
Operator for solving multifluid Helmholtz on a grid level.
Definition: CD_MFHelmholtzOp.H:43
void fillGrad(const LevelData< MFCellFAB > &a_phi) override final
Not called, I think.
Definition: CD_MFHelmholtzOp.cpp:316
RefCountedPtr< MFHelmholtzJumpBC > m_jumpBC
BC jump object. This is the one that has the stencils and can compute derivatives.
Definition: CD_MFHelmholtzOp.H:604
MFHelmholtzOp()=delete
Constructor. Must subsequently call define(...)
void AMRResidual(LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
Compute residual on this level. AMR version.
Definition: CD_MFHelmholtzOp.cpp:1048
Copier m_exchangeCopier
Copier for exchange operation.
Definition: CD_MFHelmholtzOp.H:669
void prolongIncrement(LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_correctCoarse) override final
Prolongation method.
Definition: CD_MFHelmholtzOp.cpp:939
void scale(LevelData< MFCellFAB > &a_lhs, const Real &a_scale) override final
Scale function.
Definition: CD_MFHelmholtzOp.cpp:348
MFHelmholtzOp(const MFHelmholtzOp &a_op)=delete
No copy construction allowed.
void AMROperatorNC(LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels.
Definition: CD_MFHelmholtzOp.cpp:1134
void incr(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, Real a_scale) override final
Increment function.
Definition: CD_MFHelmholtzOp.cpp:340
const RefCountedPtr< LevelData< MFFluxFAB > > & getBcoef()
Get the Helmholtz B-coefficient on faces.
Definition: CD_MFHelmholtzOp.cpp:248
Real dotProduct(const LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_2) override final
Dot product.
Definition: CD_MFHelmholtzOp.cpp:415
void applyOpNoBoundary(LevelData< MFCellFAB > &a_ans, const LevelData< MFCellFAB > &a_phi) override final
Apply operator but turn off all BCs.
Definition: CD_MFHelmholtzOp.cpp:300
MFLevelGrid m_mflg
Level grid.
Definition: CD_MFHelmholtzOp.H:619
void preCond(LevelData< MFCellFAB > &a_corr, const LevelData< MFCellFAB > &a_residual) override final
Precondition system before bottom solve.
Definition: CD_MFHelmholtzOp.cpp:533
MFLevelGrid m_mflgCoar
Coarse grid.
Definition: CD_MFHelmholtzOp.H:629
void getFlux(MFFluxFAB &a_flux, const LevelData< MFCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale) override final
Fill flux.
Definition: CD_MFHelmholtzOp.cpp:330
static constexpr int m_nComp
Number of components that we solve for.
Definition: CD_MFHelmholtzOp.H:579
void createCoarser(LevelData< MFCellFAB > &a_coarse, const LevelData< MFCellFAB > &a_fine, bool a_ghosted) override final
Create coarsened data.
Definition: CD_MFHelmholtzOp.cpp:501
void AMRProlong(LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection) override final
Prolongation onto AMR level.
Definition: CD_MFHelmholtzOp.cpp:1032
void relax(LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, int a_iterations) override final
Relaxation method. This does smoothing for the system L(correction) = residual.
Definition: CD_MFHelmholtzOp.cpp:716
void AMROperator(LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override
Apply the AMR operator, i.e. compute L(phi) in an AMR context.
Definition: CD_MFHelmholtzOp.cpp:1184
MFCoarAve m_coarAve
Coarsener.
Definition: CD_MFHelmholtzOp.H:649
void assign(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
Assignment fucntion.
Definition: CD_MFHelmholtzOp.cpp:364
int m_numPhases
Number of phases.
Definition: CD_MFHelmholtzOp.H:684
void applyOp(LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, bool a_homogeneousPhysBc) override final
Apply operator.
Definition: CD_MFHelmholtzOp.cpp:554
bool m_hasFine
True if there is a finer AMR level.
Definition: CD_MFHelmholtzOp.H:709
const RefCountedPtr< LevelData< MFBaseIVFAB > > & getBcoefIrreg()
Get the Helmholtz B-coefficient on the EB.
Definition: CD_MFHelmholtzOp.cpp:254
RefCountedPtr< LevelData< MFFluxFAB > > m_Bcoef
Helmholtz B-coefficeint.
Definition: CD_MFHelmholtzOp.H:659
void restrictResidual(LevelData< MFCellFAB > &a_resCoar, LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs) override final
Restrict residual onto coarse level.
Definition: CD_MFHelmholtzOp.cpp:912
IntVect m_ghostPhi
Number of ghost cells.
Definition: CD_MFHelmholtzOp.H:694
std::map< int, RefCountedPtr< LevelData< BaseIVFAB< Real > > > > m_dirichletBcValues
Dirichlet BC values for each phase.
Definition: CD_MFHelmholtzOp.H:614
Location::Cell m_dataLocation
Interpretation of data. Either on cell center or on cell centroid.
Definition: CD_MFHelmholtzOp.H:584
void divideByIdentityCoef(LevelData< MFCellFAB > &a_rhs) override final
Divide by the a-coefficient.
Definition: CD_MFHelmholtzOp.cpp:286
RefCountedPtr< EBHelmholtzOp > & getHelmOp(const int a_phase)
Get Helmholtz operator.
void AMROperatorNF(LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &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_MFHelmholtzOp.cpp:1098
~MFHelmholtzOp()
Destructor.
Definition: CD_MFHelmholtzOp.cpp:207
void AMRUpdateResidual(LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection) override final
Update AMR residual.
Definition: CD_MFHelmholtzOp.cpp:955
void axby(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_x, const LevelData< MFCellFAB > &a_y, const Real a_a, const Real a_b) override final
Set a_lhs = a*x + b*y.
Definition: CD_MFHelmholtzOp.cpp:672
void residual(LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs, const bool a_homogeneousPhysBc) override final
Compute residual on this level.
Definition: CD_MFHelmholtzOp.cpp:658
void relaxGSMultiColor(LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
Multi-colored gauss-seidel relaxation.
Definition: CD_MFHelmholtzOp.cpp:854
static constexpr int m_comp
Component that we solve for.
Definition: CD_MFHelmholtzOp.H:574
bool m_multifluid
Multifluid operator or not.
Definition: CD_MFHelmholtzOp.H:674
void AMRResidualNC(LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
Compute AMR residual on coarsest.
Definition: CD_MFHelmholtzOp.cpp:1082
void AMRRestrict(LevelData< MFCellFAB > &a_residualCoarse, const LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection, bool a_skip_res) override final
Restrict residual.
Definition: CD_MFHelmholtzOp.cpp:992
bool m_hasMGObjects
Has MG objects or not.
Definition: CD_MFHelmholtzOp.H:679
Vector< IntVect > m_colors
"Colors" for the multi-coloered relaxation method
Definition: CD_MFHelmholtzOp.H:594
void relaxPointJacobi(LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition: CD_MFHelmholtzOp.cpp:747
void interpolateCF(const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > *a_phiCoar, const bool a_homogeneousCF)
Do coarse-fine interpolation.
Definition: CD_MFHelmholtzOp.cpp:607
void setAcoAndBco(const RefCountedPtr< LevelData< MFCellFAB >> &a_Acoef, const RefCountedPtr< LevelData< MFFluxFAB >> &a_Bcoef, const RefCountedPtr< LevelData< MFBaseIVFAB >> &a_BcoefIrreg)
Update operators with new coefficients.
Definition: CD_MFHelmholtzOp.cpp:215
IntVect m_ghostRhs
Number of ghost cells.
Definition: CD_MFHelmholtzOp.H:699
bool m_hasMGObjcts
True if there are multigrid levels.
Definition: CD_MFHelmholtzOp.H:714
void AMRResidualNF(LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC) override final
Compute AMR residual on finest AMR level.
Definition: CD_MFHelmholtzOp.cpp:1066
void relaxGSRedBlack(LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition: CD_MFHelmholtzOp.cpp:799
void setJump(RefCountedPtr< LevelData< BaseIVFAB< Real >>> &a_jump)
Set the jump boundary condition.
Definition: CD_MFHelmholtzOp.cpp:260
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_jump
Actual BC jump in data-based format. This is the right-hand side of dphi/dn1 + dphi/dn2 = jump.
Definition: CD_MFHelmholtzOp.H:609
void setToZero(LevelData< MFCellFAB > &a_lhs) override final
Set to zero.
Definition: CD_MFHelmholtzOp.cpp:356
void createCoarsened(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, const int &a_refRat) override final
Create coarsening of data holder.
Definition: CD_MFHelmholtzOp.cpp:517
std::map< int, RefCountedPtr< EBHelmholtzOp > > m_helmOps
Helmholtz operators on each phase. Note that I'm using int rather than Phase as identifier because th...
Definition: CD_MFHelmholtzOp.H:599
Smoother m_smoother
Relaxation method.
Definition: CD_MFHelmholtzOp.H:589
void updateJumpBC(const LevelData< MFCellFAB > &a_phi, const bool a_homogeneousPhysBC)
Update the jump condition.
Definition: CD_MFHelmholtzOp.cpp:694
Vector< long long > computeOperatorLoads(LevelData< MFCellFAB > &a_phi, const int a_numApply) noexcept
Time the applyOp routine. Template parameter is std::chrono duration. E.g. std::chrono::microseconds.
Definition: CD_MFHelmholtzOpImplem.H:19
MFLevelGrid m_mflgCoFi
Coarsened version of this grid.
Definition: CD_MFHelmholtzOp.H:634
MFHelmholtzOp(const MFHelmholtzOp &&a_op)=delete
No move construction allowed.
bool m_hasCoar
True if there is a coarser AMR level.
Definition: CD_MFHelmholtzOp.H:704
void assignCopier(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, const Copier &a_copier) override final
Assign lhs.
Definition: CD_MFHelmholtzOp.cpp:372
int refToCoarser() override final
Return coarsening factor to coarser level (1 if there is no coarser level);.
Definition: CD_MFHelmholtzOp.cpp:268
MFLevelGrid m_mflgFine
Fine level grid.
Definition: CD_MFHelmholtzOp.H:624
Real norm(const LevelData< MFCellFAB > &a_lhs, int a_order) override final
Compute solution norm.
Definition: CD_MFHelmholtzOp.cpp:396
MFHelmholtzOp & operator=(const MFHelmholtzOp &a_oper)=delete
No copy assigment allowed.
void setAlphaAndBeta(const Real &a_alpha, const Real &a_beta) override final
Set alpha and beta.
Definition: CD_MFHelmholtzOp.cpp:276
RefCountedPtr< LevelData< MFBaseIVFAB > > m_BcoefIrreg
Helmholtz B-coefficeint on EB.
Definition: CD_MFHelmholtzOp.H:664
int m_refToCoar
Refinement factor to coarser AMR level.
Definition: CD_MFHelmholtzOp.H:689
const RefCountedPtr< LevelData< MFCellFAB > > & getAcoef()
Get the Helmholtz A-coefficient on cell centers.
Definition: CD_MFHelmholtzOp.cpp:242
RefCountedPtr< LevelData< MFCellFAB > > m_Acoef
Helmholtz A-coefficeint.
Definition: CD_MFHelmholtzOp.H:654
MFHelmholtzOp & operator=(const MFHelmholtzOp &&a_oper)=delete
No move assigment allowed.
MFMultigridInterpolator m_interpolator
Multi grid interpolator.
Definition: CD_MFHelmholtzOp.H:644
MFLevelGrid m_mflgCoarMG
Coarse multigrid-grid.
Definition: CD_MFHelmholtzOp.H:639
void assignLocal(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
Local assignment function.
Definition: CD_MFHelmholtzOp.cpp:380
void buildCopier(Copier &a_copier, const LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override
Build copier.
Definition: CD_MFHelmholtzOp.cpp:388
Smoother
Relaxation methods for this operator.
Definition: CD_MFHelmholtzOp.H:49
void create(LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
Create method.
Definition: CD_MFHelmholtzOp.cpp:485
void exchangeGhost(const LevelData< MFCellFAB > &a_phi) const
Perform an exchange operation, event if the data is const.
Definition: CD_MFHelmholtzOp.cpp:702
Wrapper class for holding multifluid EBLevelGrids.
Definition: CD_MFLevelGrid.H:29
Wrapper class for holding multifluid EBMultigridInterpolators.
Definition: CD_MFMultigridInterpolator.H:26
Class which wraps EBRefluxs in multiphase.
Definition: CD_MFReflux.H:26
Namespace for encapsulating various data centerings.
Definition: CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition: CD_Location.H:30