Nektar++
UnsteadyDiffusion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyDiffusion.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 diffusion solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 #include <iostream>
38 #include <iomanip>
39 
40 #include <boost/core/ignore_unused.hpp>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46  string UnsteadyDiffusion::className = GetEquationSystemFactory().RegisterCreatorFunction("UnsteadyDiffusion", UnsteadyDiffusion::create);
47 
48  UnsteadyDiffusion::UnsteadyDiffusion(
51  : UnsteadySystem(pSession, pGraph)
52  {
53  }
54 
55  /**
56  * @brief Initialisation object for the unsteady diffusion problem.
57  */
59  {
61 
62  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
63  m_session->LoadParameter("epsilon", m_epsilon, 1.0);
64 
65  m_session->MatchSolverInfo(
66  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
67 
69  {
70  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
71  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
72  }
73 
74  int npoints = m_fields[0]->GetNpoints();
75 
76  if(m_session->DefinesParameter("d00"))
77  {
79  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
80  }
81  if(m_session->DefinesParameter("d11"))
82  {
84  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
85  }
86  if(m_session->DefinesParameter("d22"))
87  {
89  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
90  }
91 
92  switch (m_projectionType)
93  {
95  {
96  std::string diffName;
97 
98  // Do not forwards transform initial condition
99  m_homoInitialFwd = false;
100 
101  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
103  CreateInstance(diffName, diffName);
104  m_diffusion->SetFluxVector(&UnsteadyDiffusion::
105  GetFluxVector, this);
106  m_diffusion->InitObject(m_session, m_fields);
107  break;
108  }
109 
112  {
113  // In case of Galerkin explicit diffusion gives an error
115  {
116  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
117  }
118  // In case of Galerkin implicit diffusion: do nothing
119  }
120  }
121 
123  {
126  }
127  else
128  {
133  }
134  }
135 
136  /**
137  * @brief Unsteady diffusion problem destructor.
138  */
140  {
141  }
142 
144  {
146  if(m_useSpecVanVisc)
147  {
148  stringstream ss;
149  ss << "SVV (cut off = " << m_sVVCutoffRatio
150  << ", coeff = " << m_sVVDiffCoeff << ")";
151  AddSummaryItem(s, "Smoothing", ss.str());
152  }
153  }
154 
155 
156  /* @brief Compute the right-hand side for the unsteady diffusion problem.
157  *
158  * @param inarray Given fields.
159  * @param outarray Calculated solution.
160  * @param time Time.
161  */
163  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
164  Array<OneD, Array<OneD, NekDouble> > &outarray,
165  const NekDouble time)
166  {
167  boost::ignore_unused(time);
168 
169  // Number of fields (variables of the problem)
170  int nVariables = inarray.size();
171 
172  // RHS computation using the new advection base class
173  m_diffusion->Diffuse(nVariables,
174  m_fields,
175  inarray,
176  outarray);
177  }
178 
179  /**
180  * @brief Compute the projection for the unsteady diffusion problem.
181  *
182  * @param inarray Given fields.
183  * @param outarray Calculated solution.
184  * @param time Time.
185  */
187  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
188  Array<OneD, Array<OneD, NekDouble> > &outarray,
189  const NekDouble time)
190  {
191  int i;
192  int nvariables = inarray.size();
193  SetBoundaryConditions(time);
194 
195  switch(m_projectionType)
196  {
198  {
199  // Just copy over array
200  int npoints = GetNpoints();
201 
202  for(i = 0; i < nvariables; ++i)
203  {
204  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
205  }
206  break;
207  }
210  {
212 
213  for(i = 0; i < nvariables; ++i)
214  {
215  m_fields[i]->FwdTrans(inarray[i], coeffs);
216  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
217  }
218  break;
219  }
220  default:
221  {
222  ASSERTL0(false, "Unknown projection scheme");
223  break;
224  }
225  }
226  }
227 
228  /**
229  * @brief Implicit solution of the unsteady diffusion problem.
230  */
232  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
233  Array<OneD, Array<OneD, NekDouble> > &outarray,
234  const NekDouble time,
235  const NekDouble lambda)
236  {
237  boost::ignore_unused(time);
238 
240 
241  int nvariables = inarray.size();
242  int npoints = m_fields[0]->GetNpoints();
243  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
244  factors[StdRegions::eFactorTau] = 1.0;
245 
246  if(m_useSpecVanVisc)
247  {
250  }
251 
252  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
253  // inarray = input: \hat{rhs} -> output: \hat{Y}
254  // outarray = output: nabla^2 \hat{Y}
255  // where \hat = modal coeffs
256  for (int i = 0; i < nvariables; ++i)
257  {
258  // Multiply 1.0/timestep/lambda
259  Vmath::Smul(npoints,
260  -factors[StdRegions::eFactorLambda],
261  inarray[i], 1,
262  outarray[i], 1);
263 
264  // Solve a system of equations with Helmholtz solver
265  m_fields[i]->HelmSolve(outarray[i],
266  m_fields[i]->UpdateCoeffs(),
267  factors,
268  m_varcoeff);
269 
270  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
271  outarray[i]);
272 
273  m_fields[i]->SetPhysState(false);
274  }
275  }
276 
277  /**
278  * @brief Return the flux vector for the unsteady diffusion problem.
279  */
281  const Array<OneD, Array<OneD, NekDouble> > &inarray,
282  const Array<OneD, Array<OneD, Array<OneD, NekDouble> > >&qfield,
283  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >&viscousTensor)
284  {
285  boost::ignore_unused(inarray);
286 
287  unsigned int nDim = qfield.size();
288  unsigned int nConvectiveFields = qfield[0].size();
289  unsigned int nPts = qfield[0][0].size();
290 
291  for (unsigned int j = 0; j < nDim; ++j)
292  {
293  for (unsigned int i = 0; i < nConvectiveFields; ++i)
294  {
295  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1,
296  viscousTensor[j][i], 1 );
297  }
298  }
299  }
300 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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 GetNpoints()
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.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
virtual 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 GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
Return the flux vector for the unsteady diffusion problem.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
StdRegions::VarCoeffMap m_varcoeff
virtual void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
virtual void v_InitObject()
Initialisation object for 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.
virtual ~UnsteadyDiffusion()
Destructor.
SolverUtils::DiffusionSharedPtr m_diffusion
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:225
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199