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