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