Nektar++
UnsteadyReactionDiffusion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyReactionDiffusion.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Unsteady reaction-diffusion solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 #include <iomanip>
37 
38 #include <boost/core/ignore_unused.hpp>
39 
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 string UnsteadyReactionDiffusion::className = GetEquationSystemFactory().
47  RegisterCreatorFunction("UnsteadyReactionDiffusion",
48  UnsteadyReactionDiffusion::create);
49 
50 UnsteadyReactionDiffusion::UnsteadyReactionDiffusion(
53  : UnsteadySystem(pSession, pGraph)
54 {
55 }
56 
57 /**
58  * @brief Initialisation object for the unsteady reaction-diffusion problem.
59  */
61 {
63 
64  ASSERTL0(m_intScheme->GetIntegrationSchemeType() == LibUtilities::eIMEX,
65  "Reaction-diffusion requires an implicit-explicit timestepping"
66  " (e.g. IMEXOrder2)");
68  "Reaction-diffusion requires use of continuous Galerkin"
69  "projection.");
70 
71  // Load diffusion parameter
72  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
73 
74  // Forcing terms
75  m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
76  m_fields, m_fields.num_elements());
77 
81 }
82 
83 /**
84  * @brief Unsteady diffusion problem destructor.
85  */
87 {
88 }
89 
90 /**
91  * @brief Compute the right-hand side for the unsteady reaction diffusion
92  * problem.
93  *
94  * @param inarray Given fields.
95  * @param outarray Calculated solution.
96  * @param time Time.
97  */
99  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
100  Array<OneD, Array<OneD, NekDouble> > &outarray,
101  const NekDouble time)
102 {
103  // RHS should be set to zero.
104  for (int i = 0; i < outarray.num_elements(); ++i)
105  {
106  Vmath::Zero(outarray[i].num_elements(), &outarray[i][0], 1);
107  }
108 
109  // Add forcing terms for reaction.
110  for (auto &x : m_forcing)
111  {
112  // set up non-linear terms
113  x->Apply(m_fields, inarray, outarray, time);
114  }
115 }
116 
117 /**
118  * @brief Compute the projection for the unsteady diffusion problem.
119  *
120  * @param inarray Given fields.
121  * @param outarray Calculated solution.
122  * @param time Time.
123  */
125  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
126  Array<OneD, Array<OneD, NekDouble> > &outarray,
127  const NekDouble time)
128 {
129  int i;
130  int nvariables = inarray.num_elements();
131  SetBoundaryConditions(time);
132 
134 
135  for(i = 0; i < nvariables; ++i)
136  {
137  m_fields[i]->FwdTrans(inarray[i], coeffs);
138  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
139  }
140 }
141 
142 /**
143  * @brief Implicit solution of the unsteady diffusion problem.
144  */
146  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
147  Array<OneD, Array<OneD, NekDouble> > &outarray,
148  const NekDouble time,
149  const NekDouble lambda)
150 {
151  boost::ignore_unused(time);
152 
154 
155  int nvariables = inarray.num_elements();
156  int npoints = m_fields[0]->GetNpoints();
157  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
158 
159  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i] inarray = input:
160  // \hat{rhs} -> output: \hat{Y} outarray = output: nabla^2 \hat{Y} where
161  // \hat = modal coeffs
162  for (int i = 0; i < nvariables; ++i)
163  {
164  // Multiply 1.0/timestep/lambda
165  Vmath::Smul(npoints, -factors[StdRegions::eFactorLambda],
166  inarray[i], 1, outarray[i], 1);
167 
168  // Solve a system of equations with Helmholtz solver
169  m_fields[i]->HelmSolve(
170  outarray[i], m_fields[i]->UpdateCoeffs(), NullFlagList, factors);
171  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
172  m_fields[i]->SetPhysState(false);
173  }
174 }
175 
176 }
virtual void v_InitObject()
Initialisation object for the unsteady reaction-diffusion problem.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Implicit solution of the unsteady diffusion problem.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection for the unsteady diffusion problem.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
Implicit Explicit General Linear Method.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:85
Base class for unsteady solvers.
double NekDouble
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side for the unsteady reaction diffusion problem.
EquationSystemFactory & GetEquationSystemFactory()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNcoeffs()
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FlagList NullFlagList
An empty flag list.