Chombo + EB + MF  3.2
MultilevelLinearOpI.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 _MULTILEVELLINEAROPI_H_
12 #define _MULTILEVELLINEAROPI_H_
13 
14 #include "AMRIO.H"
15 
16 #include "NamespaceHeader.H"
17 
18 // set a couple of defaults here (do before define so that
19 // they can be overridden before calling define)
20 template<class T>
22 {
23  // set some default values
24  m_use_multigrid_preconditioner = true;
25  m_num_mg_iterations = 1;
26  m_num_mg_smooth = 4;
27  // default is to coarsen all the way down
28  m_preCondSolverDepth = -1;
29  m_precondBottomSolverPtr = NULL;
30 }
31 
32 /// define function
33 template<class T>
34 void
36  const Vector<int>& a_refRatios,
37  const Vector<ProblemDomain>& a_domains,
38  const Vector<RealVect>& a_vectDx,
40  int a_lBase)
41 {
42  CH_TIME("MultilevelLinearOp::define");
43 
44  int numLevels = a_vectGrids.size();
45 
46  // couple of sanity checks...
47  CH_assert (a_lBase >= 0);
48  CH_assert (a_lBase < numLevels);
49  CH_assert ( a_domains.size() == numLevels);
50  // since you technically only need numLevels -1 refinement ratios...
51  CH_assert (a_refRatios.size() >= numLevels -1);
52 
53  m_lBase = a_lBase;
54 
55  m_vectGrids = a_vectGrids;
56  m_refRatios = a_refRatios;
57  m_domains = a_domains;
58  m_vectDx = a_vectDx;
59  m_vectOperators.resize(numLevels);
60 
61  // need to at least define the levels to lBase -1 (for coarser BC's)
62  int coarsestLevel = Max(0, m_lBase-1);
63  for (int level=coarsestLevel; level<numLevels; level++)
64  {
65  RefCountedPtr<AMRLevelOp<LevelData<T> > > op = RefCountedPtr<AMRLevelOp<LevelData<T> > >(a_opFactory->AMRnewOp(m_domains[level]));
66  m_vectOperators[level] = op;
67  }
68 
69  // finally, define AMRMultigrid preconditioner if required
70  if (m_use_multigrid_preconditioner)
71  {
72  // preconditioner requires a bottom smoother
74 
75  m_precondBottomSolverPtr = newPtr;
76  int bottomSolverVerbosity = 1;
77 
78  newPtr->m_verbosity = bottomSolverVerbosity;
79 
80  if (m_preCondSolverDepth >= 0)
81  {
82  m_preCondSolver.m_maxDepth = m_preCondSolverDepth;
83  }
84  m_preCondSolver.define(m_domains[0],
85  *a_opFactory,
86  m_precondBottomSolverPtr,
87  numLevels);
88 
89  // most of these have no real meaning since we're only doina a
90  // single V-cycle, rather than a full solve
91  Real solverEps = 1.0e-7;
92  int maxIterations = 1;
93  Real hang = 1.0e-7;
94  Real normThresh = 1.0e-30;
95 
96  int num_mg = 1;
97  m_preCondSolver.setSolverParameters(m_num_mg_smooth,
98  m_num_mg_smooth,
99  m_num_mg_smooth,
100  num_mg,
101  maxIterations,
102  solverEps,
103  hang,
104  normThresh);
105 
106  // set preconditioner solver to be _really_ quiet...
107  m_preCondSolver.m_verbosity = 1;
108 
109  // AMRMultiGrid::init has not yet been called
110  // (will call later, during actual solve)
111  m_isPrecondSolverInitialized = false;
112 
113  }
114 }
115 
116  ///
117 template<class T>
119 {
120  CH_TIME("MultilevelLinearOp::~MultilevelLinearOp");
121 
122  if (m_precondBottomSolverPtr != NULL)
123  {
124  delete m_precondBottomSolverPtr;
125  m_precondBottomSolverPtr = NULL;
126  }
127 }
128 
129 ///
130 /**
131  Say you are solving L(phi) = rhs. Make a_lhs = L(a_phi) - a_rhs.
132  If a_homogeneous is true, evaluate the operator using homogeneous
133  boundary conditions.
134 */
135 template<class T>
136 void
138  const Vector<LevelData<T>* >& a_phi,
139  const Vector<LevelData<T>* >& a_rhs,
140  bool a_homogeneous )
141 {
142  CH_TIME("MultilevelLinearOp::residual");
143 
144  // do this using the operators in AMRLevelOp
145  LevelData<T> dummyT;
146 
147  int numLevels = a_lhs.size();
148  CH_assert (a_phi.size() == numLevels);
149  CH_assert (a_rhs.size() == numLevels);
150  CH_assert (m_vectOperators.size() >= numLevels);
151 
152  int maxLevel = numLevels -1;
153  for (int level = m_lBase; level < numLevels; level++)
154  {
155  const LevelData<T>& thisPhi = *(a_phi[level]);
156  LevelData<T>& thisResid = *(a_lhs[level]);
157  const LevelData<T>& thisRHS = *(a_rhs[level]);
158  // theres gotta be a cleaner way to do this...
159  if (level >0)
160  {
161  const LevelData<T>& phiCrse = *a_phi[level-1];
162 
163  if (level < maxLevel)
164  {
165  const LevelData<T>& phiFine = *a_phi[level+1];
166  m_vectOperators[level]->AMRResidual(thisResid,
167  phiFine,
168  thisPhi,
169  phiCrse,
170  thisRHS,
171  a_homogeneous,
172  (m_vectOperators[level+1].operator->()));
173  }
174  else
175  {
176  // no finer level exists
177  m_vectOperators[level]->AMRResidualNF(thisResid,
178  thisPhi,
179  phiCrse,
180  thisRHS,
181  a_homogeneous);
182  } // end if no finer level
183 
184  } // end if a coarser level exists
185  else
186  {
187  // we're on level 0
188  if (level < maxLevel)
189  {
190  const LevelData<T>& phiFine = *a_phi[level+1];
191  m_vectOperators[level]->AMRResidualNC(thisResid,
192  phiFine,
193  thisPhi,
194  thisRHS,
195  a_homogeneous,
196  (m_vectOperators[level+1].operator->()));
197  }
198  else
199  {
200  // no finer level exists
201  m_vectOperators[level]->residual(thisResid,
202  thisPhi,
203  thisRHS,
204  a_homogeneous);
205  } // end if no finer level
206  } // end if on level 0
207  } // end loop over levels
208 }
209 
210 ///
211 /**
212  Given the current state of the residual the correction, apply
213  your preconditioner to a_cor.
214 */
215 template<class T>
216 void
218  const Vector<LevelData<T>* >& a_residual)
219 
220 {
221  CH_TIME("MultilevelLinearOp::preCond");
222 
223  // sanity checks
224  int numLevels = a_cor.size();
225  CH_assert (a_residual.size() == numLevels);
226  CH_assert (m_vectOperators.size() <= numLevels);
227 
228  // use AMR multigrid for a preconditioner
229  if (m_use_multigrid_preconditioner)
230  {
231  CH_TIME("MultilevelLinearOp::preCond::Multigrid");
232 
233  Vector<LevelData<T>* > localCorr;
234  Vector<LevelData<T>* > localResid;
235  create(localCorr, a_cor);
236  create(localResid, a_residual);
237 
238  // this is kinda silly, but I need to get rid of RefCountedPtr
239  // part of things
240  Vector<LevelData<T>* > tempCorr(numLevels, NULL);
241  Vector<LevelData<T>* > tempResid(numLevels, NULL);
242 
243  {
244  CH_TIME("MultilevelLinearOp::preCond::Multigrid::NewStorage");
245 
246  // in place. Otherwise, need to allocate new storage for them
247  for (int lev=0; lev<numLevels; lev++)
248  {
249  tempCorr[lev] = localCorr[lev];
250  tempResid[lev] = const_cast<LevelData<T>*>(localResid[lev]);
251  }
252  }
253 
254  int finestLevel = numLevels -1;
255  {
256  CH_TIME("MultilevelLinearOp::preCond::Multigrid::Initialize");
257 
258  if (!m_isPrecondSolverInitialized)
259  {
260  // because we're not calling solve, need to call init directly
261  m_preCondSolver.init(tempCorr, tempResid, finestLevel,
262  m_lBase);
263  // also need to set bottom solver
264  m_preCondSolver.setBottomSolver(finestLevel, m_lBase);
265  m_isPrecondSolverInitialized = true;
266  }
267  }
268 
269  for (int iter = 0; iter<m_num_mg_iterations; iter++)
270  {
271  CH_TIME("MultilevelLinearOp::preCond::Multigrid::Iterate");
272 
273  // use AMRVCycle instead of solve because we're already
274  // in residual-correction form, so we want to use
275  // homogeneous form of boundary conditions.
276 
277  // however, AMRVCycle requires that the initial correction be zero
278  // to do this properly, recompute residual before calling V-Cycle
279 
280  bool homogeneous = true;
281  residual(localResid, a_cor, a_residual, homogeneous);
282  setToZero(localCorr);
283 
284  m_preCondSolver.AMRVCycle(tempCorr, tempResid,
285  finestLevel, finestLevel,
286  m_lBase);
287 
288  // now increment a_cor with local correction
289  incr(a_cor, localCorr, 1.0);
290  }
291 
292  {
293  CH_TIME("MultilevelLinearOp::preCond::Multigrid::Clear");
294 
295  clear(localCorr);
296  clear(localResid);
297  }
298  }
299  else
300  {
301  CH_TIME("MultilevelLinearOp::preCond::LevelBased");
302 
303  // just use whatever the AMRLevelOp provides for levels finer than lBase
304  for (int level = m_lBase; level < numLevels; level++)
305  {
306  m_vectOperators[level]->preCond(*a_cor[level], *a_residual[level]);
307  }
308  }
309 }
310 
311 ///
312 /**
313  In the context of solving L(phi) = rhs, set a_lhs = L(a_phi).
314  If a_homogeneous is true,
315  evaluate the operator using homogeneous boundary conditions.
316 */
317 template<class T>
318 void
320  const Vector<LevelData<T>* >& a_phi,
321  bool a_homogeneous)
322 {
323  CH_TIME("MultilevelLinearOp::applyOp");
324 
325  int numLevels = a_phi.size();
326  int maxLevel = numLevels-1;
327  LevelData<T> dummyLDF;
328  for (int level = m_lBase; level<numLevels; level++)
329  {
330  if (level == 0)
331  {
332  if (level < maxLevel)
333  {
334  m_vectOperators[level]->AMROperatorNC(*a_lhs[level],
335  *a_phi[level+1],
336  *a_phi[level],
337  a_homogeneous,
338  m_vectOperators[level+1].operator->());
339  }
340  else
341  {
342  // this is the single-level case...
343  m_vectOperators[level]->applyOp(*a_lhs[level],
344  *a_phi[level],
345  a_homogeneous);
346  }
347  } // end if level is 0
348  else
349  {
350  // coarser level exists
351  if (level < maxLevel)
352  {
353  m_vectOperators[level]->AMROperator(*a_lhs[level],
354  *a_phi[level+1],
355  *a_phi[level],
356  *a_phi[level-1],
357  a_homogeneous,
358  m_vectOperators[level+1].operator->());
359  }
360  else
361  {
362  m_vectOperators[level]->AMROperatorNF(*a_lhs[level],
363  *a_phi[level],
364  *a_phi[level-1],
365  a_homogeneous);
366  }
367  } // end if coarser level exists
368  } // end loop over levels
369 
370 }
371 
372 ///
373 /**
374  Create data holder a_lhs that mirrors a_rhs. You do not need
375  to copy the data of a_rhs, just make a holder the same size.
376 */
377 template<class T>
378 void
380  const Vector<LevelData<T>* >& a_rhs)
381 {
382  CH_TIME("MultilevelLinearOp::create");
383 
384  a_lhs.resize(a_rhs.size());
385  // need to create an lBase-1 just in case it's needed for C/F BCs
386  int coarsestLevel = Max(0, m_lBase-1);
387  for (int level =coarsestLevel; level < a_rhs.size(); level++)
388  {
389  a_lhs[level] = new LevelData<T>;
390  m_vectOperators[level]->create(*a_lhs[level],
391  *a_rhs[level]);
392  }
393 }
394 
395 ///
396 /**
397  Clean up data holder before it goes out of scope. This is necessary
398  because create calls new.
399 */
400 template<class T>
401 void
403 {
404  CH_TIME("MultilevelLinearOp::clear");
405 
406  for (int level =0; level < a_lhs.size(); level++)
407  {
408  if (a_lhs[level] != NULL)
409  {
410  delete a_lhs[level];
411  a_lhs[level] = NULL;
412  }
413  }
414 }
415 
416 ///
417 /**
418  Set a_lhs equal to a_rhs.
419 */
420 template<class T>
421 void
423  const Vector<LevelData<T>* >& a_rhs)
424 {
425  CH_TIME("MultilevelLinearOp::assign");
426 
427  CH_assert (a_lhs.size() == a_rhs.size());
428  // only do this for lBase and finer levels
429  for (int level = m_lBase; level < a_lhs.size(); level++)
430  {
431  if (!(a_lhs[level] ==NULL) && !(a_rhs[level] == NULL))
432  {
433  m_vectOperators[level]->assign(*a_lhs[level], *a_rhs[level]);
434  }
435  }
436 }
437 
438 ///
439 /**
440  Compute and return the dot product of a_1 and a_2. In most
441  contexts, this means return the sum over all data points of a_1*a_2.
442 */
443 template<class T>
444 Real
446  const Vector<LevelData<T>* >& a_2)
447 {
448  CH_TIME("MultilevelLinearOp::dotProduct");
449 
450  // want to do this in an AMR way, so need to generate a temporary
451  Vector<LevelData<T>* > temp1, temp2;
452  create(temp1, a_1);
453  create (temp2, a_2);
454 
455  // first set to zero, then call assign, since assign only sets
456  // valid regions (that way ghost cells are set to zero)
457 
458  setToZero(temp1);
459  setToZero(temp2);
460 
461  assign(temp1, a_1);
462  assign(temp2, a_2);
463 
464  // now set covered regions to zero
465  for (int level =m_lBase; level<temp1.size()-1; level++)
466  {
467  LevelData<T>& temp1Level = *temp1[level];
468  LevelData<T>& temp2Level = *temp2[level];
469 
470  CH_assert(temp1[level]->getBoxes() == temp2[level]->getBoxes());
471  CH_assert(temp1[level+1] != NULL);
472 
473  int nRefFine = m_refRatios[level];
474  const DisjointBoxLayout& finerGrids = temp1[level+1]->getBoxes();
475  const DisjointBoxLayout& levelGrids = temp1[level]->getBoxes();
476  DataIterator levelDit = levelGrids.dataIterator();
477  LayoutIterator finerLit = finerGrids.layoutIterator();
478 
479  for (levelDit.begin(); levelDit.ok(); ++levelDit)
480  {
481  const Box& thisBox = levelGrids[levelDit];
482  for (finerLit.begin(); finerLit.ok(); ++finerLit)
483  {
484  Box testBox = finerGrids[finerLit];
485  testBox.coarsen(nRefFine);
486  testBox &= thisBox;
487  if (!testBox.isEmpty())
488  {
489  temp1Level[levelDit].setVal(0.0, testBox,
490  0, temp1Level.nComp());
491 
492  temp2Level[levelDit].setVal(0.0, testBox,
493  0, temp2Level.nComp());
494  }
495  } // end loop over finer boxes
496  } // end loop over boxes on this level
497  } // end loop over levels for setting covered regions to zero;
498 
499  // now loop over levels and call AMRLevelOp dotProduct
500  Real prod = 0.0;
501  for (int level =m_lBase; level<temp1.size(); level++)
502  {
503  Real levelProd = m_vectOperators[level]->dotProduct(*temp1[level],
504  *temp2[level]);
505 
506  // incorporate scaling for AMR dot products
507  RealVect dxLev = m_vectDx[level];
508  Real scale = D_TERM(dxLev[0], *dxLev[1], *dxLev[2]);
509  levelProd *= scale;
510  prod += levelProd;
511  }
512 
513  clear(temp1);
514  clear (temp2);
515 
516  return prod;
517 }
518 
519 ///
520 /**
521  Increment by scaled amount (a_lhs += a_scale*a_x).
522 */
523 template<class T>
524 void
526  const Vector<LevelData<T>* >& a_x,
527  Real a_scale)
528 {
529  CH_TIME("MultilevelLinearOp::incr");
530 
531  // int this case, only operate on lBase and finer
532  for (int level=m_lBase; level<a_lhs.size(); level++)
533  {
534  m_vectOperators[level]->incr(*a_lhs[level], *a_x[level], a_scale);
535  }
536 }
537 
538 ///
539 /**
540  Set input to a scaled sum (a_lhs = a_a*a_x + a_b*a_y).
541 */
542 template<class T>
543 void
545  const Vector<LevelData<T>* >& a_x,
546  const Vector<LevelData<T>* >& a_y,
547  Real a_a,
548  Real a_b)
549 {
550  CH_TIME("MultilevelLinearOp::axby");
551 
552  // only do this for lBase and finer
553  for (int level=m_lBase; level<a_lhs.size(); ++level)
554  {
555  m_vectOperators[level]->axby(*a_lhs[level],
556  *a_x[level],
557  *a_y[level],
558  a_a, a_b);
559  }
560 }
561 
562 ///
563 /**
564  Multiply the input by a given scale (a_lhs *= a_scale).
565 */
566 template<class T>
567 void
569  const Real& a_scale)
570 {
571  CH_TIME("MultilevelLinearOp::scale");
572 
573  // only do this for lBase and finer
574  for (int level =m_lBase; level<a_lhs.size(); level++)
575  {
576  m_vectOperators[level]->scale(*a_lhs[level], a_scale);
577  }
578 
579 }
580 
581 ///
582 /**
583  Return the norm of a_rhs.
584  a_ord == 0 max norm, a_ord == 1 sum(abs(a_rhs)), else, L(a_ord) norm.
585 */
586 template<class T>
587 Real
589  int a_ord)
590 {
591  CH_TIME("MultilevelLinearOp::norm");
592 
593  // now loop over levels and call AMRLevelOp AMRNorm
594  int maxLevel = a_rhs.size()-1;
595  Real thisNorm = 0.0;
596  for (int level = m_lBase; level<a_rhs.size(); level++)
597  {
598  Real normLevel;
599 
600  if (level < maxLevel)
601  {
602  // finer level exists
603  int refRatio = m_refRatios[level];
604  normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
605  *a_rhs[level+1],
606  refRatio,
607  a_ord);
608  }
609  else
610  {
611  // no finer level exists
612  LevelData<T> tempLDF;
613  int refRatio = -1;
614  normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
615  tempLDF,
616  refRatio,
617  a_ord);
618  }
619 
620  // now need to do a bit of rescaling to make things work out right
621  if (a_ord > 2)
622  {
623  MayDay::Error("MultilevelLinearOp::norm -- undefined for order > 2");
624  }
625 
626  if (a_ord == 2)
627  {
628  normLevel = normLevel*normLevel;
629  }
630 
631  if (a_ord != 0)
632  {
633  thisNorm += normLevel;
634  }
635  else if (normLevel > thisNorm)
636  {
637  thisNorm = normLevel;
638  }
639  } // end
640 
641  if (a_ord == 2)
642  {
643  thisNorm = sqrt(thisNorm);
644  }
645 
646  return thisNorm;
647 }
648 
649 ///
650 /**
651  Set a_lhs to zero.
652 */
653 template<class T>
654 void
656 {
657  CH_TIME("MultilevelLinearOp::setToZero");
658 
659  // do this for lBase -1 if necessary, to account for cases where correction
660  // lBase-1 needs to be zero
661  int coarsestLevel = Max(0, m_lBase-1);
662  for (int level = coarsestLevel; level< a_lhs.size(); level++)
663  {
664  m_vectOperators[level]->setToZero(*a_lhs[level]);
665  }
666 }
667 
668 template<class T>
669 void
671  const char* a_filename)
672 {
673 #ifdef CH_USE_HDF5
674  writeVectorLevelName(a_data, &m_refRatios, a_filename);
675 #else
676  MayDay::Warning("MultilevelLinearOp<T>::write unimplemented since CH_USE_HDF5 undefined");
677 #endif
678 }
679 
680 #include "NamespaceFooter.H"
681 
682 #endif
Box & coarsen(int refinement_ratio)
virtual void incr(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_x, Real a_scale)
Definition: MultilevelLinearOpI.H:525
virtual Real norm(const Vector< LevelData< T > * > &a_rhs, int a_ord)
Definition: MultilevelLinearOpI.H:588
A reference-counting handle class.
Definition: RefCountedPtr.H:173
virtual void write(const Vector< LevelData< T > * > *a_data, const char *a_filename)
Definition: MultilevelLinearOpI.H:670
#define CH_assert(cond)
Definition: CHArray.H:37
int m_verbosity
Definition: BiCGStabSolver.H:78
void begin()
initialize this iterator to the first Box in its Layout
Definition: LayoutIterator.H:122
virtual bool ok() const
return true if this iterator is still in its Layout
Definition: LayoutIterator.H:117
IndexTM< T, N > scale(const IndexTM< T, N > &a_p, T a_s)
Definition: IndexTMI.H:382
Definition: DataIterator.H:190
An Iterator based on a BoxLayout object.
Definition: LayoutIterator.H:35
virtual void axby(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_x, const Vector< LevelData< T > * > &a_y, Real a_a, Real a_b)
Definition: MultilevelLinearOpI.H:544
void resize(unsigned int isize)
Definition: Vector.H:346
void getBoxes(Vector< Vector< Box > > &a_boxes, Vector< int > &a_refRat, const Box &a_domain)
virtual ~MultilevelLinearOp()
Definition: MultilevelLinearOpI.H:118
virtual void setToZero(Vector< LevelData< T > * > &a_lhs)
Definition: MultilevelLinearOpI.H:655
#define CH_TIME(name)
Definition: CH_Timer.H:82
MultilevelLinearOp()
default constructor (sets defaults for preconditioner)
Definition: MultilevelLinearOpI.H:21
virtual void scale(Vector< LevelData< T > * > &a_lhs, const Real &a_scale)
Definition: MultilevelLinearOpI.H:568
new code
Definition: BoxLayoutData.H:170
double Real
Definition: REAL.H:33
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:30
virtual void clear(Vector< LevelData< T > * > &a_lhs)
Definition: MultilevelLinearOpI.H:402
size_t size() const
Definition: Vector.H:192
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 assign(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs)
Definition: MultilevelLinearOpI.H:422
virtual void applyOp(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_phi, bool a_homogeneous=false)
Definition: MultilevelLinearOpI.H:319
static void Warning(const char *const a_msg=m_nullString)
Print out message to cerr and continue.
virtual void residual(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_phi, const Vector< LevelData< T > * > &a_rhs, bool a_homogeneous=false)
compute residual using AMRLevelOps AMRResidual functionality
Definition: MultilevelLinearOpI.H:137
A Rectangular Domain on an Integer Lattice.
Definition: Box.H:469
A Real vector in SpaceDim-dimensional space.
Definition: RealVect.H:41
virtual void preCond(Vector< LevelData< T > * > &a_cor, const Vector< LevelData< T > * > &a_residual)
Apply preconditioner.
Definition: MultilevelLinearOpI.H:217
int nComp() const
Definition: BoxLayoutData.H:306
T Max(const T &a_a, const T &a_b)
Definition: Misc.H:39
virtual void define(const Vector< DisjointBoxLayout > &a_vectGrids, const Vector< int > &a_refRatios, const Vector< ProblemDomain > &a_domains, const Vector< RealVect > &a_vectDx, RefCountedPtr< AMRLevelOpFactory< LevelData< T > > > &a_opFactory, int a_lBase)
define function – opFactory should be able to define all required levels
Definition: MultilevelLinearOpI.H:35
virtual Real dotProduct(const Vector< LevelData< T > * > &a_1, const Vector< LevelData< T > * > &a_2)
Definition: MultilevelLinearOpI.H:445
void writeVectorLevelName(const Vector< LevelData< FArrayBox > *> *a_dataPtr, const Vector< int > *a_refRatios, const char *a_filename)
LayoutIterator layoutIterator() const
Iterator that processes through ALL the boxes in a BoxLayout.
Definition: AMRMultiGrid.H:233
DataIterator dataIterator() const
Parallel iterator.
bool isEmpty() const
{ Comparison Functions}
Definition: Box.H:1862
virtual void create(Vector< LevelData< T > * > &a_lhs, const Vector< LevelData< T > * > &a_rhs)
Definition: MultilevelLinearOpI.H:379
Definition: BiCGStabSolver.H:26