Chombo + EB  3.2
EBConductivityOp.H
Go to the documentation of this file.
1 #ifdef CH_LANG_CC
2 /*
3  * _______ __
4  * / ___/ / ___ __ _ / / ___
5  * / /__/ _ \/ _ \/ V \/ _ \/ _ \
6  * \___/_//_/\___/_/_/_/_.__/\___/
7  * Please refer to Copyright.txt, in Chombo's root directory.
8  */
9 #endif
10 
11 #ifndef _EBCONDUCTIVITYOP_H_
12 #define _EBCONDUCTIVITYOP_H_
13 
14 #include "REAL.H"
15 #include "Box.H"
16 #include "FArrayBox.H"
17 #include "Vector.H"
18 #include <map>
19 #include "RefCountedPtr.H"
20 
21 #include "AMRMultiGrid.H"
22 #include "EBAMRPoissonOp.H"
23 
24 #include "EBIndexSpace.H"
25 #include "EBCellFAB.H"
26 #include "EBCellFactory.H"
27 #include "EBStencil.H"
28 
29 #include "EBLevelDataOps.H"
30 #include "BaseEBBC.H"
31 #include "BaseDomainBC.H"
32 #include "CFIVS.H"
33 #include "EBFluxRegister.H"
34 #include "EBFastFR.H"
35 #include "EBMGAverage.H"
36 #include "EBMGInterp.H"
37 #include "PolyGeom.H"
38 #include "EBQuadCFInterp.H"
39 #include "EBLevelGrid.H"
40 #include "AMRTGA.H"
41 #include "AMRPoissonOp.H"
42 #include "CFRegion.H"
44 #include "NamespaceHeader.H"
45 
46 
47 #if CH_SPACEDIM==2
48 #define EBAMRPOVC_NUMSTEN 4
49 #elif CH_SPACEDIM==3
50 #define EBAMRPOVC_NUMSTEN 8
51 #else
53 {
54  THIS_WILL_ONLY_COMPILE_WHEN_CH_SPACEDIM_IS_2_OR_3;
55 }
56 #endif
57 
58 //! \class EBConductivityOp
59 //! This class implements an operator that solves the equation
60 //! (alpha a + beta div (b grad) )phi = rhs
61 //! using the AMRLevelOp interface.
62 class EBConductivityOp: public LevelTGAHelmOp<LevelData<EBCellFAB>, EBFluxFAB >
63 {
64  public:
65 
66 #ifdef CH_USE_HDF5
67  ///for debugging
68  virtual void dumpAMR(Vector<LevelData<EBCellFAB>*>& a_data, string name)
69  {
70  writeEBAMRname(&a_data, name.c_str());
71  }
72 
73  ///for debugging
74  virtual void dumpLevel(LevelData<EBCellFAB>& a_data, string name)
75  {
76  writeEBLevelname(&a_data, name.c_str());
77  }
78 #endif
79  //! Constructs a conductivity operator using the given data. This
80  //! constructor is for time-independent a and b coefficients.
81  //! If you are approaching this operator from this interface, consider backing away and using
82  //! EBConductivityOpFactory to generate these objects. Really.
83  //! Ghost cell arguments are there for caching reasons. Once you set them,
84  //! an error is thrown if you send in data that does not match.
85  //! \param a_eblgFine grid at finer level
86  //! \param a_eblg grid at this level
87  //! \param a_eblgCoar grid at coarser level
88  //! \param a_eblgCoarMG grid at intermediate multigrid level
89  //! \param a_domainBC domain boundary conditions at this level
90  //! \param a_ebBC eb boundary conditions at this level
91  //! \param a_dx grid spacing at this level
92  //! \param a_origin offset to lowest corner of the domain
93  //! \param a_refToFine refinement ratio to finer level
94  //! \param a_refToCoar refinement ratio to coarser level
95  //! \param a_hasFiner true if there is a finer AMR level, false otherwise.
96  //! \param a_hasCoarser true if there is a coarser AMR level.
97  //! \param a_hasCoarserMG true if there is a coarser MultiGrid level.
98  //! \param a_preCondIters number of iterations to do for pre-conditioning
99  //! \param a_relaxType 0 means point Jacobi, 1 is Gauss-Seidel.
100  //! \param a_acoef coefficent of identity
101  //! \param a_bcoef coefficient of gradient.
102  //! \param a_ghostCellsPhi Number of ghost cells in phi, correction
103  //! \param a_ghostCellsRhs Number of ghost cells in RHS, residual, lphi
104  EBConductivityOp(const EBLevelGrid & a_eblgFine,
105  const EBLevelGrid & a_eblg,
106  const EBLevelGrid & a_eblgCoar,
107  const EBLevelGrid & a_eblgCoarMG,
108  const RefCountedPtr<EBQuadCFInterp>& a_quadCFI,
109  const RefCountedPtr<ConductivityBaseDomainBC>& a_domainBC,
111  const Real & a_dx,
112  const Real & a_dxCoar,
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 bool& a_layoutChanged,
119  const Real& a_alpha,
120  const Real& a_beta,
121  const RefCountedPtr<LevelData<EBCellFAB> >& a_acoef,
122  const RefCountedPtr<LevelData<EBFluxFAB> >& a_bcoef,
123  const RefCountedPtr<LevelData<BaseIVFAB<Real> > >& a_bcoIrreg,
124  const IntVect& a_ghostCellsPhi,
125  const IntVect& a_ghostCellsRHS,
126  const int& a_relaxType);
127 
128 
129  //! Destructor.
131 
132  Real getSafety();
133 
134  ///
135  /**
136  This sets the data storage for the a coefficient to a different
137  object and recalculates the stuff it depends on. Use this only if you know what you're doing.
138  */
140  {
141 
142  m_acoef = a_acoef;
145  }
146 
147  //only weights by kappa. time and tga have their demands.
148  virtual void kappaScale(LevelData<EBCellFAB> & a_rhs);
149 
150  ///returns m_dx, such function is required by some LinearSolvers
151  Real dx() const
152  {
153  return m_dx;
154  }
155 
156  ///
157  /** a_residual = a_rhs - L(a_phiFine, a_phi) no coaser AMR level*/
158  void AMRResidualNC(LevelData<EBCellFAB>& a_residual,
159  const LevelData<EBCellFAB>& a_phiFine,
160  const LevelData<EBCellFAB>& a_phi,
161  const LevelData<EBCellFAB>& a_rhs,
162  bool a_homogeneousBC,
163  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
164 
165 
166  ///
167  /** apply AMR operator no coaser AMR level*/
168  void AMROperatorNC(LevelData<EBCellFAB>& a_LofPhi,
169  const LevelData<EBCellFAB>& a_phiFine,
170  const LevelData<EBCellFAB>& a_phi,
171  bool a_homogeneousBC,
172  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
173 
174  //------------------------------------
175  // Overridden methods for base classes
176  //------------------------------------
177 
178  void setAlphaAndBeta(const Real& a_alpha,
179  const Real& a_beta);
180 
181  void diagonalScale(LevelData<EBCellFAB> & a_rhs,
182  bool a_kappaWeighted);
183 
185 
186  void fillGrad(const LevelData<EBCellFAB>& a_phi);
187 
188  void getFlux(EBFluxFAB& a_flux,
189  const LevelData<EBCellFAB>& a_data,
190  const Box& a_grid,
191  const DataIndex& a_dit,
192  Real a_scale);
193 
194  void getFlux(EBFaceFAB& a_fluxCentroid,
195  const EBCellFAB& a_phi,
196  const Box& a_ghostedBox,
197  const Box& a_fabBox,
198  const ProblemDomain& a_domain,
199  const EBISBox& a_ebisBox,
200  const Real& a_dx,
201  const DataIndex& a_datInd,
202  const int& a_idir);
203 
204  // This operator may be time-dependent.
205  void setTime(Real a_oldTime, Real a_mu, Real a_dt);
206 
207  //! This is called on multigrid operators when their AMR operators
208  //! are altered.
209  void finerOperatorChanged(const MGLevelOp<LevelData<EBCellFAB> >& a_operator,
210  int a_coarseningFactor);
211 
212  //MGOp operations. no finer or coarser
213 
214  ///
215  /**
216  */
217  virtual void residual(LevelData<EBCellFAB>& a_residual,
218  const LevelData<EBCellFAB>& a_phi,
219  const LevelData<EBCellFAB>& a_rhs,
220  bool a_homogeneousPhysBC=false);
221 
222  ///
223  /**
224  */
225  virtual void preCond(LevelData<EBCellFAB>& a_opPhi,
226  const LevelData<EBCellFAB>& a_phi);
227 
228  ///
229  /**
230  This function assumes that coarse-fine boundary condtions have
231  been dealt with.
232  */
233  virtual void applyOp(LevelData<EBCellFAB>& a_opPhi,
234  const LevelData<EBCellFAB>& a_phi,
235  const LevelData<EBCellFAB>* const a_phiCoarse,
236  const bool& a_homogeneousPhysBC,
237  const bool& a_homogeneousCFBC);
238 
239  /// virtual function called by LevelTGA
240  virtual void applyOpNoBoundary(LevelData<EBCellFAB>& a_opPhi,
241  const LevelData<EBCellFAB>& a_phi);
242 
243  ///
244  /**
245  this is the linearop function. CFBC is set to homogeneous. phic is null
246  */
247  virtual void applyOp(LevelData<EBCellFAB>& a_opPhi,
248  const LevelData<EBCellFAB>& a_phi,
249  bool a_homogeneousPhysBC);
250 
251  ///
252  /**
253  */
254  virtual void create(LevelData<EBCellFAB>& a_lhs,
255  const LevelData<EBCellFAB>& a_rhs);
256 
257  ///
258  virtual void createCoarsened(LevelData<EBCellFAB>& a_lhs,
259  const LevelData<EBCellFAB>& a_rhs,
260  const int& a_refRat);
261 
262  Real
263  AMRNorm(const LevelData<EBCellFAB>& a_coarResid,
264  const LevelData<EBCellFAB>& a_fineResid,
265  const int& a_refRat,
266  const int& a_ord);
267 
268  ///
269  /**
270  */
271  virtual void assign(LevelData<EBCellFAB>& a_lhs,
272  const LevelData<EBCellFAB>& a_rhs);
273 
274  ///
275  /**
276  */
277  virtual Real dotProduct(const LevelData<EBCellFAB>& a_1,
278  const LevelData<EBCellFAB>& a_2);
279 
280  ///
281  /**
282  */
283  virtual void incr(LevelData<EBCellFAB>& a_lhs,
284  const LevelData<EBCellFAB>& a_x,
285  Real a_scale);
286 
287  ///
288  /**
289  */
290  virtual void axby(LevelData<EBCellFAB>& a_lhs,
291  const LevelData<EBCellFAB>& a_x,
292  const LevelData<EBCellFAB>& a_y,
293  Real a_a,
294  Real a_b);
295 
296  ///
297  /**
298  */
299  virtual void scale(LevelData<EBCellFAB>& a_lhs,
300  const Real& a_scale);
301 
302  ///
303  /**
304  */
305  virtual Real norm(const LevelData<EBCellFAB>& a_rhs,
306  int a_ord);
307 
308  virtual Real localMaxNorm(const LevelData<EBCellFAB>& a_rhs);
309  ///
310  /**
311  */
312  virtual void setToZero(LevelData<EBCellFAB>& a_lhs);
313 
314  ///
315  /**
316  */
317  virtual void setVal(LevelData<EBCellFAB>& a_lhs, const Real& a_value);
318 
319  ///
320  /**
321  */
322  virtual void createCoarser(LevelData<EBCellFAB>& a_coarse,
323  const LevelData<EBCellFAB>& a_fine,
324  bool a_ghosted);
325 
326  ///
327  /**
328  */
329  virtual void relax(LevelData<EBCellFAB>& a_e,
330  const LevelData<EBCellFAB>& a_residual,
331  int a_iterations);
332 
333 
334  virtual void relaxGauSai(LevelData<EBCellFAB>& a_e,
335  const LevelData<EBCellFAB>& a_residual,
336  int a_iterations);
337 
338  virtual void relaxPoiJac(LevelData<EBCellFAB>& a_e,
339  const LevelData<EBCellFAB>& a_residual,
340  int a_iterations);
341 
342 
343  virtual void relaxGSRBFast(LevelData<EBCellFAB>& a_e,
344  const LevelData<EBCellFAB>& a_residual,
345  int a_iterations);
346 
347  ///
348  /**
349  Calculate restricted residual:
350  a_resCoarse[2h] = I[h->2h] (a_rhsFine[h] - L[h](a_phiFine[h]))
351  */
352  virtual void restrictResidual(LevelData<EBCellFAB>& a_resCoarse,
353  LevelData<EBCellFAB>& a_phiFine,
354  const LevelData<EBCellFAB>& a_rhsFine);
355 
356  ///
357  /**
358  Correct the fine solution based on coarse correction:
359  a_phiThisLevel += I[2h->h] (a_correctCoarse)
360  */
361  virtual void prolongIncrement(LevelData<EBCellFAB>& a_phiThisLevel,
362  const LevelData<EBCellFAB>& a_correctCoarse);
363 
364  ///
365  /** Refinement ratio between this level and coarser level.
366  Returns 1 when there are no coarser AMRLevelOp objects */
367  virtual int refToCoarser();
368 
369  ///
370  /** Refinement ratio between this level and coarser level.
371  Returns 1 when there are no coarser AMRLevelOp objects */
372  virtual int refToFiner();
373 
374  ///
375  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
376  virtual void AMRResidual(LevelData<EBCellFAB>& a_residual,
377  const LevelData<EBCellFAB>& a_phiFine,
378  const LevelData<EBCellFAB>& a_phi,
379  const LevelData<EBCellFAB>& a_phiCoarse,
380  const LevelData<EBCellFAB>& a_rhs,
381  bool a_homogeneousBC,
382  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
383 
384  ///
385  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
386  virtual void AMRResidualNF(LevelData<EBCellFAB>& a_residual,
387  const LevelData<EBCellFAB>& a_phi,
388  const LevelData<EBCellFAB>& a_phiCoarse,
389  const LevelData<EBCellFAB>& a_rhs,
390  bool a_homogeneousBC);
391 
392 
393  ///
394  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
395  virtual void AMROperator(LevelData<EBCellFAB>& a_LofPhi,
396  const LevelData<EBCellFAB>& a_phiFine,
397  const LevelData<EBCellFAB>& a_phi,
398  const LevelData<EBCellFAB>& a_phiCoarse,
399  bool a_homogeneousBC,
400  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
401 
402  ///
403  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
404  virtual void AMROperatorNF(LevelData<EBCellFAB>& a_LofPhi,
405  const LevelData<EBCellFAB>& a_phi,
406  const LevelData<EBCellFAB>& a_phiCoarse,
407  bool a_homogeneousBC);
408 
409  ///
410  /** a_resCoarse = I[h-2h] (a_residual - L(a_correction, a_coarseCorrection)) */
411  virtual void AMRRestrict(LevelData<EBCellFAB>& a_resCoarse,
412  const LevelData<EBCellFAB>& a_residual,
413  const LevelData<EBCellFAB>& a_correction,
414  const LevelData<EBCellFAB>& a_coarseCorrection,
415  bool a_skip_res = false );
416 
417  ///
418  /** a_correction += I[2h->h](a_coarseCorrection) */
419  virtual void AMRProlong(LevelData<EBCellFAB>& a_correction,
420  const LevelData<EBCellFAB>& a_coarseCorrection);
421 
422  ///
423  /** a_residual = a_residual - L(a_correction, a_coarseCorrection) */
424  virtual void AMRUpdateResidual(LevelData<EBCellFAB>& a_residual,
425  const LevelData<EBCellFAB>& a_correction,
426  const LevelData<EBCellFAB>& a_coarseCorrection);
427 
428  void reflux(LevelData<EBCellFAB>& a_residual,
429  const LevelData<EBCellFAB>& a_phiFine,
430  const LevelData<EBCellFAB>& a_phi,
431  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
432 
433  void gsrbColor(LevelData<EBCellFAB>& a_phi,
434  const LevelData<EBCellFAB>& a_lph,
435  const LevelData<EBCellFAB>& a_rhs,
436  const IntVect& a_color);
437 
438  void getDivFStencil(VoFStencil& a_vofStencil,
439  const VolIndex& a_vof,
440  const DataIndex& a_dit);
441 
442  void getFluxStencil(VoFStencil& a_fluxStencil,
443  const FaceIndex& a_face,
444  const DataIndex& a_dit);
445 
446  void getFaceCenteredFluxStencil(VoFStencil& a_fluxStencil,
447  const FaceIndex& a_face,
448  const DataIndex& a_dit);
449 
450  void incrOpRegularDir(EBCellFAB& a_lhs,
451  const EBCellFAB& a_phi,
452  const bool& a_homogeneous,
453  const int& a_dir,
454  const DataIndex& a_datInd);
455  void applyOpIrregular(EBCellFAB& a_lhs,
456  const EBCellFAB& a_phi,
457  const bool& a_homogeneous,
458  const DataIndex& a_datInd);
459 
460  //this is needed to do wacky things like reset the a coefficient.
461  //not for the faint of heart
462  //but none but the brave deserve the fair
463  const EBLevelGrid& getEBLG() const
464  {
465  return m_eblg;
466  }
467  static void setForceNoEBCF(bool a_forceNoEBCF)
468  {
469  s_forceNoEBCF = a_forceNoEBCF;
470  }
471  ///do not call this one unless you really know what you are doing
472  void defineStencils();
473 
474  //adding suplementary methods for petsc solver
475  void getAlphaDiagWeight(LayoutData<BaseIVFAB<Real> > const*& a_alphaDiagWeight)
476  {
477  a_alphaDiagWeight = &m_alphaDiagWeight;
478  }
479 
480  void getAlphaBeta(Real& a_alpha, Real& a_beta)
481  {
482  a_alpha = m_alpha;
483  a_beta = m_beta;
484  }
486  {
487  return m_domainBC;
488  }
489 
491  {
492  return m_acoef;
493  }
494 
496  {
497  return m_bcoef;
498  }
499 
500  void getEBBCFluxStencil(LayoutData<BaseIVFAB<VoFStencil> > const*& a_ebbcFluxStencil)
501  {
502  a_ebbcFluxStencil = m_ebBC->getFluxStencil(0);
503  }
504 
505  static int s_numComps;
506  static int s_whichComp;
507 
508 protected:
509  void incrOpRegularAllDirs(Box * a_loBox,
510  Box * a_hiBox,
511  int * a_hasLo,
512  int * a_hasHi,
513  Box & a_curDblBox,
514  Box & a_curPhiBox,
515  int a_nComps,
516  BaseFab<Real> & a_curOpPhiFAB,
517  const BaseFab<Real> & a_curPhiFAB,
518  bool a_homogeneousPhysBC,
519  const DataIndex& a_dit);
520 
521  void applyDomainFlux(Box * a_loBox,
522  Box * a_hiBox,
523  int * a_hasLo,
524  int * a_hasHi,
525  Box & a_dblBox,
526  int a_nComps,
527  BaseFab<Real> & a_phiFAB,
528  bool a_homogeneousPhysBC,
529  const DataIndex& a_dit);
530 
531  void GSColorAllIrregular(EBCellFAB& a_phi,
532  const EBCellFAB& a_rhs,
533  const int& a_icolor,
534  const DataIndex& a_dit);
535 
536  static bool s_turnOffBCs;
537  static bool s_forceNoEBCF;
539  void dumpFABPoint(const EBCellFAB& a_fab, const DataIndex& a_dit, const string& a_blab);
540  void dumpLevelPoint(const LevelData<EBCellFAB>& a_res, const string& a_blab);
541 
542  virtual void calculateAlphaWeight();
543  virtual void calculateRelaxationCoefficient();
544 
545  void defineColorStencils(Box a_ideBoxLo[SpaceDim],
546  Box a_ideBoxHi[SpaceDim]);
547  //EBCF gymnastics
548  void defineEBCFStencils();
549  void getFluxEBCF(EBFaceFAB& a_flux,
550  const EBCellFAB& a_phi,
551  const Box& a_ghostedBox,
552  Vector<FaceIndex>& a_faceitEBCF,
553  Vector<VoFStencil>& a_stenEBCF);
554 
555  void getFluxRegOnly(EBFaceFAB& a_fluxCentroid,
556  const EBCellFAB& a_phi,
557  const Box& a_ghostedBox,
558  const Real& a_dx,
559  const DataIndex& a_datInd,
560  const int& a_idir);
561 
562  //stuff to make EBCF go faster
565 
566  //stuff to make relaxation go faster
567  //not for the faint of heart
573 
574 
575 
579 
581 
587 
588  // RefCountedPtr<ConductivityBaseDomainBC> m_domainBC;
589  // RefCountedPtr<ConductivityBaseEBBC> m_ebBC;
592 
596 
597  //! "Current" (time-interpolated) value of the a coefficient. For a
598  //! time-independent a coefficient, this is where the coefficient lives.
600 
601  //! "Current" (time-interpolated) value of the b coefficient. For a
602  //! time-independent a coefficient, this is where the coefficient lives.
604 
605  //! "Current" (time-interpolated) value of the b coefficient on irregular
606  //! cells.
608 
611  //weights that get multiplied by alpha
613  //weights that get multiplied by beta
617  bool m_hasEBCF;
618  bool m_hasFine;
620  bool m_hasCoar;
621 
622  //restriction object
624  //prolongation object
626 
627  //stencils for operator evaluation
629  //stencils for operator evaluation on gauss-seidel colors
630 
631  //! Multigrid relaxation coefficient
633 
634  //cache the vofiterators
635  //for irregular cell iteration (includes buffer around multivalued cells)
638  //for domain boundary conditions at ir regular cells
641 
642 
643  // Coarse-fine stencils for homogeneous CFInterp
646 
647  //flux register with finer level
649  //special mg objects for when we do not have
650  //a coarser level or when the refinement to coarse
651  //is greater than two
652  //flag for when we need special MG objects
655  //stuff below is only defined if m_hasMGObjects==true
661 
663 
664 private:
665 
666  void incrementFRCoar(EBFastFR& a_fluxReg,
667  const LevelData<EBCellFAB>& a_phiFine,
668  const LevelData<EBCellFAB>& a_phi);
669 
670  void incrementFRFine(EBFastFR& a_fluxReg,
671  const LevelData<EBCellFAB>& a_phiFine,
672  const LevelData<EBCellFAB>& a_phi,
673  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
674 
675  void getFlux(FArrayBox& a_flux,
676  const FArrayBox& a_phi,
677  const Box& a_faceBox,
678  const int& a_idir,
679  const Real& a_dx,
680  const DataIndex& a_datInd);
681 
682 
683 
684  void applyCFBCs(LevelData<EBCellFAB>& a_phi,
685  const LevelData<EBCellFAB>* const a_phiCoarse,
686  bool a_homogeneousCFBC);
687 
688  void getOpVoFStencil(VoFStencil& a_stencil,
689  const EBISBox& a_ebisbox,
690  const VolIndex& a_vof);
691 
692  void getOpVoFStencil(VoFStencil& a_stencil,
693  const int& a_dir,
694  const Vector<VolIndex>& a_allMonotoneVoFs,
695  const EBISBox& a_ebisbox,
696  const VolIndex& a_vof,
697  const bool& a_lowOrder);
698 
699 
700  void getOpFaceStencil(VoFStencil& a_stencil,
701  const Vector<VolIndex>& a_allMonotoneVofs,
702  const EBISBox& a_ebisbox,
703  const VolIndex& a_vof,
704  int a_dir,
705  const Side::LoHiSide& a_side,
706  const FaceIndex& a_face,
707  const bool& a_lowOrder);
708 
709  void levelJacobi(LevelData<EBCellFAB>& a_phi,
710  const LevelData<EBCellFAB>& a_rhs);
711 
713 
714  void applyHomogeneousCFBCs(EBCellFAB& a_phi,
715  const DataIndex& a_datInd,
716  int a_idir,
717  Side::LoHiSide a_hiorlo);
718 private:
719 
720  //! Default constructor. Creates an undefined conductivity operator.
722 
723  //copy constructor and operator= disallowed for all the usual reasons
724  EBConductivityOp(const EBConductivityOp& a_opin)
725  {
726  MayDay::Error("invalid operator");
727  }
728 
729  void operator=(const EBConductivityOp& a_opin)
730  {
731  MayDay::Error("invalid operator");
732  }
733 };
734 
735 
736 #include "NamespaceFooter.H"
737 #endif
void operator=(const EBConductivityOp &a_opin)
Definition: EBConductivityOp.H:729
LayoutData< VoFIterator > m_vofIterDomLo[CH_SPACEDIM]
Definition: EBConductivityOp.H:639
LevelData< EBCellFAB > m_relCoef
Multigrid relaxation coefficient.
Definition: EBConductivityOp.H:632
EBLevelGrid m_eblgCoarMG
Definition: EBConductivityOp.H:585
static int s_numComps
Definition: EBConductivityOp.H:505
void setAlphaAndBeta(const Real &a_alpha, const Real &a_beta)
virtual void setToZero(LevelData< EBCellFAB > &a_lhs)
virtual void restrictResidual(LevelData< EBCellFAB > &a_resCoarse, LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_rhsFine)
void fillGrad(const LevelData< EBCellFAB > &a_phi)
These functions are part of the LevelTGA interface......
LayoutData< BaseIVFAB< Real > > m_cacheEBxDomainFluxHi[EBAMRPOVC_NUMSTEN][SpaceDim]
Definition: EBConductivityOp.H:571
void levelJacobi(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
Vector< IntVect > m_colors
Definition: EBConductivityOp.H:662
void getFlux(EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale)
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
Definition: EBConductivityOp.H:612
EBFastFR-A class to encapsulate a levels worth of flux registers.
Definition: EBFastFR.H:39
bool m_hasMGObjects
Definition: EBConductivityOp.H:653
Real m_dx
Definition: EBConductivityOp.H:594
#define CH_SPACEDIM
Definition: SPACE.H:51
virtual void AMRUpdateResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
void getFaceCenteredFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit)
LayoutData< RefCountedPtr< EBStencil > > m_opEBStencil
Definition: EBConductivityOp.H:628
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:141
void THIS_IS_AN_ERROR_MESSAGE(void)
Definition: EBConductivityOp.H:52
void diagonalScale(LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted)
virtual void AMROperatorNF(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC)
void defineStencils()
do not call this one unless you really know what you are doing
Real dx() const
returns m_dx, such function is required by some LinearSolvers
Definition: EBConductivityOp.H:151
virtual void AMRProlong(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
LayoutData< CFIVS > m_loCFIVS[SpaceDim]
Definition: EBConductivityOp.H:644
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
Definition: EBConductivityOp.H:614
one dimensional dynamic array
Definition: Vector.H:53
virtual LayoutData< BaseIVFAB< VoFStencil > > * getFluxStencil(int ivar)=0
Definition: FaceIndex.H:28
Data that maintains a one-to-one mapping of T to the boxes in a BoxLayout.
Definition: BoxLayout.H:26
LayoutData< CFIVS > m_hiCFIVS[SpaceDim]
Definition: EBConductivityOp.H:645
Definition: EBISBox.H:46
void writeEBAMRname(const Vector< LevelData< EBCellFAB > * > *a_dataPtr, const char *a_filename)
LayoutData< VoFIterator > m_vofIterIrreg
Definition: EBConductivityOp.H:636
virtual void createCoarser(LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted)
void getOpFaceStencil(VoFStencil &a_stencil, const Vector< VolIndex > &a_allMonotoneVofs, const EBISBox &a_ebisbox, const VolIndex &a_vof, int a_dir, const Side::LoHiSide &a_side, const FaceIndex &a_face, const bool &a_lowOrder)
Definition: EBLevelGrid.H:30
void getFluxEBCF(EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_ghostedBox, Vector< FaceIndex > &a_faceitEBCF, Vector< VoFStencil > &a_stenEBCF)
EBFastFR m_fastFR
Definition: EBConductivityOp.H:648
virtual void dumpAMR(Vector< LevelData< EBCellFAB > *> &a_data, string name)
for debugging
Definition: EBConductivityOp.H:68
bool m_hasCoar
Definition: EBConductivityOp.H:620
void getOpVoFStencil(VoFStencil &a_stencil, const EBISBox &a_ebisbox, const VolIndex &a_vof)
virtual void applyOpNoBoundary(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
virtual function called by LevelTGA
void writeEBLevelname(const LevelData< EBCellFAB > *a_dataPtr, const char *a_filename)
void getEBBCFluxStencil(LayoutData< BaseIVFAB< VoFStencil > > const *&a_ebbcFluxStencil)
Definition: EBConductivityOp.H:500
virtual void relax(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
Definition: EBFaceFAB.H:28
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_bcoIrreg
Definition: EBConductivityOp.H:607
const RefCountedPtr< LevelData< EBCellFAB > > getAScalingCoefficients()
Definition: EBConductivityOp.H:490
Piecewise constant interpolation.
Definition: EBMGInterp.H:33
const int SpaceDim
Definition: SPACE.H:38
void getAlphaDiagWeight(LayoutData< BaseIVFAB< Real > > const *&a_alphaDiagWeight)
Definition: EBConductivityOp.H:475
Real m_dxCoar
Definition: EBConductivityOp.H:595
void divideByIdentityCoef(LevelData< EBCellFAB > &a_rhs)
Definition: AMRMultiGrid.H:39
EBMGInterp m_ebInterpMG
Definition: EBConductivityOp.H:657
DisjointBoxLayout m_dblCoarMG
Definition: EBConductivityOp.H:658
VoF-centered stencil.
Definition: Stencils.H:60
EBLevelGrid m_eblg
Definition: EBConductivityOp.H:582
virtual void createCoarsened(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat)
A EBFaceFAB-like container for edge-centered fluxes.
Definition: EBFluxFAB.H:25
virtual void scale(LevelData< EBCellFAB > &a_lhs, const Real &a_scale)
static IntVect s_ivDebug
Definition: EBConductivityOp.H:538
virtual void applyOp(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoarse, const bool &a_homogeneousPhysBC, const bool &a_homogeneousCFBC)
ProblemDomain m_domainCoarMG
Definition: EBConductivityOp.H:660
virtual int refToCoarser()
void getAlphaBeta(Real &a_alpha, Real &a_beta)
Definition: EBConductivityOp.H:480
LayoutData< VoFIterator > m_vofItIrregColorDomLo[EBAMRPOVC_NUMSTEN][SpaceDim]
Definition: EBConductivityOp.H:568
EBLevelGrid m_eblgCoarsenedFine
Definition: EBConductivityOp.H:586
const RefCountedPtr< BaseDomainBC > getDomainBC()
Definition: EBConductivityOp.H:485
void AMROperatorNC(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
virtual void residual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC=false)
EBMGInterp m_ebInterp
Definition: EBConductivityOp.H:625
~EBConductivityOp()
Destructor.
EBLevelGrid m_eblgCoar
Definition: EBConductivityOp.H:584
const char * name(const FArrayBox &a_dummySpecializationArg)
Definition: CH_HDF5.H:907
virtual int refToFiner()
Real m_alpha
Definition: EBConductivityOp.H:609
virtual void prolongIncrement(LevelData< EBCellFAB > &a_phiThisLevel, const LevelData< EBCellFAB > &a_correctCoarse)
Definition: EBCellFAB.H:29
virtual void calculateRelaxationCoefficient()
virtual void AMRResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void applyOpIrregular(EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const DataIndex &a_datInd)
void reflux(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
bool m_hasEBCF
Definition: EBConductivityOp.H:617
void dumpLevelPoint(const LevelData< EBCellFAB > &a_res, const string &a_blab)
double Real
Definition: REAL.H:33
Definition: MultiGrid.H:30
RefCountedPtr< BaseDomainBC > m_domainBC
Definition: EBConductivityOp.H:590
RefCountedPtr< BaseEBBC > m_ebBC
Definition: EBConductivityOp.H:591
LayoutData< BaseIVFAB< Real > > m_cacheEBxDomainFluxLo[EBAMRPOVC_NUMSTEN][SpaceDim]
Definition: EBConductivityOp.H:570
static void setForceNoEBCF(bool a_forceNoEBCF)
Definition: EBConductivityOp.H:467
const EBLevelGrid & getEBLG() const
Definition: EBConductivityOp.H:463
virtual void relaxGauSai(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
virtual void AMROperator(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
virtual void relaxGSRBFast(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
int m_refToCoar
Definition: EBConductivityOp.H:616
void setTime(Real a_oldTime, Real a_mu, Real a_dt)
void getFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit)
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:30
bool m_hasInterpAve
Definition: EBConductivityOp.H:619
LoHiSide
Definition: LoHiSide.H:27
Real AMRNorm(const LevelData< EBCellFAB > &a_coarResid, const LevelData< EBCellFAB > &a_fineResid, const int &a_refRat, const int &a_ord)
virtual void AMRRestrict(LevelData< EBCellFAB > &a_resCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection, bool a_skip_res=false)
EBConductivityOp(const EBConductivityOp &a_opin)
Definition: EBConductivityOp.H:724
void applyCFBCs(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoarse, bool a_homogeneousCFBC)
static int s_whichComp
Definition: EBConductivityOp.H:506
static void Error(const char *const a_msg=m_nullString, int m_exitCode=CH_DEFAULT_ERROR_CODE)
Print out message to cerr and exit with the specified exit code.
virtual void AMRResidualNF(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC)
EBISLayout m_ebislCoarMG
Definition: EBConductivityOp.H:659
EBMGAverage m_ebAverage
Definition: EBConductivityOp.H:623
LayoutData< VoFIterator > m_vofItIrregColorDomHi[EBAMRPOVC_NUMSTEN][SpaceDim]
Definition: EBConductivityOp.H:569
Definition: EBConductivityOp.H:62
EBConductivityOp()
Default constructor. Creates an undefined conductivity operator.
virtual void dumpLevel(LevelData< EBCellFAB > &a_data, string name)
for debugging
Definition: EBConductivityOp.H:74
const RefCountedPtr< LevelData< EBFluxFAB > > getBScalingCoefficients()
Definition: EBConductivityOp.H:495
virtual void preCond(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
A Rectangular Domain on an Integer Lattice.
Definition: Box.H:465
void finerOperatorChanged(const MGLevelOp< LevelData< EBCellFAB > > &a_operator, int a_coarseningFactor)
void incrementFRCoar(EBFastFR &a_fluxReg, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi)
const IntVect m_ghostCellsRHS
Definition: EBConductivityOp.H:578
virtual void create(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
LayoutData< VoFIterator > m_vofIterMulti
Definition: EBConductivityOp.H:637
Definition: DataIndex.H:112
RefCountedPtr< LevelData< EBCellFAB > > m_acoef
Definition: EBConductivityOp.H:599
void gsrbColor(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_lph, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color)
void getFluxRegOnly(EBFaceFAB &a_fluxCentroid, const EBCellFAB &a_phi, const Box &a_ghostedBox, const Real &a_dx, const DataIndex &a_datInd, const int &a_idir)
virtual void setVal(LevelData< EBCellFAB > &a_lhs, const Real &a_value)
virtual void calculateAlphaWeight()
Piecewise constant interpolation.
Definition: EBMGAverage.H:31
LayoutData< Vector< VoFStencil > > m_stencilCoar[2 *SpaceDim]
Definition: EBConductivityOp.H:564
Real m_beta
Definition: EBConductivityOp.H:610
virtual void relaxPoiJac(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
LayoutData< RefCountedPtr< EBStencil > > m_colorEBStencil[EBAMRPOVC_NUMSTEN]
Definition: EBConductivityOp.H:572
int m_relaxType
Definition: EBConductivityOp.H:576
An integer Vector in SpaceDim-dimensional space.
Definition: CHArray.H:42
bool m_hasFine
Definition: EBConductivityOp.H:618
EBMGAverage m_ebAverageMG
Definition: EBConductivityOp.H:656
Definition: FArrayBox.H:45
Definition: AMRTGA.H:159
virtual void kappaScale(LevelData< EBCellFAB > &a_rhs)
for eb only. kappa weight the rhs but do not multiply by the identity coefficient ...
Real m_dxFine
Definition: EBConductivityOp.H:593
virtual Real localMaxNorm(const LevelData< EBCellFAB > &a_rhs)
void AMRResidualNC(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
virtual void assign(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
Volume of Fluid Index.
Definition: VolIndex.H:31
static bool s_turnOffBCs
Definition: EBConductivityOp.H:536
void incrementFRFine(EBFastFR &a_fluxReg, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void defineEBCFStencils()
void incrOpRegularAllDirs(Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_curDblBox, Box &a_curPhiBox, int a_nComps, BaseFab< Real > &a_curOpPhiFAB, const BaseFab< Real > &a_curPhiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit)
RefCountedPtr< LevelData< EBFluxFAB > > m_bcoef
Definition: EBConductivityOp.H:603
Definition: EBISLayout.H:39
LayoutData< VoFIterator > m_vofIterDomHi[CH_SPACEDIM]
Definition: EBConductivityOp.H:640
void incrOpRegularDir(EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const int &a_dir, const DataIndex &a_datInd)
void getDivFStencil(VoFStencil &a_vofStencil, const VolIndex &a_vof, const DataIndex &a_dit)
virtual void incr(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, Real a_scale)
const IntVect m_ghostCellsPhi
Definition: EBConductivityOp.H:577
void GSColorAllIrregular(EBCellFAB &a_phi, const EBCellFAB &a_rhs, const int &a_icolor, const DataIndex &a_dit)
static bool s_forceNoEBCF
Definition: EBConductivityOp.H:537
virtual void axby(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, Real a_a, Real a_b)
bool m_layoutChanged
Definition: EBConductivityOp.H:654
void applyHomogeneousCFBCs(LevelData< EBCellFAB > &a_phi)
void dumpFABPoint(const EBCellFAB &a_fab, const DataIndex &a_dit, const string &a_blab)
void defineColorStencils(Box a_ideBoxLo[SpaceDim], Box a_ideBoxHi[SpaceDim])
virtual void resetACoefficient(RefCountedPtr< LevelData< EBCellFAB > > &a_acoef)
Definition: EBConductivityOp.H:139
EBLevelGrid m_eblgFine
Definition: EBConductivityOp.H:583
LayoutData< Vector< FaceIndex > > m_faceitCoar[2 *SpaceDim]
Definition: EBConductivityOp.H:563
RefCountedPtr< EBQuadCFInterp > m_quadCFIWithCoar
Definition: EBConductivityOp.H:580
int m_refToFine
Definition: EBConductivityOp.H:615
virtual Real norm(const LevelData< EBCellFAB > &a_rhs, int a_ord)
virtual Real dotProduct(const LevelData< EBCellFAB > &a_1, const LevelData< EBCellFAB > &a_2)
void applyDomainFlux(Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_dblBox, int a_nComps, BaseFab< Real > &a_phiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit)