Chombo + EB  3.0
LevelTGA.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 _LEVELTGA_H_
12 #define _LEVELTGA_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 <FArrayBox.H>
23 #include "FluxBox.H"
24 #include "LevelFluxRegister.H"
25 #include <ProblemDomain.H>
26 #include "BaseLevelTGA.H"
27 #include "BaseLevelBackwardEuler.H"
28 #include "NamespaceHeader.H"
29 
30 
31 /// Implements TGA algorithm to update solution to heat equation
32 /** The LevelTGA class implements the Runge-Kutta-based approach
33  to solving the heat equation due to Twizell, Gumel, and Arigu.
34 */
35 class LevelTGA : public BaseLevelTGA<LevelData<FArrayBox>, FluxBox, LevelFluxRegister>
36 {
37 
38 public:
39 
40  /// full constructor
41  /**
42  */
44  const Vector<int>& a_refRat,
45  const ProblemDomain& a_level0Domain,
48  :BaseLevelTGA<LevelData<FArrayBox>, FluxBox, LevelFluxRegister>(a_grids, a_refRat, a_level0Domain, a_opFact, a_solver)
49  {
50  ;
51  }
52 
53  /// destructor
54  virtual ~LevelTGA()
55  {
56  }
57 
58  /// computes time-centered diffusive term for explicit updating
59  /** compute time-centered L(phi) for use in subsequent update
60  operations. In this case, we do solve for phiNew, then
61  subtract source and old phi back out to get L(phi).
62  */
63  void computeDiffusion(LevelData<FArrayBox>& a_DiffusiveTerm,
64  LevelData<FArrayBox>& a_phiOld,
65  LevelData<FArrayBox>& a_src,
66  LevelFluxRegister* a_FineFluxRegPtr,
67  LevelFluxRegister* a_crseFluxRegPtr,
68  const LevelData<FArrayBox>* a_crsePhiOldPtr,
69  const LevelData<FArrayBox>* a_crsePhiNewPtr,
70  Real a_oldTime,
71  Real a_crseOldTime,
72  Real a_crseNewTime,
73  Real a_dt,
74  int a_level);
75 
76  /// do solve
77  /** update phi to a_phiNew
78  */
79  void updateSoln(LevelData<FArrayBox>& a_phiNew,
80  LevelData<FArrayBox>& a_phiOld,
81  LevelData<FArrayBox>& a_src,
82  LevelFluxRegister* a_fineFluxRegPtr,
83  LevelFluxRegister* a_crseFluxRegPtr,
84  const LevelData<FArrayBox>* a_crsePhiOldPtr,
85  const LevelData<FArrayBox>* a_crsePhiNewPtr,
86  Real a_oldTime,
87  Real a_crseOldTime,
88  Real a_crseNewTime,
89  Real a_dt,
90  int a_level,
91  bool a_zeroPhi = true,
92  int a_fluxStartComponent = 0);
93 
95  const DisjointBoxLayout& a_grids,
96  int a_level);
97 protected:
98 
99 };
100 
101 /// Implements BackwardEuler algorithm to update solution to heat equation
102 /** The LevelBackwardEuler class implements the Runge-Kutta-based approach
103  to solving the heat equation due to Twizell, Gumel, and Arigu.
104 */
105 class LevelBackwardEuler : public BaseLevelBackwardEuler<LevelData<FArrayBox>, FluxBox, LevelFluxRegister>
106 {
107 
108 public:
109 
110  /// full constructor
111  /**
112  */
114  const Vector<int>& a_refRat,
115  const ProblemDomain& a_level0Domain,
117  const RefCountedPtr<AMRMultiGrid<LevelData<FArrayBox> > >& a_solver)
118  :BaseLevelBackwardEuler<LevelData<FArrayBox>, FluxBox, LevelFluxRegister>(a_grids, a_refRat, a_level0Domain, a_opFact, a_solver)
119  {
120  ;
121  }
122 
123  /// destructor
125  {
126  }
127 
128  /// computes time-centered diffusive term for explicit updating
129  /** compute time-centered L(phi) for use in subsequent update
130  operations. In this case, we do solve for phiNew, then
131  subtract source and old phi back out to get L(phi).
132  */
133  void computeDiffusion(LevelData<FArrayBox>& a_DiffusiveTerm,
134  LevelData<FArrayBox>& a_phiOld,
135  LevelData<FArrayBox>& a_src,
136  LevelFluxRegister* a_FineFluxRegPtr,
137  LevelFluxRegister* a_crseFluxRegPtr,
138  const LevelData<FArrayBox>* a_crsePhiOldPtr,
139  const LevelData<FArrayBox>* a_crsePhiNewPtr,
140  Real a_oldTime,
141  Real a_crseOldTime,
142  Real a_crseNewTime,
143  Real a_dt,
144  int a_level);
145 
146  /// do solve
147  /** update phi to a_phiNew
148  */
149  void updateSoln(LevelData<FArrayBox>& a_phiNew,
150  LevelData<FArrayBox>& a_phiOld,
151  LevelData<FArrayBox>& a_src,
152  LevelFluxRegister* a_fineFluxRegPtr,
153  LevelFluxRegister* a_crseFluxRegPtr,
154  const LevelData<FArrayBox>* a_crsePhiOldPtr,
155  const LevelData<FArrayBox>* a_crsePhiNewPtr,
156  Real a_oldTime,
157  Real a_crseOldTime,
158  Real a_crseNewTime,
159  Real a_dt,
160  int a_level,
161  bool a_zeroPhi = true,
162  int a_fluxStartComponent = 0);
163 
165  const DisjointBoxLayout& a_grids,
166  int a_level);
167 protected:
168 
169 };
170 
171 #include "NamespaceFooter.H"
172 #endif
173 
174 
A reference-counting handle class.
Definition: RefCountedPtr.H:66
virtual ~LevelTGA()
destructor
Definition: LevelTGA.H:54
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:130
LevelBackwardEuler(const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< FArrayBox > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< FArrayBox > > > &a_solver)
full constructor
Definition: LevelTGA.H:113
Definition: BaseLevelBackwardEuler.H:45
LevelTGA(const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< FArrayBox > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< FArrayBox > > > &a_solver)
full constructor
Definition: LevelTGA.H:43
void setSourceGhostCells(LevelData< FArrayBox > &a_src, const DisjointBoxLayout &a_grids, int a_level)
Definition: AMRMultiGrid.H:306
void computeDiffusion(LevelData< FArrayBox > &a_DiffusiveTerm, LevelData< FArrayBox > &a_phiOld, LevelData< FArrayBox > &a_src, LevelFluxRegister *a_FineFluxRegPtr, LevelFluxRegister *a_crseFluxRegPtr, const LevelData< FArrayBox > *a_crsePhiOldPtr, const LevelData< FArrayBox > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level)
computes time-centered diffusive term for explicit updating
A FArrayBox-like container for face-centered fluxes.
Definition: FluxBox.H:22
Implements BackwardEuler algorithm to update solution to heat equation.
Definition: LevelTGA.H:105
double Real
Definition: REAL.H:33
Definition: BaseLevelTGA.H:45
A BoxLayout that has a concept of disjointedness.
Definition: DisjointBoxLayout.H:31
Implements TGA algorithm to update solution to heat equation.
Definition: LevelTGA.H:35
LevelFluxRegister-A class to encapsulate a levels worth of flux registers.
Definition: LevelFluxRegister.H:29
void updateSoln(LevelData< FArrayBox > &a_phiNew, LevelData< FArrayBox > &a_phiOld, LevelData< FArrayBox > &a_src, LevelFluxRegister *a_fineFluxRegPtr, LevelFluxRegister *a_crseFluxRegPtr, const LevelData< FArrayBox > *a_crsePhiOldPtr, const LevelData< FArrayBox > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, int a_fluxStartComponent=0)
do solve
Definition: FArrayBox.H:44
Definition: AMRMultiGrid.H:231
virtual ~LevelBackwardEuler()
destructor
Definition: LevelTGA.H:124