Chombo + EB  3.0
EBAMRPoissonOp.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 _EBAMRPOISSONOP_H_
12 #define _EBAMRPOISSONOP_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 
23 #include "EBIndexSpace.H"
24 #include "EBCellFAB.H"
25 #include "EBCellFactory.H"
26 #include "EBStencil.H"
27 
28 #include "EBLevelDataOps.H"
29 #include "BaseEBBC.H"
30 #include "BaseDomainBC.H"
31 #include "CFIVS.H"
32 #include "EBFastFR.H"
33 #include "EBMGAverage.H"
34 #include "EBMGInterp.H"
35 #include "PolyGeom.H"
36 #include "EBQuadCFInterp.H"
37 #include "EBLevelGrid.H"
38 #include "AMRTGA.H"
39 #include "EBAMRIO.H"
40 #include "AMRPoissonOp.H"
41 #include "CFRegion.H"
42 #include "NamespaceHeader.H"
43 
44 #ifdef CH_USE_PETSC
45 #include "petsc.h"
46 #include "petscmat.h"
47 #include "petscksp.h"
48 #include "petscviewer.h"
49 #endif
50 
51 #if CH_SPACEDIM==2
52 #define EBAMRPO_NUMSTEN 4
53 #elif CH_SPACEDIM==3
54 #define EBAMRPO_NUMSTEN 8
55 #else
57 {
58  THIS_WILL_ONLY_COMPILE_WHEN_CH_SPACEDIM_IS_2_OR_3;
59 }
60 #endif
61 
62 ///
63 /**
64  Operator to solve (alpha + beta lapl)phi = rhs. This follows the AMRLevelOp interface.
65 */
66 class EBAMRPoissonOp: public LevelTGAHelmOp<LevelData<EBCellFAB>, EBFluxFAB >
67 {
68 public:
69 
70  ///
71  virtual void dumpAMR(Vector<LevelData<EBCellFAB>*>& a_data, string name)
72  {
73  writeEBAMRname(&a_data, name.c_str());
74  }
75 
76  virtual void dumpLevel(LevelData<EBCellFAB>& a_data, string name)
77  {
78  writeEBLevelname(&a_data, name.c_str());
79  }
80 
81 
82  ///
83  /**
84  version that does not fill ebislcoar
85  */
86  static bool getCoarserLayouts(DisjointBoxLayout& a_dblCoar,
87  ProblemDomain& a_domainCoar,
88  const DisjointBoxLayout& a_dblFine,
89  const EBISLayout& a_ebislFine,
90  const ProblemDomain& a_domainFine,
91  int a_refToCoar,
92  const EBIndexSpace* a_ebisPtr,
93  int a_maxBoxSize,
94  bool& a_layoutChanged,
95  int a_testRef = 2);
96 
97  static Real staticMaxNorm(const LevelData<EBCellFAB>& a_rhs, const EBLevelGrid& a_eblg);
98 
99  //for tga to reset stuff
100  virtual void setAlphaAndBeta(const Real& a_alpha,
101  const Real& a_beta);
102 
103  //another tgaism
104  virtual void diagonalScale(LevelData<EBCellFAB>& a_rhs,
105  bool a_kappaWeighted = true);
107  {
108  //no acoef here.
109  }
110 
111  virtual void kappaScale(LevelData<EBCellFAB> & a_rhs)
112  {
113  //since this is a constant coefficient operator with a=1,
114  //this means the same thing as diagonal scale
115  diagonalScale(a_rhs);
116 
117  }
118 
119  ///returns m_dx, such function is required by some LinearSolvers
120  Real dx() const
121  {
122  return m_dx[0];
123  }
124 
125  ///a leveltgaism
126  virtual void fillGrad(const LevelData<EBCellFAB>& a_phi)
127  {;}
128 
129  ///
130  /**
131  dump stencil as matrix to stdout
132  */
134  {
135  public:
137  int& a_phase)
138  {
139  m_vof = a_vof;
140  m_phase = a_phase;
141  }
142 
144  {
145  MayDay::Error("No weak construction of StencilIndex class.");
146  }
147 
149  {
150  return *this;
151  }
152 
153  VolIndex vof() const
154  {
155  return m_vof;
156  }
157 
158  int phase() const
159  {
160  return m_phase;
161  }
162 
163  protected:
165  int m_phase;
166  };
167 
169  {
170  public:
171  bool operator() (const StencilIndex& a_s1, const StencilIndex& a_s2) const
172  {
173  /*
174  int p1 = a_s1.phase();
175  int p2 = a_s2.phase();
176  if (p1 == p2)
177  {
178  return (a_s1.vof()<a_s2.vof());
179  }
180  else
181  {
182  const IntVect& iv1 = a_s1.vof().gridIndex();
183  const IntVect& iv2 = a_s2.vof().gridIndex();
184  if (iv1 == iv2)
185  {
186  return (p1<p2);
187  }
188  else
189  {
190  return (iv1<iv2);
191  }
192  }
193  */
194  const IntVect& iv1 = a_s1.vof().gridIndex();
195  const IntVect& iv2 = a_s2.vof().gridIndex();
196  for (int idir=0; idir<SpaceDim; ++idir)
197  {
198  if (iv1[idir] != iv2[idir])
199  {
200  return (iv1[idir]<iv2[idir]);
201  }
202  }
203  int p1 = a_s1.phase();
204  int p2 = a_s2.phase();
205  if (p1 == p2)
206  {
207  int ci1 = a_s1.vof().cellIndex();
208  int ci2 = a_s2.vof().cellIndex();
209  return (ci1<ci2);
210  }
211  else
212  {
213  return (p1<p2);
214  }
215  }
216  };
217 
218  /*
219  Calculate stencil weights for each VoF. A mapping of each VolIndex
220  an integer is also performed. The first output is this mapping, as
221  a matrix that can be read into Matlab/Octave. The second output is
222  the stencil weight for each VoF pair, in matrix form.
223  The volume fraction is an additional output of the mapping matrix.
224  The matrices are output in the format:
225 
226  M_(RES)(imap, :) = [ gridIndex[0] gridIndex[1] (gridIndex[2]) alphaWeight kappa];
227  where 'RES' is the resolution at the EBPoissonOp's level, 'imap' is
228  the integer the VoF has been mapper to, 'alphaWeight' is the
229  weighting used to get good conditioning of the matrix, and kappa is
230  the VoF's volume fraction.
231 
232  L_(RES)(imap, jmap) = value;
233  where 'RES' is as above, 'imap' is the integer that the VoF whose
234  stencil is being computed has been mapped to, 'jmap' is the integer
235  for a VoF in that stencil. 'value' is the stencil weight.
236  */
237  void dumpStencilMatrix();
238 
239  /*
240  Calculates the stencil for VoFs at the domain boundary.
241  */
242  void getDomainFluxStencil( VoFStencil& a_stencil,
243  const VolIndex& a_vof,
244  const int a_comp,
245  const DataIndex& a_dit);
246  /*
247  A method for testing EBAMRPoissonOp::dumpStencilMatrix. Since it
248  uses calls to 'applyOp' to get the stencil weights, it is slow
249  but correct.
250  */
252 
253  virtual void getFlux(EBFluxFAB& a_flux,
254  const LevelData<EBCellFAB>& a_data,
255  const Box& a_grid,
256  const DataIndex& a_dit,
257  Real a_scale)
258  {
259  for (int idir = 0; idir < SpaceDim; idir++)
260  {
261  Box ghostedBox = a_grid;
262  ghostedBox.grow(1);
263  ghostedBox.grow(idir,-1);
264  ghostedBox &= m_eblg.getDomain();
265 
266  getFlux(a_flux[idir], a_data[a_dit], ghostedBox, a_grid,
267  m_eblg.getDomain(),
268  m_eblg.getEBISL()[a_dit], m_dx, idir);
269  }
270  }
272  {
273  return m_eblg;
274  }
276  {
277  return m_eblgCoarMG;
278  }
279  static void setOperatorTime(Real a_time)
280  {
281  s_time = a_time;
282  }
283 
284  //conforms to AMRTGA
285  virtual void setTime(Real a_time)
286  {
287  setOperatorTime(a_time);
288  }
289  ///
290  virtual ~EBAMRPoissonOp();
291 
292  ///
293  EBAMRPoissonOp();
294 
295  ///
296  /** a_residual = a_rhs - L(a_phiFine, a_phi) no coaser AMR level*/
297  void AMRResidualNC(LevelData<EBCellFAB>& a_residual,
298  const LevelData<EBCellFAB>& a_phiFine,
299  const LevelData<EBCellFAB>& a_phi,
300  const LevelData<EBCellFAB>& a_rhs,
301  bool a_homogeneousBC,
302  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
303 
304 
305  ///
306  /** apply AMR operator no coaser AMR level*/
307  void AMROperatorNC(LevelData<EBCellFAB>& a_LofPhi,
308  const LevelData<EBCellFAB>& a_phiFine,
309  const LevelData<EBCellFAB>& a_phi,
310  bool a_homogeneousBC,
311  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
312 
313  ///
314  /**
315  If you are approaching this operator from this interface, consider backing away and using
316  EBAMRPoissonOpFactory to generate these objects. Really.
317  a_eblgFine, : grid at finer level \\
318  a_eblg, : grid at this level \\
319  a_eblgCoar, : grid at coarser level \\
320  a_eblgCoarMG, : grid at intermediate multigrid level \\
321  a_domainBC, : domain boundary conditions at this level \\
322  a_ebBC: eb boundary conditions at this level \\
323  a_dx: grid spacing at this level \\
324  a_origin: offset to lowest corner of the domain \\
325  a_refToFine: refinement ratio to finer level \\
326  a_refToCoar: refinement ratio to coarser level \\
327  a_hasFiner: true if there is a finer AMR level, false otherwise. \\
328  a_hasCoarser: true if there is a coarser AMR level. \\
329  a_hasCoarserMG: true if there is a coarser MultiGrid level. \\
330  a_preCondIters: number of iterations to do for pre-conditioning \\
331  a_relaxType: 0 means point Jacobi, 1 is Gauss-Seidel. \\
332  a_alpha: coefficent of identity \\
333  a_beta: coefficient of laplacian.\\
334  a_ghostCellsPhi: Number of ghost cells in phi, correction\\
335  a_ghostCellsRhs: Number of ghost cells in RHS, residual, lphi\\
336  Ghost cell arguments are there for caching reasons. Once you set them, an error is thrown if
337  you send in data that does not match.
338  */
339 EBAMRPoissonOp(const EBLevelGrid & a_eblgFine,
340  const EBLevelGrid & a_eblg,
341  const EBLevelGrid & a_eblgCoar,
342  const EBLevelGrid & a_eblgCoarMG,
343  const RefCountedPtr<EBQuadCFInterp>& a_quadCFI,
344  const RefCountedPtr<BaseDomainBC>& a_domainBC,
345  const RefCountedPtr<BaseEBBC>& a_ebBC,
346  const RealVect& a_dx,
347  const RealVect& a_dxCoar,
348  const RealVect& a_origin,
349  const int& a_refToFine,
350  const int& a_refToCoar,
351  const bool& a_hasFine,
352  const bool& a_hasCoar,
353  const bool& a_hasMGObjects,
354  const bool& a_layoutChanged,
355  const int& a_numPreCondIters,
356  const int& a_relaxType,
357  const Real& a_alpha,
358  const Real& a_beta,
359  const IntVect& a_ghostCellsPhi,
360  const IntVect& a_ghostCellsRHS,
361  int a_testRef = 2);
362 
363  //MGOp operations. no finer or coarser
364 
365  ///
366  /**
367  */
368  virtual void residual(LevelData<EBCellFAB>& a_residual,
369  const LevelData<EBCellFAB>& a_phi,
370  const LevelData<EBCellFAB>& a_rhs,
371  bool a_homogeneousPhysBC=false);
372 
373  ///
374  /**
375  For debugging purpose, get the matrix of the op (could be change from time to time due to EBBC update) by calculating 0 - L(phi(i)) for i=0,1,...,n, each time phi(i)=1 for i-th cell and 0 for all other cells
376  */
377 
378  virtual void getOpMatrix(const LevelData<EBCellFAB>& a_phi,
379  const LevelData<EBCellFAB>& a_rhs);
380 
381  ///
382  /**
383  */
384  virtual void preCond(LevelData<EBCellFAB>& a_opPhi,
385  const LevelData<EBCellFAB>& a_phi);
386 
387  ///
388  /**
389  This function assumes that coarse-fine boundary condtions have
390  been dealt with.
391  */
392  virtual void applyOp(LevelData<EBCellFAB>& a_opPhi,
393  const LevelData<EBCellFAB>& a_phi,
394  const LevelData<EBCellFAB>* const a_phiCoarse,
395  const bool& a_homogeneousPhysBC,
396  const bool& a_homogeneousCFBC);
397 
398  /// virtual function called by LevelTGA
399  virtual void applyOpNoBoundary(LevelData<EBCellFAB>& a_opPhi,
400  const LevelData<EBCellFAB>& a_phi);
401 
402  ///
403  /**
404  get EB flux from leveldata if we have one. otherwise use m_ebBC
405  ignore input data in the case of homogeneous Phys BC
406  */
407  virtual void
408  applyOp(LevelData<EBCellFAB>& a_opPhi,
409  const LevelData<EBCellFAB>& a_phi,
410  const LevelData<EBCellFAB>* const a_phiCoar,
411  const bool& a_homogeneousPhysBC,
412  const bool& a_homogeneousCFBC,
413  const LevelData<BaseIVFAB<Real> >* const a_ebFluxBCLD);
414 
415  virtual void
416  applyOp(LevelData<EBCellFAB>& a_opPhi,
417  const LevelData<EBCellFAB>& a_phi,
418  const LevelData<EBCellFAB>* const a_phiCoar,
419  DataIterator& a_dit,
420  const bool& a_homogeneousPhysBC,
421  const bool& a_homogeneousCFBC,
422  const LevelData<BaseIVFAB<Real> >* const a_ebFluxBCLD);
423 
424  virtual void
426  const BaseFab<Real>& a_rhs,
427  const int& a_icolor,
428  const Real& a_weight,
429  const bool& a_homogeneousPhysBC,
430  const DataIndex& a_dit);
431 
432  virtual void
434  const EBCellFAB& a_rhs,
435  const int& a_icolor,
436  const bool& a_homogeneousPhysBC,
437  const DataIndex& a_dit);
438 
439  virtual void
441  const LevelData<EBCellFAB>& a_phiOld,
442  const LevelData<EBCellFAB>& a_rhs,
443  const int& a_icolor,
444  const Real& a_weight,
445  const bool& a_homogeneousPhysBC);
446 
447  virtual void
449  const LevelData<EBCellFAB>& a_phiOld,
450  const LevelData<EBCellFAB>& a_rhs,
451  const int& a_icolor,
452  const bool& a_homogeneousPhysBC);
453 
454  ///
455  /**
456  this is the linearop function. CFBC is set to homogeneous. phic is null
457  */
458  virtual void applyOp(LevelData<EBCellFAB>& a_opPhi,
459  const LevelData<EBCellFAB>& a_phi,
460  bool a_homogeneousPhysBC);
461 
462  /// no exchange of cf interp
463  void
465  const LevelData<EBCellFAB>& a_phi,
466  const LevelData<EBCellFAB>* const a_phiCoar,
467  DataIterator& a_dit,
468  const bool& a_homogeneousPhysBC,
469  const bool& a_homogeneousCFBC,
470  const LevelData<BaseIVFAB<Real> >* const a_ebFluxBCLD //only non null in multifluid
471  );
472 
473  ///
474  /**
475  */
476  virtual void create(LevelData<EBCellFAB>& a_lhs,
477  const LevelData<EBCellFAB>& a_rhs);
478 
479  ///
480  virtual void createCoarsened(LevelData<EBCellFAB>& a_lhs,
481  const LevelData<EBCellFAB>& a_rhs,
482  const int& a_refRat);
483 
484  Real
485  AMRNorm(const LevelData<EBCellFAB>& a_coarResid,
486  const LevelData<EBCellFAB>& a_fineResid,
487  const int& a_refRat,
488  const int& a_ord);
489 
490  ///
491  /**
492  */
493  virtual void assign(LevelData<EBCellFAB>& a_lhs,
494  const LevelData<EBCellFAB>& a_rhs);
495 
496  ///
497  /**
498  */
499  virtual Real dotProduct(const LevelData<EBCellFAB>& a_1,
500  const LevelData<EBCellFAB>& a_2);
501 
502  ///
503  /**
504  */
505  virtual void incr(LevelData<EBCellFAB>& a_lhs,
506  const LevelData<EBCellFAB>& a_x,
507  Real a_scale);
508 
509  ///
510  /**
511  */
512  virtual void axby(LevelData<EBCellFAB>& a_lhs,
513  const LevelData<EBCellFAB>& a_x,
514  const LevelData<EBCellFAB>& a_y,
515  Real a_a,
516  Real a_b);
517 
518  ///
519  /**
520  */
521  virtual void scale(LevelData<EBCellFAB>& a_lhs,
522  const Real& a_scale);
523 
524  ///
525  /**
526  */
527  virtual Real norm(const LevelData<EBCellFAB>& a_rhs,
528  int a_ord);
529 
530  ///
531  /**
532  */
533  virtual Real localMaxNorm(const LevelData<EBCellFAB>& a_rhs);
534 
535  ///
536  /**
537  */
538  virtual void setToZero(LevelData<EBCellFAB>& a_lhs);
539 
540  ///
541  /**
542  */
543  virtual void setVal(LevelData<EBCellFAB>& a_lhs, const Real& a_value);
544 
545  ///
546  /**
547  */
548  virtual void createCoarser(LevelData<EBCellFAB>& a_coarse,
549  const LevelData<EBCellFAB>& a_fine,
550  bool a_ghosted);
551 
552  ///
553  /**
554  */
555  virtual void relax(LevelData<EBCellFAB>& a_e,
556  const LevelData<EBCellFAB>& a_residual,
557  int a_iterations);
558 
559  ///
560  /**
561  Calculate restricted residual:
562  a_resCoarse[2h] = I[h->2h] (a_rhsFine[h] - L[h](a_phiFine[h]))
563  */
564  virtual void restrictResidual(LevelData<EBCellFAB>& a_resCoarse,
565  LevelData<EBCellFAB>& a_phiFine,
566  const LevelData<EBCellFAB>& a_rhsFine);
567 
568  ///
569  /**
570  Correct the fine solution based on coarse correction:
571  a_phiThisLevel += I[2h->h] (a_correctCoarse)
572  */
573  virtual void prolongIncrement(LevelData<EBCellFAB>& a_phiThisLevel,
574  const LevelData<EBCellFAB>& a_correctCoarse);
575 
576  ///
577  /** Refinement ratio between this level and coarser level.
578  Returns 1 when there are no coarser AMRLevelOp objects */
579  virtual int refToCoarser();
580 
581  ///
582  /** Refinement ratio between this level and coarser level.
583  Returns 1 when there are no coarser AMRLevelOp objects */
584  virtual int refToFiner();
585 
586  ///
587  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
588  virtual void AMRResidual(LevelData<EBCellFAB>& a_residual,
589  const LevelData<EBCellFAB>& a_phiFine,
590  const LevelData<EBCellFAB>& a_phi,
591  const LevelData<EBCellFAB>& a_phiCoarse,
592  const LevelData<EBCellFAB>& a_rhs,
593  bool a_homogeneousBC,
594  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
595 
596  ///
597  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
598  virtual void AMRResidualNF(LevelData<EBCellFAB>& a_residual,
599  const LevelData<EBCellFAB>& a_phi,
600  const LevelData<EBCellFAB>& a_phiCoarse,
601  const LevelData<EBCellFAB>& a_rhs,
602  bool a_homogeneousBC);
603 
604 
605  ///
606  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
607  virtual void AMROperator(LevelData<EBCellFAB>& a_LofPhi,
608  const LevelData<EBCellFAB>& a_phiFine,
609  const LevelData<EBCellFAB>& a_phi,
610  const LevelData<EBCellFAB>& a_phiCoarse,
611  bool a_homogeneousBC,
612  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
613 
614  ///
615  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
616  virtual void AMROperatorNF(LevelData<EBCellFAB>& a_LofPhi,
617  const LevelData<EBCellFAB>& a_phi,
618  const LevelData<EBCellFAB>& a_phiCoarse,
619  bool a_homogeneousBC);
620 
621  ///
622  /** a_resCoarse = I[h-2h] (a_residual - L(a_correction, a_coarseCorrection)) */
623  virtual void AMRRestrict(LevelData<EBCellFAB>& a_resCoarse,
624  const LevelData<EBCellFAB>& a_residual,
625  const LevelData<EBCellFAB>& a_correction,
626  const LevelData<EBCellFAB>& a_coarseCorrection);
627 
628  ///
629  /** a_correction += I[2h->h](a_coarseCorrection) */
630  virtual void AMRProlong(LevelData<EBCellFAB>& a_correction,
631  const LevelData<EBCellFAB>& a_coarseCorrection);
632 
633  ///
634  /** a_residual = a_residual - L(a_correction, a_coarseCorrection) */
635  virtual void AMRUpdateResidual(LevelData<EBCellFAB>& a_residual,
636  const LevelData<EBCellFAB>& a_correction,
637  const LevelData<EBCellFAB>& a_coarseCorrection);
638 
639  ///
640  /** a_residual = a_residual - L(a_correction, a_coarseCorrection)
641  used in multifluid */
642  void AMRUpdateResidual(LevelData<EBCellFAB>& a_residual,
643  const LevelData<EBCellFAB>& a_correction,
644  const LevelData<EBCellFAB>& a_coarseCorrection,
645  const LevelData<BaseIVFAB<Real> >* const a_ebFluxBCLD);
646 
647 
648  void reflux(LevelData<EBCellFAB>& a_residual,
649  const LevelData<EBCellFAB>& a_phiFine,
650  const LevelData<EBCellFAB>& a_phi,
651  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
652 
653 
654  void fast_reflux(LevelData<EBCellFAB>& a_residual,
655  const LevelData<EBCellFAB>& a_phiFine,
656  const LevelData<EBCellFAB>& a_phi,
657  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
658 
659  void setEBBC(const RefCountedPtr<BaseEBBC>& a_ebBC);
660 
661  //for debugging other operators
663  const LevelData<EBCellFAB>& a_lph,
664  const LevelData<EBCellFAB>& a_rhs,
665  const IntVect& a_color);
666 
667  //for debugging other operators
669  const LevelData<EBCellFAB>& a_rhs);
670 
672  const LevelData<EBCellFAB>& a_rhs);
673 
675  const LevelData<EBCellFAB>& a_resid,
676  const IntVect& color);
677 
678  void colorGS(LevelData<EBCellFAB>& a_phi,
679  const LevelData<EBCellFAB>& a_rhs,
680  const int& a_icolor);
681 
682  void levelGSRB(LevelData<EBCellFAB>& a_phi,
683  const LevelData<EBCellFAB>& a_rhs);
684 
685  void levelGSRB(LevelData<EBCellFAB>& a_phi,
686  const LevelData<EBCellFAB>& a_rhs,
687  const int a_color);
688 
690  const LevelData<EBCellFAB>& a_rhs);
691 
693  const LevelData<EBCellFAB>& a_phiOld,
694  const LevelData<EBCellFAB>& a_rhs,
695  const int& a_icolor);
696 
697  static void doLazyRelax(bool a_doLazyRelax);
698  static void doEBEllipticLoadBalance(bool a_doEBEllipticLoadBalance);
699  static void areaFracWeighted(bool a_areaFracWeighted);
700 
701  static void getDivFStencil(VoFStencil& a_vofStencil,
702  const VolIndex& a_vof,
703  const EBISBox& a_ebisBox,
704  const IntVectSet& a_cfivs,
705  const RealVect& a_dx,
706  bool doFaceInterp = true);
707 
708  void getDivFStencil(VoFStencil& a_vofStencil,
709  const VolIndex& a_vof,
710  const DataIndex& a_dit,
711  bool doFaceInterp);
712 
713  //CP recently added
714  void getVoFStencil(LayoutData<BaseIVFAB<VoFStencil> > const*& a_vofStencil)
715  {
716  //I'm defining const* to make sure the content in m_opStencil is not modified
717  a_vofStencil = &m_opStencil;
718  }
719 
720  void getAlphaDiagWeight(LayoutData<BaseIVFAB<Real> > const*& a_alphaDiagWeight)
721  {
722  a_alphaDiagWeight = &m_alphaDiagWeight;
723  }
724 
725  void getAlphaBeta(Real& a_alpha, Real& a_beta)
726  {
727  a_alpha = m_alpha;
728  a_beta = m_beta;
729  }
730  //CP recently added
732  {
733  //I'm defining const* to make sure the content in m_opStencil is not modified
734  return m_domainBC;
735  }
736 
737  static void getFluxStencil(VoFStencil& a_fluxStencil,
738  const FaceIndex& a_face,
739  const EBISBox& a_ebisBox,
740  const IntVectSet& a_cfivs,
741  const RealVect& a_dx,
742  bool a_doFaceInterp);
743 
744  static void getFaceCenteredFluxStencil(VoFStencil& a_fluxStencil,
745  const FaceIndex& a_face,
746  const RealVect& a_dx);
747 
748  void applyCFBCs(LevelData<EBCellFAB>& a_phi,
749  const LevelData<EBCellFAB>* const a_phiCoarse,
750  bool a_homogeneousCFBC,
751  bool a_doOnlyRegularInterp = false);
752 
753  ///
754  /**
755  This function computes: a_lhs = (1/diagonal)*a_rhs
756  It is used to initialize the preconditioner, and by
757  MFPoissonOp::levelJacobi.
758  Consider using one of the level Gauss-Seidel methods
759  instead of monkeying with this.
760  */
762  const LevelData<EBCellFAB>& a_rhs);
763 
764 protected:
765 
766  void defineMGObjects(const EBLevelGrid& a_eblgCoarMG);
767  void defineWithCoarser(const EBLevelGrid& a_eblgCoar, const int& a_refToCoar);
768  void defineWithFiner(const EBLevelGrid& a_eblgFine,
769  const int& a_refToFine);
770 
771  static bool s_turnOffBCs;
773  static bool s_areaFracWeighted;
774  void defineStencils();
775  void defineEBCFStencils();
776 
778 
781 
783 
789 
792 
796 
804  static Real s_time;
805  static bool s_doLazyRelax;
807  static bool s_doTrimEdges;
811  bool m_hasFine;
813  bool m_hasCoar;
816 
817  bool m_hasEBCF;
818 
820  //restriction object
822  //prolongation object
824 
825  //stencils for operator evaluation
830  //stencils for operator evaluation on gauss-seidel colors
835  //weights that get multiplied by alpha
837  //weights that get multiplied by beta
839  //constant for using in place of weight in EBStencil->apply for inhom dom bcs
841 
842  //cache the irreg vofiterator
845 
848 
853  // LayoutData<BaseIVFAB<Real> > m_cacheIrregFluxDomLo[EBAMRPO_NUMSTEN][SpaceDim];
854  // LayoutData<BaseIVFAB<Real> > m_cacheIrregFluxDomHi[EBAMRPO_NUMSTEN][SpaceDim];
855 
856  // Coarse-fine stencils for homogeneous CFInterp
859 
860  //flux register with finer level
862 
863  //stuff to make EBCF go faster
866 
867 
868  //special mg objects for when we do not have
869  //a coarser level or when the refinement to coarse
870  //is greater than two
871  //flag for when we need special MG objects
874 
876 
877 
878  //stuff below is only defined if m_hasMGObjects==true
884 
885 private:
886  //internally useful but not for general consumption
887  //lots of hidden assumptions and all that
888 
889  void fast_incrementFRCoar(const LevelData<EBCellFAB>& a_phiFine,
890  const LevelData<EBCellFAB>& a_phi);
891 
892  void fast_incrementFRFine(const LevelData<EBCellFAB>& a_phiFine,
893  const LevelData<EBCellFAB>& a_phi,
894  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
895 
896  ///this one gets called by base level tga (called by the public getFlux)
897  void getFlux(EBFaceFAB& a_flux,
898  const EBCellFAB& a_phi,
899  const Box& a_ghostedBox,
900  const Box& a_fabBox,
901  const ProblemDomain& a_domainBox,
902  const EBISBox& a_ebisBox,
903  const RealVect& a_dx,
904  const int& a_idir);
905 
906  ///this one is internal (called by refluxing)
907  void getFluxEBCF(EBFaceFAB& a_flux,
908  const EBCellFAB& a_phi,
909  const Box& a_ghostedBox,
910  Vector<FaceIndex>& a_faceitEBCF,
911  Vector<VoFStencil>& a_ebcfsten,
912  const RealVect& a_dx);
913 
914  ///this one is internal (called by refluxing)
915  void getFluxRegO(EBFaceFAB& a_flux,
916  const EBCellFAB& a_phi,
917  const Box& a_ghostedBox,
918  const RealVect& a_dx);
919 
920  void getOpVoFStencil(VoFStencil& a_stencil,
921  const EBISBox& a_ebisbox,
922  const VolIndex& a_vof);
923 
924  void getOpVoFStencil(VoFStencil& a_stencil,
925  const int& a_dir,
926  const Vector<VolIndex>& a_allMonotoneVoFs,
927  const EBISBox& a_ebisbox,
928  const VolIndex& a_vof,
929  const bool& a_lowOrder);
930 
931 
932  void getOpFaceStencil(VoFStencil& a_stencil,
933  const Vector<VolIndex>& a_allMonotoneVofs,
934  const EBISBox& a_ebisbox,
935  const VolIndex& a_vof,
936  int a_dir,
937  const Side::LoHiSide& a_side,
938  const FaceIndex& a_face,
939  const bool& a_lowOrder);
940 
941  void levelJacobi(LevelData<EBCellFAB>& a_phi,
942  const LevelData<EBCellFAB>& a_rhs);
943 
945 
946  void applyHomogeneousCFBCs(EBCellFAB& a_phi,
947  const DataIndex& a_datInd,
948  int a_idir,
949  Side::LoHiSide a_hiorlo);
950 
951  Real getRadius(const FaceIndex& a_face, const RealVect& a_centroid);
952 
953  ///
954  /** applyOp() on the regular cells for all directions
955  opPhi comes in holding alpha*phi. this adds on beta*lapl phi*/
956  void applyOpRegularAllDirs(Box * a_loBox,
957  Box * a_hiBox,
958  int * a_hasLo,
959  int * a_hasHi,
960  Box & a_curOpPhiBox,
961  Box & a_curPhiBox,
962  int a_nComps,
963  BaseFab<Real> & a_curOpPhiFAB,
964  const BaseFab<Real> & a_curPhiFAB,
965  bool a_homogeneousPhysBC,
966  const DataIndex& a_dit,
967  const Real& a_beta);
968 
969  void applyDomainFlux(Box * a_loBox,
970  Box * a_hiBox,
971  int * a_hasLo,
972  int * a_hasHi,
973  Box & a_curPhiBox,
974  int a_nComps,
975  BaseFab<Real> & a_phiFAB,
976  bool a_homogeneousPhysBC,
977  const DataIndex& a_dit,
978  const Real& a_beta);
979 
980 private:
981  //copy constructor and operator= disallowed for all the usual reasons
983  {
984  MayDay::Error("invalid operator");
985  }
986 
987  void operator=(const EBAMRPoissonOp& a_opin)
988  {
989  MayDay::Error("invalid operator");
990  }
991 };
992 
993 
994 #include "NamespaceFooter.H"
995 #endif
bool m_hasCoar
Definition: EBAMRPoissonOp.H:813
EBMGAverage m_ebAverage
Definition: EBAMRPoissonOp.H:821
Vector< IntVect > m_colors
Definition: EBAMRPoissonOp.H:875
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)
LayoutData< VoFIterator > m_vofItIrregDomHi[SpaceDim]
Definition: EBAMRPoissonOp.H:847
virtual void relax(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
LayoutData< EBCellFAB > m_cacheInhomDomBCLo[SpaceDim]
Definition: EBAMRPoissonOp.H:851
int m_testRef
Definition: EBAMRPoissonOp.H:777
virtual void create(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
RealVect m_dxFine
Definition: EBAMRPoissonOp.H:793
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)
static void setOperatorTime(Real a_time)
Definition: EBAMRPoissonOp.H:279
EBFastFR-A class to encapsulate a levels worth of flux registers.
Definition: EBFastFR.H:39
RealVect m_origin
Definition: EBAMRPoissonOp.H:808
EBISLayout getEBISL() const
Definition: EBLevelGrid.H:93
An irregular domain on an integer lattice.
Definition: IntVectSet.H:44
virtual void prolongIncrement(LevelData< EBCellFAB > &a_phiThisLevel, const LevelData< EBCellFAB > &a_correctCoarse)
bool m_hasFine
Definition: EBAMRPoissonOp.H:811
static Real staticMaxNorm(const LevelData< EBCellFAB > &a_rhs, const EBLevelGrid &a_eblg)
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:130
virtual void setToZero(LevelData< EBCellFAB > &a_lhs)
virtual void createCoarsened(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat)
void defineWithFiner(const EBLevelGrid &a_eblgFine, const int &a_refToFine)
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
Definition: EBAMRPoissonOp.H:838
Definition: EBIndexSpace.H:260
virtual void kappaScale(LevelData< EBCellFAB > &a_rhs)
for eb only. kappa weight the rhs but do not multiply by the identity coefficient ...
Definition: EBAMRPoissonOp.H:111
void getDomainFluxStencil(VoFStencil &a_stencil, const VolIndex &a_vof, const int a_comp, const DataIndex &a_dit)
EBFastFR m_fastFR
Definition: EBAMRPoissonOp.H:861
LayoutData< RefCountedPtr< EBStencil > > m_opEBStencilInhomDomHi[SpaceDim]
Definition: EBAMRPoissonOp.H:829
void defineWithCoarser(const EBLevelGrid &a_eblgCoar, const int &a_refToCoar)
Real getRadius(const FaceIndex &a_face, const RealVect &a_centroid)
virtual void GSColorAllIrregular(EBCellFAB &a_phi, const EBCellFAB &a_rhs, const int &a_icolor, const bool &a_homogeneousPhysBC, const DataIndex &a_dit)
void getFluxEBCF(EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_ghostedBox, Vector< FaceIndex > &a_faceitEBCF, Vector< VoFStencil > &a_ebcfsten, const RealVect &a_dx)
this one is internal (called by refluxing)
one dimensional dynamic array
Definition: Vector.H:52
virtual Real localMaxNorm(const LevelData< EBCellFAB > &a_rhs)
VolIndex m_vof
Definition: EBAMRPoissonOp.H:164
Definition: FaceIndex.H:28
Data that maintains a one-to-one mapping of T to the boxes in a BoxLayout.
Definition: LayoutData.H:46
const ProblemDomain & getDomain() const
Definition: EBLevelGrid.H:132
virtual void restrictResidual(LevelData< EBCellFAB > &a_resCoarse, LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_rhsFine)
Definition: EBAMRPoissonOp.H:168
int m_refToCoar
Definition: EBAMRPoissonOp.H:810
void applyOpNoCFBCs(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoar, DataIterator &a_dit, const bool &a_homogeneousPhysBC, const bool &a_homogeneousCFBC, const LevelData< BaseIVFAB< Real > > *const a_ebFluxBCLD)
no exchange of cf interp
A strange but true thing to make copying from one boxlayoutdata to another fast.
Definition: Copier.H:137
const IntVect m_ghostCellsRHS
Definition: EBAMRPoissonOp.H:780
LayoutData< RefCountedPtr< EBStencil > > m_opEBStencilInhomDomLo[SpaceDim]
Definition: EBAMRPoissonOp.H:828
Definition: EBISBox.H:46
static void doEBEllipticLoadBalance(bool a_doEBEllipticLoadBalance)
void writeEBAMRname(const Vector< LevelData< EBCellFAB > * > *a_dataPtr, const char *a_filename)
const IntVect & gridIndex() const
Definition: VolIndex.H:119
LayoutData< CFIVS > m_loCFIVS[SpaceDim]
Definition: EBAMRPoissonOp.H:857
LayoutData< BaseIVFAB< VoFStencil > > m_opStencil
Definition: EBAMRPoissonOp.H:826
Definition: EBLevelGrid.H:30
virtual Real dotProduct(const LevelData< EBCellFAB > &a_1, const LevelData< EBCellFAB > &a_2)
LayoutData< VoFIterator > m_vofItIrregColor[EBAMRPO_NUMSTEN]
Definition: EBAMRPoissonOp.H:844
EBISLayout m_ebislCoarMG
Definition: EBAMRPoissonOp.H:882
EBLevelGrid getEBLG()
Definition: EBAMRPoissonOp.H:271
static bool s_doInconsistentRelax
Definition: EBAMRPoissonOp.H:806
Definition: DataIterator.H:140
virtual void dumpAMR(Vector< LevelData< EBCellFAB > *> &a_data, string name)
Definition: EBAMRPoissonOp.H:71
LayoutData< RefCountedPtr< EBStencil > > m_opEBStencil
Definition: EBAMRPoissonOp.H:827
void writeEBLevelname(const LevelData< EBCellFAB > *a_dataPtr, const char *a_filename)
static bool s_turnOffBCs
Definition: EBAMRPoissonOp.H:771
int m_relaxType
Definition: EBAMRPoissonOp.H:815
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)
Definition: EBFaceFAB.H:28
EBMGInterp m_ebInterp
Definition: EBAMRPoissonOp.H:823
void reflux(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void levelMultiColorGS(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
LayoutData< RefCountedPtr< EBStencil > > m_colorEBStencilDomHi[EBAMRPO_NUMSTEN][SpaceDim]
Definition: EBAMRPoissonOp.H:833
Copier m_exchangeCopier
Definition: EBAMRPoissonOp.H:819
virtual void GSColorAllRegularClone(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiOld, const LevelData< EBCellFAB > &a_rhs, const int &a_icolor, const Real &a_weight, const bool &a_homogeneousPhysBC)
void fast_incrementFRFine(const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
Real AMRNorm(const LevelData< EBCellFAB > &a_coarResid, const LevelData< EBCellFAB > &a_fineResid, const int &a_refRat, const int &a_ord)
int m_phase
Definition: EBAMRPoissonOp.H:165
Piecewise constant interpolation.
Definition: EBMGInterp.H:33
LayoutData< RefCountedPtr< EBStencil > > m_invDiagEBStencil
Definition: EBAMRPoissonOp.H:834
void operator=(const EBAMRPoissonOp &a_opin)
Definition: EBAMRPoissonOp.H:987
int m_refToFine
Definition: EBAMRPoissonOp.H:809
const int SpaceDim
Definition: SPACE.H:39
static bool s_areaFracWeighted
Definition: EBAMRPoissonOp.H:773
virtual void scale(LevelData< EBCellFAB > &a_lhs, const Real &a_scale)
Definition: AMRMultiGrid.H:35
static Real s_time
Definition: EBAMRPoissonOp.H:804
virtual int refToCoarser()
virtual void assign(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
VoF-centered stencil.
Definition: Stencils.H:59
RealVect m_invDx
Definition: EBAMRPoissonOp.H:797
RealVect m_invDx2
Definition: EBAMRPoissonOp.H:798
RealVect m_dx
Definition: EBAMRPoissonOp.H:794
void getAlphaDiagWeight(LayoutData< BaseIVFAB< Real > > const *&a_alphaDiagWeight)
Definition: EBAMRPoissonOp.H:720
Definition: EBAMRPoissonOp.H:66
A EBFaceFAB-like container for edge-centered fluxes.
Definition: EBFluxFAB.H:25
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)
EBLevelGrid m_eblgCoarMG
Definition: EBAMRPoissonOp.H:787
virtual ~EBAMRPoissonOp()
Real m_dxScale
Definition: EBAMRPoissonOp.H:799
void getInvDiagRHS(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
static void getFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const EBISBox &a_ebisBox, const IntVectSet &a_cfivs, const RealVect &a_dx, bool a_doFaceInterp)
virtual void AMRRestrict(LevelData< EBCellFAB > &a_resCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
static bool getCoarserLayouts(DisjointBoxLayout &a_dblCoar, ProblemDomain &a_domainCoar, const DisjointBoxLayout &a_dblFine, const EBISLayout &a_ebislFine, const ProblemDomain &a_domainFine, int a_refToCoar, const EBIndexSpace *a_ebisPtr, int a_maxBoxSize, bool &a_layoutChanged, int a_testRef=2)
virtual int refToFiner()
virtual void AMRProlong(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
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)
StencilIndex()
Definition: EBAMRPoissonOp.H:143
const char * name(const FArrayBox &a_dummySpecializationArg)
Definition: CH_HDF5.H:741
void getAlphaBeta(Real &a_alpha, Real &a_beta)
Definition: EBAMRPoissonOp.H:725
Definition: EBCellFAB.H:29
Real m_beta
Definition: EBAMRPoissonOp.H:802
LayoutData< VoFIterator > m_vofItIrreg
Definition: EBAMRPoissonOp.H:843
void dumpStencilMatrix()
void fast_incrementFRCoar(const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi)
void defineEBCFStencils()
virtual void setTime(Real a_time)
Definition: EBAMRPoissonOp.H:285
EBMGInterp m_ebInterpMG
Definition: EBAMRPoissonOp.H:880
double Real
Definition: REAL.H:33
EBAMRPoissonOp(const EBAMRPoissonOp &a_opin)
Definition: EBAMRPoissonOp.H:982
static void doLazyRelax(bool a_doLazyRelax)
EBLevelGrid m_eblgCoar
Definition: EBAMRPoissonOp.H:786
bool m_hasMGObjects
Definition: EBAMRPoissonOp.H:872
void defineStencils()
Real m_bCoef
Definition: EBAMRPoissonOp.H:803
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:31
LayoutData< VoFIterator > m_vofItIrregColorDomHi[EBAMRPO_NUMSTEN][SpaceDim]
Definition: EBAMRPoissonOp.H:850
RefCountedPtr< BaseDomainBC > m_domainBC
Definition: EBAMRPoissonOp.H:790
LoHiSide
Definition: LoHiSide.H:27
static bool s_doEBEllipticLoadBalance
Definition: EBAMRPoissonOp.H:772
static void areaFracWeighted(bool a_areaFracWeighted)
void defineMGObjects(const EBLevelGrid &a_eblgCoarMG)
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 AMROperatorNF(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC)
EBLevelGrid getEBLGCoarMG()
Definition: EBAMRPoissonOp.H:275
virtual void fillGrad(const LevelData< EBCellFAB > &a_phi)
a leveltgaism
Definition: EBAMRPoissonOp.H:126
RefCountedPtr< EBQuadCFInterp > m_quadCFIWithCoar
Definition: EBAMRPoissonOp.H:782
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
Definition: EBAMRPoissonOp.H:836
bool m_hasInterpAve
Definition: EBAMRPoissonOp.H:812
static bool s_doTrimEdges
Definition: EBAMRPoissonOp.H:807
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 Real norm(const LevelData< EBCellFAB > &a_rhs, int a_ord)
virtual void applyOpNoBoundary(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
virtual function called by LevelTGA
virtual void preCond(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
LayoutData< Vector< FaceIndex > > m_faceitCoar[2 *SpaceDim]
Definition: EBAMRPoissonOp.H:864
RefCountedPtr< BaseEBBC > m_ebBC
Definition: EBAMRPoissonOp.H:791
void levelGSRB(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
virtual void setAlphaAndBeta(const Real &a_alpha, const Real &a_beta)
void setEBBC(const RefCountedPtr< BaseEBBC > &a_ebBC)
LayoutData< RefCountedPtr< EBStencil > > m_colorEBStencilDomLo[EBAMRPO_NUMSTEN][SpaceDim]
Definition: EBAMRPoissonOp.H:832
virtual void setVal(LevelData< EBCellFAB > &a_lhs, const Real &a_value)
bool m_layoutChanged
Definition: EBAMRPoissonOp.H:873
Real m_alpha
Definition: EBAMRPoissonOp.H:800
A Rectangular Domain on an Integer Lattice.
Definition: Box.H:465
A Real vector in SpaceDim-dimensional space.
Definition: RealVect.H:41
const RefCountedPtr< BaseDomainBC > getDomainBC()
Definition: EBAMRPoissonOp.H:731
void getOpVoFStencil(VoFStencil &a_stencil, const EBISBox &a_ebisbox, const VolIndex &a_vof)
void AMROperatorNC(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
LayoutData< EBCellFAB > m_cacheInhomDomBCHi[SpaceDim]
Definition: EBAMRPoissonOp.H:852
LayoutData< Vector< VoFStencil > > m_stencilCoar[2 *SpaceDim]
Definition: EBAMRPoissonOp.H:865
void fast_reflux(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void slowGSRBColor(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_lph, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color)
virtual void AMRUpdateResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
Definition: DataIndex.H:112
bool m_hasEBCF
Definition: EBAMRPoissonOp.H:817
Real dx() const
returns m_dx, such function is required by some LinearSolvers
Definition: EBAMRPoissonOp.H:120
int cellIndex() const
Definition: VolIndex.H:127
LayoutData< VoFIterator > m_vofItIrregDomLo[SpaceDim]
Definition: EBAMRPoissonOp.H:846
Piecewise constant interpolation.
Definition: EBMGAverage.H:31
Real m_aCoef
Definition: EBAMRPoissonOp.H:801
StencilIndex & operator=(const StencilIndex &a_sin)
Definition: EBAMRPoissonOp.H:148
void dumpReferenceStencilMatrix()
EBLevelGrid m_eblg
Definition: EBAMRPoissonOp.H:784
void colorGS(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, const int &a_icolor)
LayoutData< VoFIterator > m_vofItIrregColorDomLo[EBAMRPO_NUMSTEN][SpaceDim]
Definition: EBAMRPoissonOp.H:849
An integer Vector in SpaceDim-dimensional space.
Definition: CHArray.H:42
ProblemDomain m_domainCoarMG
Definition: EBAMRPoissonOp.H:883
Definition: AMRTGA.H:159
void colorGSClone(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiOld, const LevelData< EBCellFAB > &a_rhs, const int &a_icolor)
LayoutData< CFIVS > m_hiCFIVS[SpaceDim]
Definition: EBAMRPoissonOp.H:858
void THIS_IS_AN_ERROR_MESSAGE(void)
Definition: EBAMRPoissonOp.H:56
Volume of Fluid Index.
Definition: VolIndex.H:31
void applyCFBCs(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoarse, bool a_homogeneousCFBC, bool a_doOnlyRegularInterp=false)
void levelMultiColorGSClone(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
Box & grow(int i)
grow functions
Definition: Box.H:2268
EBMGAverage m_ebAverageMG
Definition: EBAMRPoissonOp.H:879
EBLevelGrid m_eblgCoarsenedFine
Definition: EBAMRPoissonOp.H:788
Definition: EBISLayout.H:39
virtual void getFlux(EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale)
Definition: EBAMRPoissonOp.H:253
virtual void diagonalScale(LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted=true)
static bool s_doLazyRelax
Definition: EBAMRPoissonOp.H:805
void applyOpRegularAllDirs(Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_curOpPhiBox, Box &a_curPhiBox, int a_nComps, BaseFab< Real > &a_curOpPhiFAB, const BaseFab< Real > &a_curPhiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit, const Real &a_beta)
void levelSlowRelax(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
virtual void residual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC=false)
StencilIndex(VolIndex &a_vof, int &a_phase)
Definition: EBAMRPoissonOp.H:136
const IntVect m_ghostCellsPhi
Definition: EBAMRPoissonOp.H:779
virtual void axby(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, Real a_a, Real a_b)
virtual void GSColorAllIrregularClone(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiOld, const LevelData< EBCellFAB > &a_rhs, const int &a_icolor, const bool &a_homogeneousPhysBC)
virtual void GSColorAllRegular(BaseFab< Real > &a_phi, const BaseFab< Real > &a_rhs, const int &a_icolor, const Real &a_weight, const bool &a_homogeneousPhysBC, const DataIndex &a_dit)
void getFluxRegO(EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_ghostedBox, const RealVect &a_dx)
this one is internal (called by refluxing)
void applyDomainFlux(Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_curPhiBox, int a_nComps, BaseFab< Real > &a_phiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit, const Real &a_beta)
virtual void dumpLevel(LevelData< EBCellFAB > &a_data, string name)
Definition: EBAMRPoissonOp.H:76
VolIndex vof() const
Definition: EBAMRPoissonOp.H:153
static void getFaceCenteredFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const RealVect &a_dx)
virtual void incr(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, Real a_scale)
void levelJacobi(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
int m_numPreCondIters
Definition: EBAMRPoissonOp.H:814
LayoutData< BaseIVFAB< Real > > m_one
Definition: EBAMRPoissonOp.H:840
virtual void getOpMatrix(const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
void getVoFStencil(LayoutData< BaseIVFAB< VoFStencil > > const *&a_vofStencil)
Definition: EBAMRPoissonOp.H:714
EBLevelGrid m_eblgFine
Definition: EBAMRPoissonOp.H:785
Definition: EBAMRPoissonOp.H:133
virtual void divideByIdentityCoef(LevelData< EBCellFAB > &a_rhs)
Definition: EBAMRPoissonOp.H:106
static void getDivFStencil(VoFStencil &a_vofStencil, const VolIndex &a_vof, const EBISBox &a_ebisBox, const IntVectSet &a_cfivs, const RealVect &a_dx, bool doFaceInterp=true)
DisjointBoxLayout m_dblCoarMG
Definition: EBAMRPoissonOp.H:881
LayoutData< RefCountedPtr< EBStencil > > m_colorEBStencil[EBAMRPO_NUMSTEN]
Definition: EBAMRPoissonOp.H:831
RealVect m_dxCoar
Definition: EBAMRPoissonOp.H:795
void applyHomogeneousCFBCs(LevelData< EBCellFAB > &a_phi)
int phase() const
Definition: EBAMRPoissonOp.H:158
virtual void createCoarser(LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted)