Chombo + EB  3.0
EBViscousTensorOp.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 _EBVISCOUSTENSOROP_H_
12 #define _EBVISCOUSTENSOROP_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 
27 #include "EBLevelDataOps.H"
28 #include "BaseEBBC.H"
29 #include "AMRTGA.H"
30 #include "BaseDomainBC.H"
31 #include "CFIVS.H"
32 #include "EBFastFR.H"
33 #include "EBMGAverage.H"
34 #include "EBMGInterp.H"
35 #include "EBTensorCFInterp.H"
36 #include "EBCoarsen.H"
37 #include "PolyGeom.H"
38 #include "EBQuadCFInterp.H"
39 #include "EBLevelGrid.H"
40 #include "EBAMRIO.H"
41 #include "DivergenceStencil.H"
42 
43 #include "NamespaceHeader.H"
44 
45 // Forward declaration.
46 class EBViscousTensorOp;
47 
48 //! \class EBViscousTensorOp
49 //! This class implements an operator that solves the equation
50 //! alpha a I + (divF) = alpha a I + beta*div(eta(grad B + grad B^T) + lambda* I div B ) = rhs
51 //! where beta, lambda, and eta are incorporated into the flux F. It uses
52 //! the AMRLevelOp interface.
53 //class EBViscousTensorOp: public AMRLevelOp<LevelData<EBCellFAB> >
54 class EBViscousTensorOp: public LevelTGAHelmOp<LevelData<EBCellFAB>, EBFluxFAB >
55 {
56  public:
57 
58  static int s_whichLev;
59  static int s_step;
60 
61  //! Constructs a viscous tensor operator using the given data. This
62  //! constructor is for time-independent a and b coefficients.
63  EBViscousTensorOp(const EBLevelGrid & a_eblgFine,
64  const EBLevelGrid & a_eblg,
65  const EBLevelGrid & a_eblgCoar,
66  const EBLevelGrid & a_eblgCoarMG,
67  const Real& a_alpha,
68  const Real& a_beta,
69  const RefCountedPtr<LevelData<EBCellFAB> >& a_acoef,
70  const RefCountedPtr<LevelData<EBFluxFAB> >& a_eta,
71  const RefCountedPtr<LevelData<EBFluxFAB> >& a_lambda,
72  const RefCountedPtr<LevelData<BaseIVFAB<Real> > >& a_etaIrreg,
73  const RefCountedPtr<LevelData<BaseIVFAB<Real> > >& a_lambdaIrreg,
74  const Real& a_dx,
75  const Real& a_dxCoar,
76  const int& a_refToFine,
77  const int& a_refToCoar,
78  const RefCountedPtr<ViscousBaseDomainBC>& a_domainBC,
79  const RefCountedPtr<ViscousBaseEBBC>& a_ebBC,
80  const bool& a_hasMGObjects,
81  const IntVect& a_ghostCellsPhi,
82  const IntVect& a_ghostCellsRHS);
83 
84  //! Constructs a viscous tensor operator using the given data. This
85  //! constructor is for a time-dependent a coefficient.
86  EBViscousTensorOp(const EBLevelGrid & a_eblgFine,
87  const EBLevelGrid & a_eblg,
88  const EBLevelGrid & a_eblgCoar,
89  const EBLevelGrid & a_eblgCoarMG,
90  const Real& a_alpha,
91  const Real& a_beta,
92  const RefCountedPtr<LevelData<EBCellFAB> >& a_acoef0,
93  const RefCountedPtr<LevelData<EBCellFAB> >& a_acoef1,
94  const RefCountedPtr<LevelData<EBCellFAB> >& a_acoef,
95  const RefCountedPtr<LevelData<EBFluxFAB> >& a_eta,
96  const RefCountedPtr<LevelData<EBFluxFAB> >& a_lambda,
97  const RefCountedPtr<LevelData<BaseIVFAB<Real> > >& a_etaIrreg,
98  const RefCountedPtr<LevelData<BaseIVFAB<Real> > >& a_lambdaIrreg,
99  const Real& a_dx,
100  const Real& a_dxCoar,
101  const int& a_refToFine,
102  const int& a_refToCoar,
103  const RefCountedPtr<ViscousBaseDomainBC>& a_domainBC,
104  const RefCountedPtr<ViscousBaseEBBC>& a_ebBC,
105  const bool& a_hasMGObjects,
106  const IntVect& a_ghostCellsPhi,
107  const IntVect& a_ghostCellsRHS);
108 
109  //! Destructor.
110  virtual ~EBViscousTensorOp();
111 
112  ///for tga to reset stuff
113  virtual void setAlphaAndBeta(const Real& a_alpha,
114  const Real& a_beta);
115 
116 #ifdef CH_USE_HDF5
117  virtual void outputLevel(LevelData<EBCellFAB>& a_rhs, string& a_name)
118  {
119  char filename[1024];
120  sprintf(filename, "%s.ebvto.step%d.lev.%d.hdf5",a_name.c_str(), s_step, s_whichLev );
121  writeEBLevelname(&a_rhs, filename);
122  }
123 
124 
125  virtual void outputAMR(Vector< LevelData<EBCellFAB>* >& a_rhs, string& a_name)
126  {
127  char filename[1024];
128  sprintf(filename, "%s.ebvto.step%d.lev.%d.hdf5", a_name.c_str(), s_step, s_whichLev );
129  writeEBAMRname(&a_rhs, filename);
130  }
131 #endif
132 
133  ///it's tga's world---we just live in it.
134  virtual void kappaScale(LevelData<EBCellFAB> & a_rhs)
135  {
137  }
138 
139  /// compute (tau dot grad u) (for energy equation)
141  const LevelData<EBCellFAB>& a_gradU,
142  int a_level);
143 
144  /// compute volfrac(sigma dot grad u) (for energy equation)
145  void getKappaDivSigmaU(LevelData<EBCellFAB> & a_divSigmaU,
146  const LevelData<EBCellFAB>& a_velocity,
147  const LevelData<EBCellFAB>* a_veloCoar,
148  int a_level);
149 
150  //average surrounding coeffcients onto faces (simple averages)
152  LevelData<EBCellFAB>& a_lambdaCell);
153 
154  //get stress at cell centers given gradient of velocity
155  void getCCSigma(LevelData<EBCellFAB> & a_source,
156  const LevelData<EBCellFAB>& a_gradU,
157  const LevelData<EBCellFAB>& a_eta,
158  const LevelData<EBCellFAB>& a_lambda);
159  ///another tgaism
160  virtual void diagonalScale(LevelData<EBCellFAB> & a_rhs,
161  bool a_kappaWeighted)
162  {
163  if (a_kappaWeighted)
165  for (DataIterator dit = m_eblg.getDBL().dataIterator(); dit.ok(); ++dit)
166  {
167  for (int idir = 0; idir < SpaceDim; idir++)
168  {
169  int isrc = 0; int idst = idir; int inco = 1;
170  a_rhs[dit()].mult((*m_acoef)[dit()], isrc, idst, inco);
171  }
172  }
173  }
174 
176  {
177  for (DataIterator dit = m_eblg.getDBL().dataIterator(); dit.ok(); ++dit)
178  {
179  for (int idir = 0; idir < SpaceDim; idir++)
180  {
181  int isrc = 0; int idst = idir; int inco = 1;
182  a_rhs[dit()].divide((*m_acoef)[dit()], isrc, idst, inco);
183  }
184  }
185  }
186 
187  ///a leveltgaism
188  virtual void fillGrad(const LevelData<EBCellFAB>& a_phi);
189 
190  ///another leveltgaism
191  virtual void getFlux(EBFluxFAB& a_flux,
192  const LevelData<EBCellFAB>& a_data,
193  const Box& a_grid,
194  const DataIndex& a_dit,
195  Real a_scale);
196 
197  void getVelDotSigma(LevelData<EBFluxFAB> & a_velDotSigma,
198  const LevelData<EBFluxFAB>& a_vel,
199  const LevelData<EBFluxFAB>& a_sigma);
200 
201  void getFlux(EBFaceFAB& a_fluxCentroid,
202  const EBCellFAB& a_phi,
203  const Box& a_ghostedBox,
204  const Box& a_fabBox,
205  const ProblemDomain& a_domain,
206  const EBISBox& a_ebisBox,
207  const Real& a_dx,
208  const DataIndex& a_datInd,
209  const int& a_idir);
210 
211  // This operator may be time-dependent.
212  void setTime(Real a_oldTime, Real a_mu, Real a_dt);
213 
214  /// apply operator without any boundary or coarse-fine boundary conditions and no finer level
216  const LevelData<EBCellFAB>& a_phi)
217  {
218  s_turnOffBCs = true;
219  applyOp(a_opPhi, a_phi, true);
220  s_turnOffBCs = false;
221  }
222 
223  //! This is called on multigrid operators when their AMR operators
224  //! are altered.
225  void finerOperatorChanged(const MGLevelOp<LevelData<EBCellFAB> >& a_operator,
226  int a_coarseningFactor);
227 
228  //functions used by the wider world
229 
230  ///
231  /** a_residual = a_rhs - L(a_phiFine, a_phi) no coaser AMR level*/
232  void AMRResidualNC(LevelData<EBCellFAB>& a_residual,
233  const LevelData<EBCellFAB>& a_phiFine,
234  const LevelData<EBCellFAB>& a_phi,
235  const LevelData<EBCellFAB>& a_rhs,
236  bool a_homogeneousBC,
237  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
238 
239 
240  ///
241  /** apply AMR operator no coaser AMR level*/
242  void AMROperatorNC(LevelData<EBCellFAB>& a_LofPhi,
243  const LevelData<EBCellFAB>& a_phiFine,
244  const LevelData<EBCellFAB>& a_phi,
245  bool a_homogeneousBC,
246  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
247 
248  ///
249  /**
250  */
251  virtual void residual(LevelData<EBCellFAB>& a_residual,
252  const LevelData<EBCellFAB>& a_phi,
253  const LevelData<EBCellFAB>& a_rhs,
254  bool a_homogeneousPhysBC=false);
255 
256  ///
257  /**
258  */
259  virtual void preCond(LevelData<EBCellFAB>& a_opPhi,
260  const LevelData<EBCellFAB>& a_phi);
261 
262 
263  ///
264  /**
265  this is the linearop function. CFBC is set to homogeneous. phic is null
266  */
267  virtual void applyOp(LevelData<EBCellFAB>& a_opPhi,
268  const LevelData<EBCellFAB>& a_phi,
269  bool a_homogeneousPhysBC);
270 
271  ///
272  /**
273  */
274  virtual void create(LevelData<EBCellFAB>& a_lhs,
275  const LevelData<EBCellFAB>& a_rhs);
276 
277  ///
278  virtual void createCoarsened(LevelData<EBCellFAB>& a_lhs,
279  const LevelData<EBCellFAB>& a_rhs,
280  const int& a_refRat);
281 
282  ///
283  Real
284  AMRNorm(const LevelData<EBCellFAB>& a_coarResid,
285  const LevelData<EBCellFAB>& a_fineResid,
286  const int& a_refRat,
287  const int& a_ord);
288 
289 
290  ///
291  /**
292  */
293  virtual void assign(LevelData<EBCellFAB>& a_lhs,
294  const LevelData<EBCellFAB>& a_rhs);
295 
296  ///
297  /**
298  */
299  virtual Real dotProduct(const LevelData<EBCellFAB>& a_1,
300  const LevelData<EBCellFAB>& a_2);
301 
302  ///
303  /**
304  */
305  virtual void incr(LevelData<EBCellFAB>& a_lhs,
306  const LevelData<EBCellFAB>& a_x,
307  Real a_scale);
308 
309  ///
310  /**
311  */
312  virtual void axby(LevelData<EBCellFAB>& a_lhs,
313  const LevelData<EBCellFAB>& a_x,
314  const LevelData<EBCellFAB>& a_y,
315  Real a_a,
316  Real a_b);
317 
318  ///
319  /**
320  */
321  virtual void scale(LevelData<EBCellFAB>& a_lhs,
322  const Real& a_scale);
323 
324  ///
325  /**
326  */
327  virtual Real norm(const LevelData<EBCellFAB>& a_rhs,
328  int a_ord);
329 
330  ///
331  /**
332  */
333  virtual void setToZero(LevelData<EBCellFAB>& a_lhs);
334 
335  ///
336  /**
337  */
338  virtual void setVal(LevelData<EBCellFAB>& a_lhs, const Real& a_value);
339 
340  ///
341  /**
342  */
343  virtual void createCoarser(LevelData<EBCellFAB>& a_coarse,
344  const LevelData<EBCellFAB>& a_fine,
345  bool a_ghosted);
346 
347  ///
348  /**
349  */
350  virtual void relax(LevelData<EBCellFAB>& a_e,
351  const LevelData<EBCellFAB>& a_residual,
352  int a_iterations);
353 
354  ///
355  /**
356  Calculate restricted residual:
357  a_resCoarse[2h] = I[h->2h] (a_rhsFine[h] - L[h](a_phiFine[h]))
358  */
359  virtual void restrictResidual(LevelData<EBCellFAB>& a_resCoarse,
360  LevelData<EBCellFAB>& a_phiFine,
361  const LevelData<EBCellFAB>& a_rhsFine);
362 
363  ///
364  /**
365  Correct the fine solution based on coarse correction:
366  a_phiThisLevel += I[2h->h] (a_correctCoarse)
367  */
368  virtual void prolongIncrement(LevelData<EBCellFAB>& a_phiThisLevel,
369  const LevelData<EBCellFAB>& a_correctCoarse);
370 
371  ///
372  /** Refinement ratio between this level and coarser level.
373  Returns 1 when there are no coarser AMRLevelOp objects */
374  virtual int refToCoarser();
375 
376  ///
377  /** Refinement ratio between this level and coarser level.
378  Returns 1 when there are no coarser AMRLevelOp objects */
379  virtual int refToFiner();
380 
381  ///
382  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
383  virtual void AMRResidual(LevelData<EBCellFAB>& a_residual,
384  const LevelData<EBCellFAB>& a_phiFine,
385  const LevelData<EBCellFAB>& a_phi,
386  const LevelData<EBCellFAB>& a_phiCoarse,
387  const LevelData<EBCellFAB>& a_rhs,
388  bool a_homogeneousBC,
389  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
390 
391  ///
392  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
393  virtual void AMRResidualNF(LevelData<EBCellFAB>& a_residual,
394  const LevelData<EBCellFAB>& a_phi,
395  const LevelData<EBCellFAB>& a_phiCoarse,
396  const LevelData<EBCellFAB>& a_rhs,
397  bool a_homogeneousBC);
398 
399 
400  ///
401  /** a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse) */
402  virtual void AMROperator(LevelData<EBCellFAB>& a_LofPhi,
403  const LevelData<EBCellFAB>& a_phiFine,
404  const LevelData<EBCellFAB>& a_phi,
405  const LevelData<EBCellFAB>& a_phiCoarse,
406  bool a_homogeneousBC,
407  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
408 
409  ///
410  /** a_residual = a_rhs - L(a_phi, a_phiCoarse) */
411  virtual void AMROperatorNF(LevelData<EBCellFAB>& a_LofPhi,
412  const LevelData<EBCellFAB>& a_phi,
413  const LevelData<EBCellFAB>& a_phiCoarse,
414  bool a_homogeneousBC);
415 
416  ///
417  /** a_resCoarse = I[h-2h] (a_residual - L(a_correction, a_coarseCorrection)) */
418  virtual void AMRRestrict(LevelData<EBCellFAB>& a_resCoarse,
419  const LevelData<EBCellFAB>& a_residual,
420  const LevelData<EBCellFAB>& a_correction,
421  const LevelData<EBCellFAB>& a_coarseCorrection);
422 
423  ///
424  /** a_correction += I[2h->h](a_coarseCorrection) */
425  virtual void AMRProlong(LevelData<EBCellFAB>& a_correction,
426  const LevelData<EBCellFAB>& a_coarseCorrection);
427 
428  ///
429  /** a_residual = a_residual - L(a_correction, a_coarseCorrection) */
430  virtual void AMRUpdateResidual(LevelData<EBCellFAB>& a_residual,
431  const LevelData<EBCellFAB>& a_correction,
432  const LevelData<EBCellFAB>& a_coarseCorrection);
433 
434  ///
435  /** a_residual = a_residual - L(a_correction, a_coarseCorrection) */
436  static void
437  getDivergenceStencil(VoFStencil& a_divStencil,
438  const FaceIndex& a_face,
439  const DataIndex& a_dit,
440  const Real & a_dx,
441  const EBLevelGrid& a_eblg);
442 
443  ///
444  /** a_residual = a_residual - L(a_correction, a_coarseCorrection) */
445  static void
446  getGradientStencil(VoFStencil& a_gradStencil,
447  int a_ivar,
448  int a_diffDir,
449  const FaceIndex& a_face,
450  const DataIndex& a_dit,
451  const Real & a_dx,
452  const EBLevelGrid& a_eblg);
453 
454  ///
455  static void doLazyRelax(bool a_doLazyRelax)
456  {
457  s_doLazyRelax = a_doLazyRelax;
458  }
459 
460  //! (Re)define the stencils for the given coefficients.
461  void defineStencils();
462 
463  protected:
464 
465  //internal functions not usable without significant context
466  // void fillGrad(const LevelData<EBCellFAB>& a_phi);
467 
468  void incrOpRegularDir(EBCellFAB& a_lhs,
469  const EBCellFAB& a_phi,
470  const bool& a_homogeneous,
471  const int& a_dir,
472  const DataIndex& a_datInd);
473 
474  void applyOpIrregular(EBCellFAB& a_lhs,
475  const EBCellFAB& a_phi,
476  const bool& a_homogeneous,
477  const DataIndex& a_datInd);
478 
479  void getFlux(FArrayBox& a_flux,
480  const FArrayBox& a_phi,
481  const FArrayBox& a_gradPhi,
482  const Box& a_faceBox,
483  const int& a_idir,
484  const DataIndex& a_datInd);
485  void
486  cfinterp(const LevelData<EBCellFAB>& a_phi,
487  const LevelData<EBCellFAB>& a_phiCoarse);
488 
489  void reflux(const LevelData<EBCellFAB>& a_phiFine,
490  const LevelData<EBCellFAB>& a_phi,
491  LevelData<EBCellFAB>& residual,
492  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
493 
495 
496  void getVoFStencil(VoFStencil& a_vofStencil,
497  const VolIndex& a_vof,
498  const DataIndex& a_dit,
499  int a_ivar);
500 
501  void getFluxStencil(VoFStencil& a_fluxStencil,
502  const FaceIndex& a_face,
503  const DataIndex& a_dit,
504  int a_ivar);
505 
506  void getDivergenceStencil(VoFStencil& a_gradStencil,
507  const FaceIndex& a_face,
508  const DataIndex& a_dit);
509 
510  void getGradientStencil(VoFStencil& a_gradStencil,
511  int a_ivar,
512  int a_diffDir,
513  const FaceIndex& a_face,
514  const DataIndex& a_dit);
515 
516 
517  void gsrbColor(LevelData<EBCellFAB>& a_phi,
518  const LevelData<EBCellFAB>& a_lph,
519  const LevelData<EBCellFAB>& a_rhs,
520  const IntVect& a_color);
521 
522  void getFaceCenteredFluxStencil(VoFStencil& a_fluxStencil,
523  const FaceIndex& a_face,
524  const DataIndex& a_dit,
525  int a_ivar);
526 
527  void cellGrad(EBCellFAB& a_gradPhi,
528  const EBCellFAB& a_phi,
529  const Box& a_grid,
530  const DataIndex& a_datInd);
531 
532  void
534  const LevelData<EBCellFAB>& a_phiFine,
535  const LevelData<EBCellFAB>& a_phi);
536 
537  //simple averaging
538  void
540  const LevelData<EBFluxFAB>& a_faceCoef,
541  const LevelData<BaseIVFAB<Real> >& a_irregCoef);
542 
543  void
545  const LevelData<EBCellFAB>& a_phiFine,
546  const LevelData<EBCellFAB>& a_phi,
547  AMRLevelOp<LevelData<EBCellFAB> >* a_finerOp);
548 
549  //member data, bloated though it may be
550  static bool s_turnOffBCs;
551  static bool s_doLazyRelax;
552  //locations of coarse fine interfaces
556 
558  //grid information
563 
564  //operator coefficients
567  //weights that get multiplied by alpha
569  //weights that get multiplied by beta
571 
572  //! "Current" (time-interpolated) value of the a coefficient. For a
573  //! time-independent a coefficient, this is where the coefficient lives.
575 
576  //! The a coefficient at the beginning of the time step. This is NULL if
577  //! the conductivity operator is time independent.
579 
580  //! The a coefficient at the end of the time step. This is NULL if
581  //! the conductivity operator is time independent.
583 
588 
589 
590  //flux register with finer level
592 
593  //grid spacing
594  //dxcoar not implied by ref ratio--we could be in some MG level where dx > dxcoar
595  //needed for homogeneous cf interp
598  Real getSafety();
599  void calculateAlphaWeight();
601 
602  //refinement ratios and whether this object has multigrid objects
603  bool m_hasFine;
604  bool m_hasCoar;
608 
609  //ghost cell information needed for ebstencil optimizations
612 
613  //stencil to apply irregular operator
615 
616  //relaxation coefficent
618 
619  //gradient of solution at cell centers
621 
624  //for domain boundary conditions at ir regular cells
627 
628  //for inhomogeneous cf interpolation
630 
631  //averaging operators
634 
635  //interpolation operators
638 
639  //boundary condition operators
643 private:
644  ///weak construction is bad
646  {
647  MayDay::Error("invalid operator");
648  }
649 
650  //copy constructor and operator= disallowed for all the usual reasons
652  {
653  MayDay::Error("invalid operator");
654  }
655 
656  void operator=(const EBViscousTensorOp& a_opin)
657  {
658  MayDay::Error("invalid operator");
659  }
660 };
661 
662 
663 #include "NamespaceFooter.H"
664 #endif
static void doLazyRelax(bool a_doLazyRelax)
Definition: EBViscousTensorOp.H:455
RefCountedPtr< LevelData< EBFluxFAB > > m_lambda
Definition: EBViscousTensorOp.H:585
virtual void setToZero(LevelData< EBCellFAB > &a_lhs)
RefCountedPtr< LevelData< EBCellFAB > > m_acoef0
Definition: EBViscousTensorOp.H:578
IntVect m_ghostCellsRHS
Definition: EBViscousTensorOp.H:611
EBFastFR m_fastFR
Definition: EBViscousTensorOp.H:591
EBMGAverage m_ebAverage
Definition: EBViscousTensorOp.H:632
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
Definition: EBViscousTensorOp.H:568
virtual void scale(LevelData< EBCellFAB > &a_lhs, const Real &a_scale)
virtual void kappaScale(LevelData< EBCellFAB > &a_rhs)
it&#39;s tga&#39;s world—we just live in it.
Definition: EBViscousTensorOp.H:134
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)
virtual void createCoarsened(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat)
virtual void restrictResidual(LevelData< EBCellFAB > &a_resCoarse, LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_rhsFine)
A reference-counting handle class.
Definition: RefCountedPtr.H:66
EBFastFR-A class to encapsulate a levels worth of flux registers.
Definition: EBFastFR.H:39
void cfinterp(const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse)
#define CH_SPACEDIM
Definition: SPACE.H:52
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_etaIrreg
Definition: EBViscousTensorOp.H:586
LevelData< EBCellFAB > m_grad
Definition: EBViscousTensorOp.H:620
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:130
EBMGAverage m_ebAverageMG
Definition: EBViscousTensorOp.H:633
virtual void AMROperatorNF(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC)
static int s_whichLev
Definition: EBViscousTensorOp.H:58
EBViscousTensorOp(const EBViscousTensorOp &a_opin)
Definition: EBViscousTensorOp.H:651
void reflux(const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, LevelData< EBCellFAB > &residual, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
Real m_dxCoar
Definition: EBViscousTensorOp.H:597
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
Definition: EBViscousTensorOp.H:570
virtual void AMRProlong(LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
static void kappaWeight(LevelData< EBCellFAB > &a_data)
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)
static void getDivergenceStencil(VoFStencil &a_divStencil, const FaceIndex &a_face, const DataIndex &a_dit, const Real &a_dx, const EBLevelGrid &a_eblg)
one dimensional dynamic array
Definition: Vector.H:52
RefCountedPtr< LevelData< EBCellFAB > > m_acoef1
Definition: EBViscousTensorOp.H:582
void incrementFRCoar(EBFastFR &a_fr, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi)
EBMGInterp m_ebInterpMG
Definition: EBViscousTensorOp.H:637
Definition: FaceIndex.H:28
EBMGInterp m_ebInterp
Definition: EBViscousTensorOp.H:636
virtual void AMRUpdateResidual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
virtual int refToFiner()
RefCountedPtr< LevelData< EBCellFAB > > m_acoef
Definition: EBViscousTensorOp.H:574
static int s_step
Definition: EBViscousTensorOp.H:59
void averageToCells(LevelData< EBCellFAB > &a_cellCoef, const LevelData< EBFluxFAB > &a_faceCoef, const LevelData< BaseIVFAB< Real > > &a_irregCoef)
Definition: EBISBox.H:46
void getCCSigma(LevelData< EBCellFAB > &a_source, const LevelData< EBCellFAB > &a_gradU, const LevelData< EBCellFAB > &a_eta, const LevelData< EBCellFAB > &a_lambda)
void writeEBAMRname(const Vector< LevelData< EBCellFAB > * > *a_dataPtr, const char *a_filename)
void setTime(Real a_oldTime, Real a_mu, Real a_dt)
Definition: EBLevelGrid.H:30
LayoutData< RefCountedPtr< DivergenceStencil > > m_divergenceStencil
Definition: EBViscousTensorOp.H:555
DisjointBoxLayout getDBL() const
Definition: EBLevelGrid.H:86
virtual void assign(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
virtual bool ok() const
return true if this iterator is still in its Layout
Definition: LayoutIterator.H:110
void getFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit, int a_ivar)
Definition: DataIterator.H:140
EBViscousTensorOp()
weak construction is bad
Definition: EBViscousTensorOp.H:645
void writeEBLevelname(const LevelData< EBCellFAB > *a_dataPtr, const char *a_filename)
void getCellCenteredCoefficients(LevelData< EBCellFAB > &a_etaCell, LevelData< EBCellFAB > &a_lambdaCell)
Definition: EBFaceFAB.H:28
void getShearStressDotGradU(LevelData< EBCellFAB > &a_source, const LevelData< EBCellFAB > &a_gradU, int a_level)
compute (tau dot grad u) (for energy equation)
virtual Real norm(const LevelData< EBCellFAB > &a_rhs, int a_ord)
virtual Real dotProduct(const LevelData< EBCellFAB > &a_1, const LevelData< EBCellFAB > &a_2)
RefCountedPtr< EBTensorCFInterp > m_interpWithCoarser
Definition: EBViscousTensorOp.H:629
LayoutData< CFIVS > m_hiCFIVS[CH_SPACEDIM]
Definition: EBViscousTensorOp.H:554
void getVelDotSigma(LevelData< EBFluxFAB > &a_velDotSigma, const LevelData< EBFluxFAB > &a_vel, const LevelData< EBFluxFAB > &a_sigma)
virtual void divideByIdentityCoef(LevelData< EBCellFAB > &a_rhs)
Definition: EBViscousTensorOp.H:175
Piecewise constant interpolation.
Definition: EBMGInterp.H:33
int m_refToCoar
Definition: EBViscousTensorOp.H:606
const int SpaceDim
Definition: SPACE.H:39
Definition: AMRMultiGrid.H:35
VoF-centered stencil.
Definition: Stencils.H:59
static void getGradientStencil(VoFStencil &a_gradStencil, int a_ivar, int a_diffDir, const FaceIndex &a_face, const DataIndex &a_dit, const Real &a_dx, const EBLevelGrid &a_eblg)
static bool s_turnOffBCs
Definition: EBViscousTensorOp.H:550
Real m_dx
Definition: EBViscousTensorOp.H:596
virtual void getFlux(EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale)
another leveltgaism
LayoutData< VoFIterator > m_vofIterMulti
Definition: EBViscousTensorOp.H:623
A EBFaceFAB-like container for edge-centered fluxes.
Definition: EBFluxFAB.H:25
virtual void outputLevel(LevelData< EBCellFAB > &a_rhs, string &a_name)
Definition: EBViscousTensorOp.H:117
void finerOperatorChanged(const MGLevelOp< LevelData< EBCellFAB > > &a_operator, int a_coarseningFactor)
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_lambdaIrreg
Definition: EBViscousTensorOp.H:587
IntVect m_ghostCellsPhi
Definition: EBViscousTensorOp.H:610
virtual void axby(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, Real a_a, Real a_b)
void calculateRelaxationCoefficient()
LayoutData< CFIVS > m_loCFIVS[CH_SPACEDIM]
Definition: EBViscousTensorOp.H:553
virtual void applyOpNoBoundary(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
apply operator without any boundary or coarse-fine boundary conditions and no finer level ...
Definition: EBViscousTensorOp.H:215
virtual void setAlphaAndBeta(const Real &a_alpha, const Real &a_beta)
for tga to reset stuff
Definition: EBCellFAB.H:29
virtual void create(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
EBLevelGrid m_eblgFine
Definition: EBViscousTensorOp.H:559
double Real
Definition: REAL.H:33
Definition: MultiGrid.H:30
IntVect m_ghostPhi
Definition: EBViscousTensorOp.H:557
EBLevelGrid m_eblg
Definition: EBViscousTensorOp.H:560
virtual void incr(LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, Real a_scale)
virtual void prolongIncrement(LevelData< EBCellFAB > &a_phiThisLevel, const LevelData< EBCellFAB > &a_correctCoarse)
void getVoFStencil(VoFStencil &a_vofStencil, const VolIndex &a_vof, const DataIndex &a_dit, int a_ivar)
void gsrbColor(LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_lph, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color)
Real m_beta
Definition: EBViscousTensorOp.H:566
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)
RefCountedPtr< ViscousBaseDomainBC > m_domainBC
Definition: EBViscousTensorOp.H:640
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)
void getFaceCenteredFluxStencil(VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit, int a_ivar)
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.
Real AMRNorm(const LevelData< EBCellFAB > &a_coarResid, const LevelData< EBCellFAB > &a_fineResid, const int &a_refRat, const int &a_ord)
Vector< IntVect > m_colors
Definition: EBViscousTensorOp.H:642
void applyOpIrregular(EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const DataIndex &a_datInd)
virtual void diagonalScale(LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted)
another tgaism
Definition: EBViscousTensorOp.H:160
void defineStencils()
(Re)define the stencils for the given coefficients.
void incrOpRegularDir(EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const int &a_dir, const DataIndex &a_datInd)
A Rectangular Domain on an Integer Lattice.
Definition: Box.H:465
void homogeneousCFInterp(LevelData< EBCellFAB > &a_phi)
void getKappaDivSigmaU(LevelData< EBCellFAB > &a_divSigmaU, const LevelData< EBCellFAB > &a_velocity, const LevelData< EBCellFAB > *a_veloCoar, int a_level)
compute volfrac(sigma dot grad u) (for energy equation)
Definition: EBViscousTensorOp.H:54
Definition: DataIndex.H:112
IntVect m_ghostRHS
Definition: EBViscousTensorOp.H:557
LayoutData< VoFIterator > m_vofIterIrreg
Definition: EBViscousTensorOp.H:622
virtual ~EBViscousTensorOp()
Destructor.
bool m_hasCoar
Definition: EBViscousTensorOp.H:604
Piecewise constant interpolation.
Definition: EBMGAverage.H:31
EBLevelGrid m_eblgCoar
Definition: EBViscousTensorOp.H:561
int m_refToFine
Definition: EBViscousTensorOp.H:605
virtual void relax(LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
void cellGrad(EBCellFAB &a_gradPhi, const EBCellFAB &a_phi, const Box &a_grid, const DataIndex &a_datInd)
An integer Vector in SpaceDim-dimensional space.
Definition: CHArray.H:42
virtual void residual(LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC=false)
Definition: FArrayBox.H:44
Definition: AMRTGA.H:159
void AMROperatorNC(LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
RefCountedPtr< ViscousBaseEBBC > m_ebBC
Definition: EBViscousTensorOp.H:641
Volume of Fluid Index.
Definition: VolIndex.H:31
EBLevelGrid m_eblgCoarMG
Definition: EBViscousTensorOp.H:562
virtual int refToCoarser()
virtual void preCond(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
virtual void fillGrad(const LevelData< EBCellFAB > &a_phi)
a leveltgaism
RefCountedPtr< LevelData< EBFluxFAB > > m_eta
Definition: EBViscousTensorOp.H:584
void operator=(const EBViscousTensorOp &a_opin)
Definition: EBViscousTensorOp.H:656
virtual void applyOp(LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousPhysBC)
virtual void outputAMR(Vector< LevelData< EBCellFAB > * > &a_rhs, string &a_name)
Definition: EBViscousTensorOp.H:125
void incrementFRFine(EBFastFR &a_fr, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
DataIterator dataIterator() const
Parallel iterator.
Real m_alpha
Definition: EBViscousTensorOp.H:565
void calculateAlphaWeight()
LayoutData< VoFIterator > m_vofIterDomHi[CH_SPACEDIM]
Definition: EBViscousTensorOp.H:626
virtual void createCoarser(LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted)
virtual void AMRRestrict(LevelData< EBCellFAB > &a_resCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
LevelData< EBCellFAB > m_relCoef
Definition: EBViscousTensorOp.H:617
bool m_hasFine
Definition: EBViscousTensorOp.H:603
static bool s_doLazyRelax
Definition: EBViscousTensorOp.H:551
LayoutData< VoFIterator > m_vofIterDomLo[CH_SPACEDIM]
Definition: EBViscousTensorOp.H:625
LayoutData< RefCountedPtr< EBStencil > > m_opEBStencil[CH_SPACEDIM]
Definition: EBViscousTensorOp.H:614
bool m_hasMGObjects
Definition: EBViscousTensorOp.H:607
virtual void setVal(LevelData< EBCellFAB > &a_lhs, const Real &a_value)