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 <iomanip>
36 #include <iostream>
37 
38 #include <boost/core/ignore_unused.hpp>
39 
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47 string UnsteadyReactionDiffusion::className =
49  "UnsteadyReactionDiffusion", UnsteadyReactionDiffusion::create);
50 
51 UnsteadyReactionDiffusion::UnsteadyReactionDiffusion(
54  : UnsteadySystem(pSession, pGraph)
55 {
56 }
57 
58 /**
59  * @brief Initialisation object for the unsteady reaction-diffusion problem.
60  */
62 {
63  UnsteadySystem::v_InitObject(DeclareFields);
64 
65  ASSERTL0(m_intScheme->GetIntegrationSchemeType() == LibUtilities::eIMEX,
66  "Reaction-diffusion requires an implicit-explicit timestepping"
67  " (e.g. IMEXOrder2)");
69  "Reaction-diffusion requires use of continuous Galerkin"
70  "projection.");
71 
72  // Load diffusion parameter
73  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
74 
75  // Forcing terms
76  m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
77  m_fields, m_fields.size());
78 
82  this);
83 }
84 
85 /**
86  * @brief Unsteady diffusion problem destructor.
87  */
89 {
90 }
91 
92 /**
93  * @brief Compute the right-hand side for the unsteady reaction diffusion
94  * problem.
95  *
96  * @param inarray Given fields.
97  * @param outarray Calculated solution.
98  * @param time Time.
99  */
101  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
102  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
103 {
104  // RHS should be set to zero.
105  for (int i = 0; i < outarray.size(); ++i)
106  {
107  Vmath::Zero(outarray[i].size(), &outarray[i][0], 1);
108  }
109 
110  // Add forcing terms for reaction.
111  for (auto &x : m_forcing)
112  {
113  // set up non-linear terms
114  x->Apply(m_fields, inarray, outarray, time);
115  }
116 }
117 
118 /**
119  * @brief Compute the projection for the unsteady diffusion problem.
120  *
121  * @param inarray Given fields.
122  * @param outarray Calculated solution.
123  * @param time Time.
124  */
126  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
127  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
128 {
129  int i;
130  int nvariables = inarray.size();
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(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, const NekDouble time,
148  const NekDouble lambda)
149 {
150  boost::ignore_unused(time);
151 
153 
154  int nvariables = inarray.size();
155  int npoints = m_fields[0]->GetNpoints();
156  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
157 
158  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i] inarray = input:
159  // \hat{rhs} -> output: \hat{Y} outarray = output: nabla^2 \hat{Y} where
160  // \hat = modal coeffs
161  for (int i = 0; i < nvariables; ++i)
162  {
163  // Multiply 1.0/timestep/lambda
164  Vmath::Smul(npoints, -factors[StdRegions::eFactorLambda], inarray[i], 1,
165  outarray[i], 1);
166 
167  // Solve a system of equations with Helmholtz solver
168  m_fields[i]->HelmSolve(outarray[i], m_fields[i]->UpdateCoeffs(),
169  factors);
170  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
171  m_fields[i]->SetPhysState(false);
172  }
173 }
174 
175 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
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:120
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
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 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.
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.
virtual void v_InitObject(bool DeclareFields=true) override
Initialisation object for the unsteady reaction-diffusion problem.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eIMEX
Implicit Explicit General Linear Method.
EquationSystemFactory & GetEquationSystemFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492