Chombo + EB  3.2
EBLevelTGA.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 _EBLEVELTGA_H_
12 #define _EBLEVELTGA_H_
13 
14 #include <iostream>
15 #include <math.h>
16 #include "SPACE.H"
17 #include <stdlib.h>
18 #include "REAL.H"
19 #include "Box.H"
20 #include "DisjointBoxLayout.H"
21 #include "LevelData.H"
22 #include "EBCellFAB.H"
23 #include "EBLevelGrid.H"
24 #include "EBFluxFAB.H"
25 #include "EBFluxRegister.H"
26 #include "ProblemDomain.H"
27 #include "BaseLevelTGA.H"
28 #include "BaseLevelBackwardEuler.H"
29 #include "BaseLevelCrankNicolson.H"
30 #include "NamespaceHeader.H"
31 
32 
33 //! \class EBLevelTGA
34 //! This class implements an imlplicit integrator for a diffusion
35 //! equation using an approach devised by Twizzell, Gumel, and Arigu.
36 class EBLevelTGA : public BaseLevelTGA<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>
37 {
38 
39 public:
40 
41  //! Create a new TGA level integrator.
42  //! \param a_grids The disjoint box layout on which the level integrator is defined.
43  //! \param a_refRat The refinement ratios for the boxes.
44  //! \param a_level0Domain The coarsest grid level defining the problem domain.
45  //! \param a_opFact A factory that generates level diffusion operators.
46  //! \param a_solver A multigrid solver.
48  const Vector<int>& a_refRat,
49  const ProblemDomain& a_level0Domain,
52  :BaseLevelTGA<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>(a_grids, a_refRat, a_level0Domain, a_opFact, a_solver)
53  {
54  m_isEBLGSet = false;
55  }
56 
57  //! Destructor
58  virtual ~EBLevelTGA()
59  {
60  }
61 
62  //! Computes the time-centered diffusive term for explicit updating.
63  //! This integrates the diffusion equation for a new value of a quantity phi
64  //! and then computes the diffusion term using a finite difference stencil
65  //! in time.
66  //! \param a_diffusionTerm The term computed by this method.
67  //! \param a_phiOld The value of phi at the beginning of the time step,
68  //! at the current grid level.
69  //! \param a_src The source term in the diffusion equation at the current
70  //! grid level.
71  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
72  //! adjacent grid level (or NULL if there is none).
73  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
74  //! adjacent grid level (or NULL if there is none).
75  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
76  //! old value of phi on the coarser adjacent grid
77  //! level (or NULL if there is none).
78  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
79  //! updated value of phi on the coarser adjacent grid
80  //! level (or NULL if there is none).
81  //! \param a_oldTime The time at the beginning of the integration step on
82  //! the current grid level.
83  //! \param a_crseOldTime The time at the beginning of the integration step
84  //! on the coarser adjacent grid level.
85  //! \param a_crseNewTime The time at the end of the integration step
86  //! on the coarser adjacent grid level.
87  //! \param a_dt The time step on the current grid level.
88  //! \param a_level The current grid level.
89  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
90  //! be set to zero. Otherwise, \a a_phiOld will be used.
91  void computeDiffusion(LevelData<EBCellFAB>& a_DiffusiveTerm,
92  LevelData<EBCellFAB>& a_phiOld,
93  LevelData<EBCellFAB>& a_src,
94  EBFluxRegister* a_FineFluxRegPtr,
95  EBFluxRegister* a_crseFluxRegPtr,
96  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
97  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
98  Real a_oldTime,
99  Real a_crseOldTime,
100  Real a_crseNewTime,
101  Real a_dt,
102  int a_level,
103  bool a_zeroPhi = true,
104  bool a_rhsAlreadyKappaWeighted = false
105  );
106 
107  //! Integrates the diffusion equation, storing the result in \a a_phiNew.
108  //! \param a_phiNew The new value of phi at the current grid level.
109  //! \param a_phiOld The value of phi at the beginning of the time step,
110  //! at the current grid level.
111  //! \param a_src The source term in the diffusion equation at the current
112  //! grid level.
113  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
114  //! adjacent grid level (or NULL if there is none).
115  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
116  //! adjacent grid level (or NULL if there is none).
117  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
118  //! old value of phi on the coarser adjacent grid
119  //! level (or NULL if there is none).
120  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
121  //! updated value of phi on the coarser adjacent grid
122  //! level (or NULL if there is none).
123  //! \param a_oldTime The time at the beginning of the integration step on
124  //! the current grid level.
125  //! \param a_crseOldTime The time at the beginning of the integration step
126  //! on the coarser adjacent grid level.
127  //! \param a_crseNewTime The time at the end of the integration step
128  //! on the coarser adjacent grid level.
129  //! \param a_dt The time step on the current grid level.
130  //! \param a_level The current grid level.
131  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
132  //! be set to zero. Otherwise, \a a_phiOld will be used.
133  void updateSoln(LevelData<EBCellFAB>& a_phiNew,
134  LevelData<EBCellFAB>& a_phiOld,
135  LevelData<EBCellFAB>& a_src,
136  EBFluxRegister* a_fineFluxRegPtr,
137  EBFluxRegister* a_crseFluxRegPtr,
138  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
139  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
140  Real a_oldTime,
141  Real a_crseOldTime,
142  Real a_crseNewTime,
143  Real a_dt,
144  int a_level,
145  bool a_zeroPhi = true,
146  bool a_rhsAlreadyKappaWeighted = false,
147  int a_fluxStartComponent = 0);
148 
150  const DisjointBoxLayout& a_grids,
151  int a_lev);
152 
154  {
155  m_isEBLGSet = true;
156  m_eblg = a_eblg;
157  }
158 
159 protected:
162 
163 
164 };
165 
166 
167 //! \class EBLevelBackwardEuler
168 //! This class implements an imlplicit integrator for a diffusion
169 //! equation using an approach devised by Twizzell, Gumel, and Arigu.
170 class EBLevelBackwardEuler : public BaseLevelBackwardEuler<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>
171 {
172 
173 public:
174 
175  //! Create a new BackwardEuler level integrator.
176  //! \param a_grids The disjoint box layout on which the level integrator is defined.
177  //! \param a_refRat The refinement ratios for the boxes.
178  //! \param a_level0Domain The coarsest grid level defining the problem domain.
179  //! \param a_opFact A factory that generates level diffusion operators.
180  //! \param a_solver A multigrid solver.
182  const Vector<int>& a_refRat,
183  const ProblemDomain& a_level0Domain,
185  const RefCountedPtr<AMRMultiGrid<LevelData<EBCellFAB> > >& a_solver)
186  :BaseLevelBackwardEuler<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>(a_grids, a_refRat, a_level0Domain, a_opFact, a_solver)
187  {
188  m_isEBLGSet = false;
189  }
190 
191  //! Destructor
193  {
194  }
195 
196  //! Computes the time-centered diffusive term for explicit updating.
197  //! This integrates the diffusion equation for a new value of a quantity phi
198  //! and then computes the diffusion term using a finite difference stencil
199  //! in time.
200  //! \param a_diffusionTerm The term computed by this method.
201  //! \param a_phiOld The value of phi at the beginning of the time step,
202  //! at the current grid level.
203  //! \param a_src The source term in the diffusion equation at the current
204  //! grid level.
205  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
206  //! adjacent grid level (or NULL if there is none).
207  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
208  //! adjacent grid level (or NULL if there is none).
209  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
210  //! old value of phi on the coarser adjacent grid
211  //! level (or NULL if there is none).
212  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
213  //! updated value of phi on the coarser adjacent grid
214  //! level (or NULL if there is none).
215  //! \param a_oldTime The time at the beginning of the integration step on
216  //! the current grid level.
217  //! \param a_crseOldTime The time at the beginning of the integration step
218  //! on the coarser adjacent grid level.
219  //! \param a_crseNewTime The time at the end of the integration step
220  //! on the coarser adjacent grid level.
221  //! \param a_dt The time step on the current grid level.
222  //! \param a_level The current grid level.
223  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
224  //! be set to zero. Otherwise, \a a_phiOld will be used.
225  void computeDiffusion(LevelData<EBCellFAB>& a_DiffusiveTerm,
226  LevelData<EBCellFAB>& a_phiOld,
227  LevelData<EBCellFAB>& a_src,
228  EBFluxRegister* a_FineFluxRegPtr,
229  EBFluxRegister* a_crseFluxRegPtr,
230  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
231  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
232  Real a_oldTime,
233  Real a_crseOldTime,
234  Real a_crseNewTime,
235  Real a_dt,
236  int a_level,
237  bool a_zeroPhi = true,
238  bool a_rhsAlreadyKappaWeighted = false
239  );
240 
241  //! Integrates the diffusion equation, storing the result in \a a_phiNew.
242  //! \param a_phiNew The new value of phi at the current grid level.
243  //! \param a_phiOld The value of phi at the beginning of the time step,
244  //! at the current grid level.
245  //! \param a_src The source term in the diffusion equation at the current
246  //! grid level.
247  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
248  //! adjacent grid level (or NULL if there is none).
249  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
250  //! adjacent grid level (or NULL if there is none).
251  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
252  //! old value of phi on the coarser adjacent grid
253  //! level (or NULL if there is none).
254  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
255  //! updated value of phi on the coarser adjacent grid
256  //! level (or NULL if there is none).
257  //! \param a_oldTime The time at the beginning of the integration step on
258  //! the current grid level.
259  //! \param a_crseOldTime The time at the beginning of the integration step
260  //! on the coarser adjacent grid level.
261  //! \param a_crseNewTime The time at the end of the integration step
262  //! on the coarser adjacent grid level.
263  //! \param a_dt The time step on the current grid level.
264  //! \param a_level The current grid level.
265  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
266  //! be set to zero. Otherwise, \a a_phiOld will be used.
267  void updateSoln(LevelData<EBCellFAB>& a_phiNew,
268  LevelData<EBCellFAB>& a_phiOld,
269  LevelData<EBCellFAB>& a_src,
270  EBFluxRegister* a_fineFluxRegPtr,
271  EBFluxRegister* a_crseFluxRegPtr,
272  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
273  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
274  Real a_oldTime,
275  Real a_crseOldTime,
276  Real a_crseNewTime,
277  Real a_dt,
278  int a_level,
279  bool a_zeroPhi = true,
280  bool a_rhsAlreadyKappaWeighted = false,
281  int a_fluxStartComponent = 0);
282 
283 
285  {
286  m_isEBLGSet = true;
287  m_eblg = a_eblg;
288  }
289 
290 protected:
293 
294 
295 };
296 
297 //! \class EBLevelCrankNicolson
298 //! This class implements an imlplicit integrator for a diffusion
299 //! equation using an approach devised by Twizzell, Gumel, and Arigu.
300 class EBLevelCrankNicolson : public BaseLevelCrankNicolson<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>
301 {
302 
303 public:
304 
305  //! Create a new CrankNicolson level integrator.
306  //! \param a_grids The disjoint box layout on which the level integrator is defined.
307  //! \param a_refRat The refinement ratios for the boxes.
308  //! \param a_level0Domain The coarsest grid level defining the problem domain.
309  //! \param a_opFact A factory that generates level diffusion operators.
310  //! \param a_solver A multigrid solver.
312  const Vector<int>& a_refRat,
313  const ProblemDomain& a_level0Domain,
315  const RefCountedPtr<AMRMultiGrid<LevelData<EBCellFAB> > >& a_solver)
316  :BaseLevelCrankNicolson<LevelData<EBCellFAB>, EBFluxFAB, EBFluxRegister>(a_grids, a_refRat, a_level0Domain, a_opFact, a_solver)
317  {
318  m_isEBLGSet = false;
319  }
320 
321  //! Destructor
323  {
324  }
325 
326  //! Computes the time-centered diffusive term for explicit updating.
327  //! This integrates the diffusion equation for a new value of a quantity phi
328  //! and then computes the diffusion term using a finite difference stencil
329  //! in time.
330  //! \param a_diffusionTerm The term computed by this method.
331  //! \param a_phiOld The value of phi at the beginning of the time step,
332  //! at the current grid level.
333  //! \param a_src The source term in the diffusion equation at the current
334  //! grid level.
335  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
336  //! adjacent grid level (or NULL if there is none).
337  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
338  //! adjacent grid level (or NULL if there is none).
339  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
340  //! old value of phi on the coarser adjacent grid
341  //! level (or NULL if there is none).
342  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
343  //! updated value of phi on the coarser adjacent grid
344  //! level (or NULL if there is none).
345  //! \param a_oldTime The time at the beginning of the integration step on
346  //! the current grid level.
347  //! \param a_crseOldTime The time at the beginning of the integration step
348  //! on the coarser adjacent grid level.
349  //! \param a_crseNewTime The time at the end of the integration step
350  //! on the coarser adjacent grid level.
351  //! \param a_dt The time step on the current grid level.
352  //! \param a_level The current grid level.
353  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
354  //! be set to zero. Otherwise, \a a_phiOld will be used.
355  void computeDiffusion(LevelData<EBCellFAB>& a_DiffusiveTerm,
356  LevelData<EBCellFAB>& a_phiOld,
357  LevelData<EBCellFAB>& a_src,
358  EBFluxRegister* a_FineFluxRegPtr,
359  EBFluxRegister* a_crseFluxRegPtr,
360  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
361  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
362  Real a_oldTime,
363  Real a_crseOldTime,
364  Real a_crseNewTime,
365  Real a_dt,
366  int a_level,
367  bool a_zeroPhi = true,
368  bool a_rhsAlreadyKappaWeighted = false
369  );
370 
371  //! Integrates the diffusion equation, storing the result in \a a_phiNew.
372  //! \param a_phiNew The new value of phi at the current grid level.
373  //! \param a_phiOld The value of phi at the beginning of the time step,
374  //! at the current grid level.
375  //! \param a_src The source term in the diffusion equation at the current
376  //! grid level.
377  //! \param a_fineFluxRegPtr A pointer to the flux register at the finer
378  //! adjacent grid level (or NULL if there is none).
379  //! \param a_crseFluxRegPtr A pointer to the flux register at the coarser
380  //! adjacent grid level (or NULL if there is none).
381  //! \param a_crsePhiOldPtr A pointer to the LevelData object holding the
382  //! old value of phi on the coarser adjacent grid
383  //! level (or NULL if there is none).
384  //! \param a_crsePhiNewPtr A pointer to the LevelData object holding the
385  //! updated value of phi on the coarser adjacent grid
386  //! level (or NULL if there is none).
387  //! \param a_oldTime The time at the beginning of the integration step on
388  //! the current grid level.
389  //! \param a_crseOldTime The time at the beginning of the integration step
390  //! on the coarser adjacent grid level.
391  //! \param a_crseNewTime The time at the end of the integration step
392  //! on the coarser adjacent grid level.
393  //! \param a_dt The time step on the current grid level.
394  //! \param a_level The current grid level.
395  //! \param a_zeroPhi If this flag is true, the initial estimate of phi will
396  //! be set to zero. Otherwise, \a a_phiOld will be used.
397  void updateSoln(LevelData<EBCellFAB>& a_phiNew,
398  LevelData<EBCellFAB>& a_phiOld,
399  LevelData<EBCellFAB>& a_src,
400  EBFluxRegister* a_fineFluxRegPtr,
401  EBFluxRegister* a_crseFluxRegPtr,
402  const LevelData<EBCellFAB>* a_crsePhiOldPtr,
403  const LevelData<EBCellFAB>* a_crsePhiNewPtr,
404  Real a_oldTime,
405  Real a_crseOldTime,
406  Real a_crseNewTime,
407  Real a_dt,
408  int a_level,
409  bool a_zeroPhi = true,
410  bool a_rhsAlreadyKappaWeighted = false,
411  int a_fluxStartComponent = 0);
412 
413 
415  {
416  m_isEBLGSet = true;
417  m_eblg = a_eblg;
418  }
419 
420 protected:
423 
424 
425 };
426 
427 #include "NamespaceFooter.H"
428 #endif
429 
430 
EBLevelBackwardEuler(const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
Definition: EBLevelTGA.H:181
A reference-counting handle class.
Definition: RefCountedPtr.H:173
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:141
Definition: EBLevelTGA.H:300
void setEBLG(Vector< EBLevelGrid > &a_eblg)
Definition: EBLevelTGA.H:153
virtual ~EBLevelBackwardEuler()
Destructor.
Definition: EBLevelTGA.H:192
EBLevelTGA(const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
Definition: EBLevelTGA.H:47
Vector< EBLevelGrid > m_eblg
Definition: EBLevelTGA.H:292
Definition: BaseLevelBackwardEuler.H:44
void setEBLG(Vector< EBLevelGrid > &a_eblg)
Definition: EBLevelTGA.H:284
bool m_isEBLGSet
Definition: EBLevelTGA.H:421
bool m_isEBLGSet
Definition: EBLevelTGA.H:291
Definition: AMRMultiGrid.H:308
Definition: EBLevelTGA.H:170
bool m_isEBLGSet
Definition: EBLevelTGA.H:160
A EBFaceFAB-like container for edge-centered fluxes.
Definition: EBFluxFAB.H:25
virtual ~EBLevelCrankNicolson()
Destructor.
Definition: EBLevelTGA.H:322
Definition: EBCellFAB.H:29
double Real
Definition: REAL.H:33
Definition: BaseLevelTGA.H:44
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:30
void updateSoln(LevelData< EBCellFAB > &a_phiNew, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, EBFluxRegister *a_fineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false, int a_fluxStartComponent=0)
EBLevelCrankNicolson(const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
Definition: EBLevelTGA.H:311
Vector< EBLevelGrid > m_eblg
Definition: EBLevelTGA.H:161
void setSourceGhostCells(LevelData< EBCellFAB > &a_src, const DisjointBoxLayout &a_grids, int a_lev)
void setEBLG(Vector< EBLevelGrid > &a_eblg)
Definition: EBLevelTGA.H:414
Vector< EBLevelGrid > m_eblg
Definition: EBLevelTGA.H:422
virtual ~EBLevelTGA()
Destructor.
Definition: EBLevelTGA.H:58
Definition: AMRMultiGrid.H:233
EBFluxRegister-A class to encapsulate a levels worth of flux registers.
Definition: EBFluxRegister.H:37
void computeDiffusion(LevelData< EBCellFAB > &a_DiffusiveTerm, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, EBFluxRegister *a_FineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false)
Definition: BaseLevelCrankNicolson.H:44
Definition: EBLevelTGA.H:36