Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Unsteady diffusion solve routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <iostream>
38 #include <iomanip>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44  string UnsteadyDiffusion::className = GetEquationSystemFactory().RegisterCreatorFunction("UnsteadyDiffusion", UnsteadyDiffusion::create);
45 
46  UnsteadyDiffusion::UnsteadyDiffusion(
48  : UnsteadySystem(pSession)
49  {
50  }
51 
52  /**
53  * @brief Initialisation object for the unsteady diffusion problem.
54  */
56  {
58 
59  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
60  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
61 
62  m_session->MatchSolverInfo(
63  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
64 
66  {
67  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
68  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
69  }
70 
71  int npoints = m_fields[0]->GetNpoints();
72 
73  if(m_session->DefinesParameter("d00"))
74  {
76  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
77  }
78  if(m_session->DefinesParameter("d11"))
79  {
81  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
82  }
83  if(m_session->DefinesParameter("d22"))
84  {
86  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
87  }
88 
89  switch (m_projectionType)
90  {
92  {
93  std::string diffName;
94 
95  // Do not forwards transform initial condition
96  m_homoInitialFwd = false;
97 
98  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
100  CreateInstance(diffName, diffName);
101  m_diffusion->SetFluxVector(&UnsteadyDiffusion::
102  GetFluxVector, this);
103  m_diffusion->InitObject(m_session, m_fields);
104  break;
105  }
106 
109  {
110  // In case of Galerkin explicit diffusion gives an error
112  {
113  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
114  }
115  // In case of Galerkin implicit diffusion: do nothing
116  }
117  }
118 
119 
121  {
124  }
125  else
126  {
129  }
130  }
131 
132  /**
133  * @brief Unsteady diffusion problem destructor.
134  */
136  {
137  }
138 
140  {
142  if(m_useSpecVanVisc)
143  {
144  stringstream ss;
145  ss << "SVV (cut off = " << m_sVVCutoffRatio
146  << ", coeff = " << m_sVVDiffCoeff << ")";
147  AddSummaryItem(s, "Smoothing", ss.str());
148  }
149  }
150 
151 
152  /* @brief Compute the right-hand side for the unsteady diffusion problem.
153  *
154  * @param inarray Given fields.
155  * @param outarray Calculated solution.
156  * @param time Time.
157  */
159  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
160  Array<OneD, Array<OneD, NekDouble> > &outarray,
161  const NekDouble time)
162  {
163  // Number of fields (variables of the problem)
164  int nVariables = inarray.num_elements();
165 
166  // RHS computation using the new advection base class
167  m_diffusion->Diffuse(nVariables,
168  m_fields,
169  inarray,
170  outarray);
171  }
172 
173  /**
174  * @brief Compute the projection for the unsteady diffusion problem.
175  *
176  * @param inarray Given fields.
177  * @param outarray Calculated solution.
178  * @param time Time.
179  */
181  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
182  Array<OneD, Array<OneD, NekDouble> > &outarray,
183  const NekDouble time)
184  {
185  int i;
186  int nvariables = inarray.num_elements();
187  SetBoundaryConditions(time);
188 
189  switch(m_projectionType)
190  {
192  {
193  // Just copy over array
194  int npoints = GetNpoints();
195 
196  for(i = 0; i < nvariables; ++i)
197  {
198  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
199  }
200  break;
201  }
204  {
206 
207  for(i = 0; i < nvariables; ++i)
208  {
209  m_fields[i]->FwdTrans(inarray[i], coeffs);
210  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
211  }
212  break;
213  }
214  default:
215  {
216  ASSERTL0(false, "Unknown projection scheme");
217  break;
218  }
219  }
220  }
221 
222  /**
223  * @brief Implicit solution of the unsteady diffusion problem.
224  */
226  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
227  Array<OneD, Array<OneD, NekDouble> > &outarray,
228  const NekDouble time,
229  const NekDouble lambda)
230  {
232 
233  int nvariables = inarray.num_elements();
234  int npoints = m_fields[0]->GetNpoints();
235  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
236  factors[StdRegions::eFactorTau] = 1.0;
237 
238  if(m_useSpecVanVisc)
239  {
242  }
243 
244  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
245  // inarray = input: \hat{rhs} -> output: \hat{Y}
246  // outarray = output: nabla^2 \hat{Y}
247  // where \hat = modal coeffs
248  for (int i = 0; i < nvariables; ++i)
249  {
250  // Multiply 1.0/timestep/lambda
251  Vmath::Smul(npoints,
252  -factors[StdRegions::eFactorLambda],
253  inarray[i], 1,
254  m_fields[i]->UpdatePhys(), 1);
255 
256  // Solve a system of equations with Helmholtz solver
257  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(),
258  m_fields[i]->UpdateCoeffs(),
259  NullFlagList,
260  factors,
261  m_varcoeff);
262 
263  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
264  m_fields[i]->UpdatePhys());
265 
266  m_fields[i]->SetPhysState(false);
267 
268  // The solution is Y[i]
269  outarray[i] = m_fields[i]->GetPhys();
270  }
271  }
272 
273  /**
274  * @brief Return the flux vector for the unsteady diffusion problem.
275  */
277  const int i,
278  const int j,
279  const Array<OneD, Array<OneD, NekDouble> > &physfield,
280  Array<OneD, Array<OneD, NekDouble> > &derivatives,
282  {
283  for(int k = 0; k < flux.num_elements(); ++k)
284  {
285  Vmath::Zero(GetNpoints(), flux[k], 1);
286  }
287  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
288  }
289 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void GetFluxVector(const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
Return the flux vector for the unsteady diffusion problem.
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:47
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:42
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
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.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
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:213
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:50
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.
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()
SolverUtils::DiffusionSharedPtr m_diffusion
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
virtual void v_InitObject()
Initialisation object for the unsteady diffusion problem.
StdRegions::VarCoeffMap m_varcoeff
static FlagList NullFlagList
An empty flag list.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215