Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 namespace Nektar
41 {
43 
46  : UnsteadySystem(pSession)
47  {
48  }
49 
50  /**
51  * @brief Initialisation object for the unsteady diffusion problem.
52  */
54  {
55  UnsteadySystem::v_InitObject();
56 
57  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
58  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
59 
60  m_session->MatchSolverInfo(
61  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
62 
64  {
65  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
66  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
67  }
68 
69  int npoints = m_fields[0]->GetNpoints();
70 
71  if(m_session->DefinesParameter("d00"))
72  {
74  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
75  }
76  if(m_session->DefinesParameter("d11"))
77  {
79  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
80  }
81  if(m_session->DefinesParameter("d22"))
82  {
84  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
85  }
86 
87  switch (m_projectionType)
88  {
90  {
91  std::string diffName;
92 
93  // Do not forwards transform initial condition
94  m_homoInitialFwd = false;
95 
96  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
98  CreateInstance(diffName, diffName);
99  m_diffusion->SetFluxVector(&UnsteadyDiffusion::
100  GetFluxVector, this);
101  m_diffusion->InitObject(m_session, m_fields);
102  break;
103  }
104 
107  {
108  // In case of Galerkin explicit diffusion gives an error
110  {
111  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
112  }
113  // In case of Galerkin implicit diffusion: do nothing
114  }
115  }
116 
117 
119  {
122  }
123  else
124  {
127  }
128  }
129 
130  /**
131  * @brief Unsteady diffusion problem destructor.
132  */
134  {
135  }
136 
138  {
139  UnsteadySystem::v_GenerateSummary(s);
140  if(m_useSpecVanVisc)
141  {
142  stringstream ss;
143  ss << "SVV (cut off = " << m_sVVCutoffRatio
144  << ", coeff = " << m_sVVDiffCoeff << ")";
145  AddSummaryItem(s, "Smoothing", ss.str());
146  }
147  }
148 
149 
150  /* @brief Compute the right-hand side for the unsteady diffusion problem.
151  *
152  * @param inarray Given fields.
153  * @param outarray Calculated solution.
154  * @param time Time.
155  */
157  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
158  Array<OneD, Array<OneD, NekDouble> > &outarray,
159  const NekDouble time)
160  {
161  // Number of fields (variables of the problem)
162  int nVariables = inarray.num_elements();
163 
164  // RHS computation using the new advection base class
165  m_diffusion->Diffuse(nVariables,
166  m_fields,
167  inarray,
168  outarray);
169  }
170 
171  /**
172  * @brief Compute the projection for the unsteady diffusion problem.
173  *
174  * @param inarray Given fields.
175  * @param outarray Calculated solution.
176  * @param time Time.
177  */
179  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
180  Array<OneD, Array<OneD, NekDouble> > &outarray,
181  const NekDouble time)
182  {
183  int i;
184  int nvariables = inarray.num_elements();
185  SetBoundaryConditions(time);
186 
187  switch(m_projectionType)
188  {
190  {
191  // Just copy over array
192  int npoints = GetNpoints();
193 
194  for(i = 0; i < nvariables; ++i)
195  {
196  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
197  }
198  break;
199  }
202  {
204 
205  for(i = 0; i < nvariables; ++i)
206  {
207  m_fields[i]->FwdTrans(inarray[i], coeffs);
208  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
209  }
210  break;
211  }
212  default:
213  {
214  ASSERTL0(false, "Unknown projection scheme");
215  break;
216  }
217  }
218  }
219 
220  /**
221  * @brief Implicit solution of the unsteady diffusion problem.
222  */
224  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
225  Array<OneD, Array<OneD, NekDouble> > &outarray,
226  const NekDouble time,
227  const NekDouble lambda)
228  {
230 
231  int nvariables = inarray.num_elements();
232  int npoints = m_fields[0]->GetNpoints();
233  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
234  factors[StdRegions::eFactorTau] = 1.0;
235 
236  if(m_useSpecVanVisc)
237  {
240  }
241 
242  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
243  // inarray = input: \hat{rhs} -> output: \hat{Y}
244  // outarray = output: nabla^2 \hat{Y}
245  // where \hat = modal coeffs
246  for (int i = 0; i < nvariables; ++i)
247  {
248  // Multiply 1.0/timestep/lambda
249  Vmath::Smul(npoints,
250  -factors[StdRegions::eFactorLambda],
251  inarray[i], 1,
252  m_fields[i]->UpdatePhys(), 1);
253 
254  // Solve a system of equations with Helmholtz solver
255  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(),
256  m_fields[i]->UpdateCoeffs(),
257  NullFlagList,
258  factors,
259  m_varcoeff);
260 
261  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
262  m_fields[i]->UpdatePhys());
263 
264  m_fields[i]->SetPhysState(false);
265 
266  // The solution is Y[i]
267  outarray[i] = m_fields[i]->GetPhys();
268  }
269  }
270 
271  /**
272  * @brief Return the flux vector for the unsteady diffusion problem.
273  */
275  const int i,
276  const int j,
277  const Array<OneD, Array<OneD, NekDouble> > &physfield,
278  Array<OneD, Array<OneD, NekDouble> > &derivatives,
280  {
281  for(int k = 0; k < flux.num_elements(); ++k)
282  {
283  Vmath::Zero(GetNpoints(), flux[k], 1);
284  }
285  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
286  }
287 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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.
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:42
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: CellModel.h:51
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
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
UnsteadyDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession)
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.
static std::string className
Name of class.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
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:199
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.
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:359
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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