chombo-discharge
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 
40 class EBHelmholtzOp : public LevelTGAHelmOp<LevelData<EBCellFAB>, EBFluxFAB>
41 {
42 public:
46  enum class Smoother
47  {
48  NoRelax,
49  PointJacobi,
50  GauSaiRedBlack,
51  GauSaiMultiColor,
52  };
53 
57  EBHelmholtzOp() = delete;
58 
63  EBHelmholtzOp(const EBHelmholtzOp& a_other) = delete;
64 
69  EBHelmholtzOp(const EBHelmholtzOp&& a_other) = delete;
70 
100  EBHelmholtzOp(const Location::Cell a_dataLocation,
101  const EBLevelGrid& a_eblgFine,
102  const EBLevelGrid& a_eblg,
103  const EBLevelGrid& a_eblgCoFi,
104  const EBLevelGrid& a_eblgCoar,
105  const EBLevelGrid& a_eblgCoarMG,
106  const RefCountedPtr<EBMultigridInterpolator>& a_interpolator,
107  const RefCountedPtr<EBReflux>& a_fluxReg,
108  const RefCountedPtr<EBCoarAve>& a_coarAve,
109  const RefCountedPtr<EBHelmholtzDomainBC>& a_domainBC,
110  const RefCountedPtr<EBHelmholtzEBBC>& a_ebBC,
111  const RealVect& a_probLo,
112  const Real& a_dx,
113  const int& a_refToFine,
114  const int& a_refToCoar,
115  const bool& a_hasFine,
116  const bool& a_hasCoar,
117  const bool& a_hasMGObjects,
118  const Real& a_alpha,
119  const Real& a_beta,
120  const RefCountedPtr<LevelData<EBCellFAB>>& a_Acoef,
121  const RefCountedPtr<LevelData<EBFluxFAB>>& a_Bcoef,
122  const RefCountedPtr<LevelData<BaseIVFAB<Real>>>& a_BcoIrreg,
123  const IntVect& a_ghostCellsPhi,
124  const IntVect& a_ghostCellsRHS,
125  const Smoother& a_smoother);
126 
130  virtual ~EBHelmholtzOp();
131 
137  operator=(const EBHelmholtzOp& a_oper) = delete;
138 
144  operator=(const EBHelmholtzOp&& a_oper) = delete;
145 
149  void
150  turnOffCFInterp();
151 
155  void
156  turnOnCFInterp();
157 
162  void
163  turnOffExchange();
164 
169  void
170  turnOnExchange();
171 
176  void
178 
183  void
185 
193  void
194  setAcoAndBco(const RefCountedPtr<LevelData<EBCellFAB>>& a_Acoef,
195  const RefCountedPtr<LevelData<EBFluxFAB>>& a_Bcoef,
196  const RefCountedPtr<LevelData<BaseIVFAB<Real>>>& a_BcoefIrreg);
197 
202  const RefCountedPtr<LevelData<EBCellFAB>>&
203  getAcoef();
204 
209  const RefCountedPtr<LevelData<EBFluxFAB>>&
210  getBcoef();
211 
216  const RefCountedPtr<LevelData<BaseIVFAB<Real>>>&
217  getBcoefIrreg();
218 
225  void
226  coarsenCell(LevelData<EBCellFAB>& a_phi, const LevelData<EBCellFAB>& a_phiFine);
227 
235  void
236  coarsenFlux(LevelData<EBFluxFAB>& a_flux, const LevelData<EBFluxFAB>& a_fineFlux);
237 
249  void
250  pointJacobiKernel(EBCellFAB& a_Lcorr,
251  EBCellFAB& a_corr,
252  const EBCellFAB& a_resid,
253  const EBCellFAB& a_Acoef,
254  const EBFluxFAB& a_Bcoef,
255  const BaseIVFAB<Real>& a_BcoefIrreg,
256  const Box& a_cellBox,
257  const DataIndex& a_dit) const noexcept;
258 
271  void
272  gauSaiRedBlackKernel(EBCellFAB& a_Lcorr,
273  EBCellFAB& a_corr,
274  const EBCellFAB& a_resid,
275  const EBCellFAB& a_Acoef,
276  const EBFluxFAB& a_Bcoef,
277  const BaseIVFAB<Real>& a_BcoefIrreg,
278  const Box& a_cellBox,
279  const DataIndex& a_dit,
280  const int& a_redBlack) const noexcept;
281 
294  void
295  gauSaiMultiColorKernel(EBCellFAB& a_Lcorr,
296  EBCellFAB& a_corr,
297  const EBCellFAB& a_resid,
298  const EBCellFAB& a_Acoef,
299  const EBFluxFAB& a_Bcoef,
300  const BaseIVFAB<Real>& a_BcoefIrreg,
301  const Box& a_cellBox,
302  const DataIndex& a_dit,
303  const IntVect& a_color) const noexcept;
304 
312  void
313  residual(LevelData<EBCellFAB>& a_residual,
314  const LevelData<EBCellFAB>& a_phi,
315  const LevelData<EBCellFAB>& a_rhs,
316  const bool a_homogeneousPhysBc) override;
317 
324  void
325  preCond(LevelData<EBCellFAB>& a_corr, const LevelData<EBCellFAB>& a_residual) override final;
326 
333  void
334  interpolateCF(LevelData<EBCellFAB>& a_phiFine, const LevelData<EBCellFAB>* a_phiCoar, const bool a_homogeneousCFBC);
335 
340  void
341  homogeneousCFInterp(LevelData<EBCellFAB>& a_phi);
342 
348  void
349  inhomogeneousCFInterp(LevelData<EBCellFAB>& a_phi, const LevelData<EBCellFAB>& a_phiCoar);
350 
358  void
359  gauSaiMultiColor(LevelData<EBCellFAB>& a_phi,
360  const LevelData<EBCellFAB>& a_Lphi,
361  const LevelData<EBCellFAB>& a_rhs,
362  const IntVect& a_color) const;
363 
373  void
374  applyOp(LevelData<EBCellFAB>& a_Lphi,
375  const LevelData<EBCellFAB>& a_phi,
376  const LevelData<EBCellFAB>* const a_phiCoar,
377  const bool a_homogeneousPhysBC,
378  const bool a_homogeneousCFBC);
379 
387  void
388  applyOp(LevelData<EBCellFAB>& a_Lphi, const LevelData<EBCellFAB>& a_phi, bool a_homogeneousPhysBc) override final;
389 
398  void
399  applyOp(EBCellFAB& a_Lphi,
400  EBCellFAB& a_phi,
401  const EBCellFAB& a_Acoef,
402  const EBFluxFAB& a_Bcoef,
403  const BaseIVFAB<Real>& a_BcoefIrreg,
404  const Box& a_cellBox,
405  const DataIndex& a_dit,
406  const bool a_homogeneousPhysBC) const noexcept;
407 
416  void
417  applyOpRegular(EBCellFAB& a_Lphi,
418  EBCellFAB& a_phi,
419  const EBCellFAB& a_Acoef,
420  const EBFluxFAB& a_Bcoef,
421  const BaseIVFAB<Real>& a_BcoefIrreg,
422  const Box& a_cellBox,
423  const DataIndex& a_dit,
424  const bool a_homogeneousPhysBC) const noexcept;
425 
436  void
437  applyDomainFlux(EBCellFAB& a_phi,
438  const EBFluxFAB& a_Bcoef,
439  const Box& a_cellBox,
440  const DataIndex& a_dit,
441  const bool a_homogeneousPhysBc) const noexcept;
442 
451  void
452  fillDomainFlux(EBFluxFAB& a_flux, const EBCellFAB& a_phi, const Box& a_cellBox, const DataIndex& a_dit);
453 
462  void
463  applyOpIrregular(EBCellFAB& a_Lphi,
464  const EBCellFAB& a_phi,
465  const EBCellFAB& a_Acoef,
466  const EBFluxFAB& a_Bcoef,
467  const BaseIVFAB<Real>& a_BcoefIrreg,
468  const BaseIVFAB<Real>& a_alphaDiagWeight,
469  const Box& a_cellBox,
470  const DataIndex& a_dit,
471  const bool a_homogeneousPhysBC) const noexcept;
472 
478  void
479  create(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override final;
480 
485  void
486  assign(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override final;
487 
496  void
497  assignCopier(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs, const Copier& a_copier) override final;
498 
504  void
505  assignLocal(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override final;
506 
513  void
514  buildCopier(Copier& a_copier, const LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override;
515 
519  Real
520  dotProduct(const LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs) override final;
521 
528  void
529  incr(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs, const Real a_scale) override final;
530 
539  void
540  axby(LevelData<EBCellFAB>& a_lhs,
541  const LevelData<EBCellFAB>& a_x,
542  const LevelData<EBCellFAB>& a_y,
543  const Real a_a,
544  const Real a_b) override final;
545 
551  void
552  scale(LevelData<EBCellFAB>& a_lhs, const Real& a_scale) override final;
553 
559  Real
560  norm(const LevelData<EBCellFAB>& a_rhs, const int a_order) override final;
561 
566  void
567  setToZero(LevelData<EBCellFAB>& a_lhs) override final;
568 
575  void
576  createCoarser(LevelData<EBCellFAB>& a_coarse, const LevelData<EBCellFAB>& a_fine, bool a_ghosted) override final;
577 
584  void
585  relax(LevelData<EBCellFAB>& a_correction, const LevelData<EBCellFAB>& a_residual, int a_iterations) override final;
586 
593  void
594  restrictResidual(LevelData<EBCellFAB>& a_resCoar,
595  LevelData<EBCellFAB>& a_phi,
596  const LevelData<EBCellFAB>& a_rhs) override final;
597 
603  void
604  prolongIncrement(LevelData<EBCellFAB>& a_phi, const LevelData<EBCellFAB>& a_correctCoarse) override final;
605 
609  int
610  refToCoarser() override final;
611 
622  void
623  AMROperator(LevelData<EBCellFAB>& a_Lphi,
624  const LevelData<EBCellFAB>& a_phiFine,
625  const LevelData<EBCellFAB>& a_phi,
626  const LevelData<EBCellFAB>& a_phiCoar,
627  const bool a_homogeneousPhysBC,
628  AMRLevelOp<LevelData<EBCellFAB>>* a_finerOp) override final;
629 
641  void
642  refluxFreeAMROperator(LevelData<EBCellFAB>& a_Lphi,
643  const LevelData<EBCellFAB>& a_phiFine,
644  const LevelData<EBCellFAB>& a_phi,
645  const LevelData<EBCellFAB>& a_phiCoar,
646  const bool a_homogeneousPhysBC,
647  AMRLevelOp<LevelData<EBCellFAB>>* a_finerOp);
648 
659  void
660  AMRResidual(LevelData<EBCellFAB>& a_residual,
661  const LevelData<EBCellFAB>& a_phiFine,
662  const LevelData<EBCellFAB>& a_phi,
663  const LevelData<EBCellFAB>& a_phiCoar,
664  const LevelData<EBCellFAB>& a_rhs,
665  bool a_homogeneousPhysBC,
666  AMRLevelOp<LevelData<EBCellFAB>>* a_finerOp) override final;
667 
676  void
677  AMRResidualNF(LevelData<EBCellFAB>& a_residual,
678  const LevelData<EBCellFAB>& a_phi,
679  const LevelData<EBCellFAB>& a_phiCoar,
680  const LevelData<EBCellFAB>& a_rhs,
681  bool a_homogeneousPhysBC) override final;
682 
692  void
693  AMRResidualNC(LevelData<EBCellFAB>& a_residual,
694  const LevelData<EBCellFAB>& a_phiFine,
695  const LevelData<EBCellFAB>& a_phi,
696  const LevelData<EBCellFAB>& a_rhs,
697  bool a_homogeneousPhysBC,
698  AMRLevelOp<LevelData<EBCellFAB>>* a_finerOp) override final;
699 
708  void
709  AMROperatorNF(LevelData<EBCellFAB>& a_Lphi,
710  const LevelData<EBCellFAB>& a_phi,
711  const LevelData<EBCellFAB>& a_phiCoar,
712  bool a_homogeneousPhysBC) override final;
713 
722  void
723  AMROperatorNC(LevelData<EBCellFAB>& a_Lphi,
724  const LevelData<EBCellFAB>& a_phi,
725  const LevelData<EBCellFAB>& a_phiCoar,
726  bool a_homogeneousPhysBC,
727  AMRLevelOp<LevelData<EBCellFAB>>* a_finerOp) override final;
728 
737  void
738  AMRRestrict(LevelData<EBCellFAB>& a_residualCoarse,
739  const LevelData<EBCellFAB>& a_residual,
740  const LevelData<EBCellFAB>& a_correction,
741  const LevelData<EBCellFAB>& a_coarseCorrection,
742  bool a_skip_res) override final;
743 
749  void
750  AMRProlong(LevelData<EBCellFAB>& a_correction, const LevelData<EBCellFAB>& a_coarseCorrection) override final;
751 
758  void
759  AMRUpdateResidual(LevelData<EBCellFAB>& a_residual,
760  const LevelData<EBCellFAB>& a_correction,
761  const LevelData<EBCellFAB>& a_coarseCorrection) override final;
762 
768  void
769  setAlphaAndBeta(const Real& a_alpha, const Real& a_beta) override final;
770 
777  void
778  createCoarsened(LevelData<EBCellFAB>& a_lhs, const LevelData<EBCellFAB>& a_rhs, const int& a_refRat) override final;
779 
785  void
786  diagonalScale(LevelData<EBCellFAB>& a_rhs, bool a_kappaWeighted) override final;
787 
792  void
793  divideByIdentityCoef(LevelData<EBCellFAB>& a_rhs) override final;
794 
798  void
799  applyOpNoBoundary(LevelData<EBCellFAB>& a_ans, const LevelData<EBCellFAB>& a_phi) override final;
800 
805  void
806  fillGrad(const LevelData<EBCellFAB>& a_phi) override final;
807 
816  void
817  getFlux(EBFluxFAB& a_flux,
818  const LevelData<EBCellFAB>& a_data,
819  const Box& a_grid,
820  const DataIndex& a_dit,
821  Real a_scale) override final;
822 
826  void
827  allocateFlux() const noexcept;
828 
832  void
833  deallocateFlux() const noexcept;
834 
838  LevelData<EBFluxFAB>&
839  getFlux() const;
840 
841 protected:
845  static constexpr int m_nComp = 1;
846 
850  static constexpr int m_comp = 0;
851 
856 
860  Interval m_interval;
861 
866 
871 
876 
880  bool m_profile;
881 
886 
890  bool m_hasFine;
891 
895  bool m_hasCoar;
896 
901 
906 
912 
918 
924 
928  IntVect m_ghostPhi;
929 
933  IntVect m_ghostRhs;
934 
938  Real m_alpha;
939 
943  Real m_beta;
944 
948  Real m_dx;
949 
953  RealVect m_vecDx;
954 
958  RealVect m_probLo;
959 
964 
969 
973  EBLevelGrid m_eblgFine;
974 
978  EBLevelGrid m_eblg;
979 
983  EBLevelGrid m_eblgCoFi;
984 
988  EBLevelGrid m_eblgCoar;
989 
993  EBLevelGrid m_eblgCoarMG;
994 
999 
1004 
1009 
1014 
1019 
1023  RefCountedPtr<EBHelmholtzEBBC> m_ebBc;
1024 
1029 
1033  RefCountedPtr<EBReflux> m_fluxReg;
1034 
1038  RefCountedPtr<EBCoarAve> m_coarAve;
1039 
1043  RefCountedPtr<LevelData<EBCellFAB>> m_Acoef;
1044 
1048  RefCountedPtr<LevelData<EBFluxFAB>> m_Bcoef;
1049 
1053  RefCountedPtr<LevelData<BaseIVFAB<Real>>> m_BcoefIrreg;
1054 
1058  LevelData<EBCellFAB> m_relCoef;
1059 
1063  mutable LevelData<EBFluxFAB>* m_flux;
1064 
1068  LayoutData<BaseIFFAB<Real>> m_interpolant[SpaceDim];
1069 
1073  LayoutData<BaseIFFAB<FaceStencil>> m_interpStencil[SpaceDim];
1074 
1078  LayoutData<BaseIFFAB<VoFStencil>> m_centroidFluxStencil[SpaceDim];
1079 
1085  LayoutData<BaseIVFAB<VoFStencil>> m_relaxStencils;
1086 
1091  LayoutData<RefCountedPtr<VCAggStencil>> m_aggRelaxStencil;
1092 
1096  LayoutData<BaseIVFAB<Real>> m_alphaDiagWeight;
1097 
1101  LayoutData<BaseIVFAB<Real>> m_betaDiagWeight;
1102 
1106  mutable LayoutData<VoFIterator> m_vofIterIrreg;
1107 
1111  mutable LayoutData<VoFIterator> m_vofIterMulti;
1112 
1116  mutable LayoutData<VoFIterator> m_vofIterStenc;
1117 
1121  mutable LayoutData<VoFIterator> m_vofIterDomLo[SpaceDim];
1122 
1126  mutable LayoutData<VoFIterator> m_vofIterDomHi[SpaceDim];
1127 
1131  std::map<std::pair<int, Side::LoHiSide>, Box> m_sideBox;
1132 
1136  Vector<IntVect> m_colors;
1137 
1144  void
1145  relaxPointJacobi(LevelData<EBCellFAB>& a_correction, const LevelData<EBCellFAB>& a_residual, const int a_iterations);
1146 
1153  void
1154  relaxGSRedBlack(LevelData<EBCellFAB>& a_correction, const LevelData<EBCellFAB>& a_residual, const int a_iterations);
1155 
1162  void
1163  relaxGSMultiColor(LevelData<EBCellFAB>& a_correction, const LevelData<EBCellFAB>& a_residual, const int a_iterations);
1164 
1168  void
1170 
1174  void
1176 
1180  void
1181  makeAggStencil();
1182 
1186  void
1187  defineStencils();
1188 
1196  VoFStencil
1197  getFaceCenterFluxStencil(const FaceIndex& a_face, const DataIndex& a_dit) const;
1198 
1206  VoFStencil
1207  getFaceCentroidFluxStencil(const FaceIndex& a_face, const DataIndex& a_dit) const;
1208 
1217  void
1218  computeFlux(EBFaceFAB& a_fluxCentroid,
1219  const EBCellFAB& a_phi,
1220  const Box& a_cellBox,
1221  const DataIndex& a_dit,
1222  const int a_dir);
1223 
1232  void
1233  computeFaceCenteredFlux(EBFaceFAB& a_fluxCenter,
1234  const EBCellFAB& a_phi,
1235  const Box& a_cellBox,
1236  const DataIndex& a_dit,
1237  const int a_dir);
1238 
1247  void
1248  computeFaceCentroidFlux(EBFaceFAB& a_flux,
1249  const EBCellFAB& a_phi,
1250  const Box& a_cellBox,
1251  const DataIndex& a_dit,
1252  const int a_dir);
1253 
1259  void
1260  computeFlux(const LevelData<EBCellFAB>& a_phi);
1261 
1269  void
1270  reflux(LevelData<EBCellFAB>& a_Lphi,
1271  const LevelData<EBCellFAB>& a_phiFine,
1272  const LevelData<EBCellFAB>& a_phi,
1273  AMRLevelOp<LevelData<EBCellFAB>>& a_finerOp);
1274 };
1275 
1276 #include <CD_NamespaceFooter.H>
1277 
1278 #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:2091
void relaxGSMultiColor(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Multi-colored gauss-seidel relaxation.
Definition: CD_EBHelmholtzOp.cpp:1789
LevelData< EBFluxFAB > & getFlux() const
Returns m_flux. This is used in refluxing routines.
Definition: CD_EBHelmholtzOp.cpp:261
static constexpr int m_nComp
Number of components that we solve for (always one..)
Definition: CD_EBHelmholtzOp.H:845
RefCountedPtr< EBCoarAve > m_coarAve
Conservative coarsener.
Definition: CD_EBHelmholtzOp.H:1038
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:158
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:625
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:518
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:850
void turnOnCFInterp()
Turn on BCs.
Definition: CD_EBHelmholtzOp.cpp:205
EBLevelGrid m_eblgCoFi
Coarsened of m_eblg.
Definition: CD_EBHelmholtzOp.H:983
const RefCountedPtr< LevelData< BaseIVFAB< Real > > > & getBcoefIrreg()
Get the Helmholtz B-coefficient on the EB.
Definition: CD_EBHelmholtzOp.cpp:186
void AMRProlong(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection) override final
Prolongation onto AMR level.
Definition: CD_EBHelmholtzOp.cpp:1075
Smoother m_smoother
Relaxation method.
Definition: CD_EBHelmholtzOp.H:865
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:1185
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:1085
EBHelmholtzOp & operator=(const EBHelmholtzOp &&a_oper)=delete
No move assigment allowed.
bool m_hasCoar
True if there's a coarser level.
Definition: CD_EBHelmholtzOp.H:895
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:1015
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:2214
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:1239
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:1111
void allocateFlux() const noexcept
Allocate m_flux.
Definition: CD_EBHelmholtzOp.cpp:245
void computeDiagWeight()
Calculate the weight of the diagonal term.
Definition: CD_EBHelmholtzOp.cpp:1902
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:637
EBHelmholtzOp()=delete
Disallowed default constructor.
virtual ~EBHelmholtzOp()
Dtor.
Definition: CD_EBHelmholtzOp.cpp:191
Copier m_exchangeCopierFine
Pre-built exchange copier.
Definition: CD_EBHelmholtzOp.H:968
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:1581
void turnOffCoarsening()
Turn off coarsening operation.
Definition: CD_EBHelmholtzOp.cpp:229
void relaxPointJacobi(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition: CD_EBHelmholtzOp.cpp:1615
RefCountedPtr< EBHelmholtzEBBC > m_ebBc
Domain bc object.
Definition: CD_EBHelmholtzOp.H:1023
void makeAggStencil()
Compute aggregated stencils.
Definition: CD_EBHelmholtzOp.cpp:2021
void homogeneousCFInterp(LevelData< EBCellFAB > &a_phi)
Do homogeneous coarse-fine interpolation.
Definition: CD_EBHelmholtzOp.cpp:1530
EBLevelGrid m_eblg
Grid.
Definition: CD_EBHelmholtzOp.H:978
LayoutData< RefCountedPtr< VCAggStencil > > m_aggRelaxStencil
For making irregular stencil applications go faster.
Definition: CD_EBHelmholtzOp.H:1091
EBMGRestrict m_restrictOp
Restriction operator for AMR levels.
Definition: CD_EBHelmholtzOp.H:998
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:801
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:1655
bool m_refluxFree
Use reflux-free formulation or not.
Definition: CD_EBHelmholtzOp.H:875
RefCountedPtr< EBHelmholtzDomainBC > m_domainBc
Domain bc object.
Definition: CD_EBHelmholtzOp.H:1018
Real dotProduct(const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Compute the dot product??
Definition: CD_EBHelmholtzOp.cpp:550
RefCountedPtr< EBReflux > m_fluxReg
Flux register.
Definition: CD_EBHelmholtzOp.H:1033
void preCond(LevelData< EBCellFAB > &a_corr, const LevelData< EBCellFAB > &a_residual) override final
Precondition system before bottom solve.
Definition: CD_EBHelmholtzOp.cpp:510
bool m_doInterpCF
Do coarse-fine interpolation or not.
Definition: CD_EBHelmholtzOp.H:911
int m_refToFine
Refinement factor to fine level.
Definition: CD_EBHelmholtzOp.H:905
Smoother
Relaxation method for the operators.
Definition: CD_EBHelmholtzOp.H:47
LevelData< EBFluxFAB > * m_flux
For holding fluxes.
Definition: CD_EBHelmholtzOp.H:1063
bool m_hasMGObjects
True if there is a multigrid level below this operator.
Definition: CD_EBHelmholtzOp.H:885
Real m_dx
Grid resolution;.
Definition: CD_EBHelmholtzOp.H:948
bool m_profile
Profile the operator.
Definition: CD_EBHelmholtzOp.H:880
void turnOffExchange()
Turn off exchange operation.
Definition: CD_EBHelmholtzOp.cpp:213
void relaxGSRedBlack(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
Jacobi relaxation.
Definition: CD_EBHelmholtzOp.cpp:1681
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
Weights of diagonal alpha terms.
Definition: CD_EBHelmholtzOp.H:1096
void fillGrad(const LevelData< EBCellFAB > &a_phi) override final
Not called, I think.
Definition: CD_EBHelmholtzOp.cpp:1510
VoFStencil getFaceCenterFluxStencil(const FaceIndex &a_face, const DataIndex &a_dit) const
Get the face-centered flux stencil.
Definition: CD_EBHelmholtzOp.cpp:2072
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:729
void divideByIdentityCoef(LevelData< EBCellFAB > &a_rhs) override final
Divide by the a-coefficient.
Definition: CD_EBHelmholtzOp.cpp:1484
void assign(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Assign data.
Definition: CD_EBHelmholtzOp.cpp:526
RealVect m_probLo
Lower-left corner of domain.
Definition: CD_EBHelmholtzOp.H:958
RefCountedPtr< EBMultigridInterpolator > m_interpolator
Interpolator object.
Definition: CD_EBHelmholtzOp.H:1028
void computeRelaxationCoefficient()
Calculate relaxation coefficient.
Definition: CD_EBHelmholtzOp.cpp:1957
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:2186
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:973
EBMGProlong m_prolongOp
Prolongation operator for AMR levels.
Definition: CD_EBHelmholtzOp.H:1008
IntVect m_ghostRhs
Ghost cells for rhs (note, the operator rhs)
Definition: CD_EBHelmholtzOp.H:933
int refToCoarser() override final
Return coarsening factor to coarser level (1 if there is no coarser level);.
Definition: CD_EBHelmholtzOp.cpp:754
LayoutData< BaseIFFAB< VoFStencil > > m_centroidFluxStencil[SpaceDim]
Face centroid flux stencil. Defined on all faces connecting one or more irregular vofs.
Definition: CD_EBHelmholtzOp.H:1078
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:760
Interval m_interval
Interval.
Definition: CD_EBHelmholtzOp.H:860
EBHelmholtzOp & operator=(const EBHelmholtzOp &a_oper)=delete
No copy assigment allowed.
void applyOpNoBoundary(LevelData< EBCellFAB > &a_ans, const LevelData< EBCellFAB > &a_phi) override final
Apply operator but turn off all BCs.
Definition: CD_EBHelmholtzOp.cpp:1500
void assignLocal(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
Local assignment function.
Definition: CD_EBHelmholtzOp.cpp:542
bool m_doExchange
Turn on/off exchange operation.
Definition: CD_EBHelmholtzOp.H:917
EBMGProlong m_prolongOpMG
Prolongation operator for MG levels.
Definition: CD_EBHelmholtzOp.H:1013
Real norm(const LevelData< EBCellFAB > &a_rhs, const int a_order) override final
Compute norm of data.
Definition: CD_EBHelmholtzOp.cpp:645
LayoutData< BaseIFFAB< FaceStencil > > m_interpStencil[SpaceDim]
Face centroid interpolation stencil.
Definition: CD_EBHelmholtzOp.H:1073
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:497
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:1326
EBHelmholtzOp(const EBHelmholtzOp &&a_other)=delete
Disallowed move constructor.
RefCountedPtr< LevelData< EBFluxFAB > > m_Bcoef
B-coefficient in Helmholtz equation.
Definition: CD_EBHelmholtzOp.H:1048
void deallocateFlux() const noexcept
Deallocate m_flux.
Definition: CD_EBHelmholtzOp.cpp:253
Timer m_timer
Timer so user can profile.
Definition: CD_EBHelmholtzOp.H:855
Location::Cell m_dataLocation
Data centering.
Definition: CD_EBHelmholtzOp.H:870
void createCoarser(LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted) override final
Create coarsened data.
Definition: CD_EBHelmholtzOp.cpp:707
EBLevelGrid m_eblgCoar
Coarse level grid (if the operator has a coarse level)
Definition: CD_EBHelmholtzOp.H:988
IntVect m_ghostPhi
Ghost cells for phi.
Definition: CD_EBHelmholtzOp.H:928
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:2167
EBLevelGrid m_eblgFine
Fine level grid (if the operator has a fine level)
Definition: CD_EBHelmholtzOp.H:973
void turnOffCFInterp()
Turn off BCs.
Definition: CD_EBHelmholtzOp.cpp:197
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:2268
int m_refToCoar
Refinement factor to coarse level.
Definition: CD_EBHelmholtzOp.H:900
Real m_alpha
Alpha-coefficient.
Definition: CD_EBHelmholtzOp.H:938
LevelData< EBCellFAB > m_relCoef
Relaxation coefficient.
Definition: CD_EBHelmholtzOp.H:1058
void diagonalScale(LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted) override final
Do diagonal scaling.
Definition: CD_EBHelmholtzOp.cpp:1462
void setToZero(LevelData< EBCellFAB > &a_lhs) override final
Set data to zero.
Definition: CD_EBHelmholtzOp.cpp:699
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:483
Copier m_exchangeCopier
Pre-built exchange copier.
Definition: CD_EBHelmholtzOp.H:963
void interpolateCF(LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > *a_phiCoar, const bool a_homogeneousCFBC)
Apply coarse-fine boundary conditions.
Definition: CD_EBHelmholtzOp.cpp:1554
void coarsenCell(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiFine)
Coarsen data from fine to coar level.
Definition: CD_EBHelmholtzOp.cpp:2306
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:1728
LayoutData< BaseIFFAB< Real > > m_interpolant[SpaceDim]
Interpolant for when we want centroid fluxes.
Definition: CD_EBHelmholtzOp.H:1068
void turnOnExchange()
Turn on exchange operation.
Definition: CD_EBHelmholtzOp.cpp:221
void inhomogeneousCFInterp(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar)
Inhomogeneous coarse-fine interpolation.
Definition: CD_EBHelmholtzOp.cpp:1542
void coarsenFlux(LevelData< EBFluxFAB > &a_flux, const LevelData< EBFluxFAB > &a_fineFlux)
Coarsen fluxes on the fine level onto this level.
Definition: CD_EBHelmholtzOp.cpp:2314
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:1030
bool m_hasFine
True if there's a finer level.
Definition: CD_EBHelmholtzOp.H:890
bool m_doCoarsen
Turn on/off exchange operation.
Definition: CD_EBHelmholtzOp.H:923
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:984
void assignCopier(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Copier &a_copier) override final
Assign lhs.
Definition: CD_EBHelmholtzOp.cpp:534
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:1835
std::map< std::pair< int, Side::LoHiSide >, Box > m_sideBox
Domain boxes on each side.
Definition: CD_EBHelmholtzOp.H:1131
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:1116
const RefCountedPtr< LevelData< EBCellFAB > > & getAcoef()
Get the Helmholtz A-coefficient on cell centers.
Definition: CD_EBHelmholtzOp.cpp:174
RealVect m_vecDx
Vector resolution.
Definition: CD_EBHelmholtzOp.H:953
void incr(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Real a_scale) override final
Increment operator.
Definition: CD_EBHelmholtzOp.cpp:617
LayoutData< VoFIterator > m_vofIterMulti
VoFIterator for "multi-cells".
Definition: CD_EBHelmholtzOp.H:1111
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
Weights of diagonal beta terms.
Definition: CD_EBHelmholtzOp.H:1101
EBHelmholtzOp(const EBHelmholtzOp &a_other)=delete
Disallowed copy constructor.
Vector< IntVect > m_colors
"Colors" for the relaxation methods
Definition: CD_EBHelmholtzOp.H:1136
const RefCountedPtr< LevelData< EBFluxFAB > > & getBcoef()
Get the Helmholtz B-coefficient on faces.
Definition: CD_EBHelmholtzOp.cpp:180
RefCountedPtr< LevelData< EBCellFAB > > m_Acoef
A-coefficient in Helmholtz equation.
Definition: CD_EBHelmholtzOp.H:1043
void turnOnCoarsening()
Turn on coarsening operation.
Definition: CD_EBHelmholtzOp.cpp:237
LayoutData< VoFIterator > m_vofIterDomLo[SpaceDim]
VoF iterators for lo domain side.
Definition: CD_EBHelmholtzOp.H:1121
LayoutData< VoFIterator > m_vofIterIrreg
VoFIterator for irregular cells.
Definition: CD_EBHelmholtzOp.H:1106
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:1047
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:997
EBMGRestrict m_restrictOpMG
Restriction operator if this is an MG level.
Definition: CD_EBHelmholtzOp.H:1003
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:717
Real m_beta
Beta-coefficient.
Definition: CD_EBHelmholtzOp.H:943
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:1380
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:1083
LayoutData< VoFIterator > m_vofIterDomHi[SpaceDim]
VoF iterators for hi domain side.
Definition: CD_EBHelmholtzOp.H:1126
void buildCopier(Copier &a_copier, const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override
Build copier.
Definition: CD_EBHelmholtzOp.cpp:2322
void defineStencils()
Define stencils.
Definition: CD_EBHelmholtzOp.cpp:269
void prolongIncrement(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_correctCoarse) override final
Prolongation method.
Definition: CD_EBHelmholtzOp.cpp:746
EBLevelGrid m_eblgCoarMG
Coarser grids (multigrid level)
Definition: CD_EBHelmholtzOp.H:993
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_BcoefIrreg
B-coefficient in Helmholtz equation, but on EB faces.
Definition: CD_EBHelmholtzOp.H:1053
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
Namespace for encapsulating various data centerings.
Definition: CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition: CD_Location.H:30