Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
UnsteadyAdvectionDiffusion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyAdvectionDiffusion.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 advection-diffusion solve routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 
39 
40 namespace Nektar
41 {
44  "UnsteadyAdvectionDiffusion",
46 
49  : UnsteadySystem(pSession),
50  AdvectionSystem(pSession)
51  {
52  m_planeNumber = 0;
53  }
54 
55  /**
56  * @brief Initialisation object for the unsteady linear advection
57  * diffusion equation.
58  */
60  {
62 
63  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
64  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
65 
66  // Define Velocity fields
67  m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
68  std::vector<std::string> vel;
69  vel.push_back("Vx");
70  vel.push_back("Vy");
71  vel.push_back("Vz");
72  vel.resize(m_spacedim);
73 
74  EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
75 
76  m_session->MatchSolverInfo(
77  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
78 
80  {
81  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
82  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
83  }
84 
85  // Type of advection and diffusion classes to be used
86  switch(m_projectionType)
87  {
88  // Discontinuous field
90  {
91  // Do not forwards transform initial condition
92  m_homoInitialFwd = false;
93 
94  // Advection term
95  string advName;
96  string riemName;
97  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
99  CreateInstance(advName, advName);
100  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
101  GetFluxVectorAdv, this);
102  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
104  CreateInstance(riemName);
105  m_riemannSolver->SetScalar("Vn", &UnsteadyAdvectionDiffusion::
106  GetNormalVelocity, this);
107  m_advObject->SetRiemannSolver(m_riemannSolver);
108  m_advObject->InitObject (m_session, m_fields);
109 
110  // Diffusion term
111  std::string diffName;
112  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
114  CreateInstance(diffName, diffName);
115  m_diffusion->SetFluxVector(&UnsteadyAdvectionDiffusion::
116  GetFluxVectorDiff, this);
117  m_diffusion->InitObject(m_session, m_fields);
118  break;
119  }
120  // Continuous field
123  {
124  // Advection term
125  std::string advName;
126  m_session->LoadSolverInfo("AdvectionType", advName,
127  "NonConservative");
129  CreateInstance(advName, advName);
130  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
131  GetFluxVectorAdv, this);
132 
133  // In case of Galerkin explicit diffusion gives an error
135  {
136  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
137  }
138  // In case of Galerkin implicit diffusion: do nothing
139  break;
140  }
141  default:
142  {
143  ASSERTL0(false, "Unsupported projection type.");
144  break;
145  }
146  }
147 
150 
152  m_explicitDiffusion == 1)
153  {
155  }
156  }
157 
158  /**
159  * @brief Unsteady linear advection diffusion equation destructor.
160  */
162  {
163  }
164 
165  /**
166  * @brief Get the normal velocity for the unsteady linear advection
167  * diffusion equation.
168  */
170  {
171  // Number of trace (interface) points
172  int i;
173  int nTracePts = GetTraceNpoints();
174 
175  // Auxiliary variable to compute the normal velocity
176  Array<OneD, NekDouble> tmp(nTracePts);
177  m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
178 
179  // Reset the normal velocity
180  Vmath::Zero(nTracePts, m_traceVn, 1);
181 
182  for (i = 0; i < m_velocity.num_elements(); ++i)
183  {
184  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
185 
186  Vmath::Vvtvp(nTracePts,
187  m_traceNormals[i], 1,
188  tmp, 1,
189  m_traceVn, 1,
190  m_traceVn, 1);
191  }
192 
193  return m_traceVn;
194  }
195 
196  /**
197  * @brief Compute the right-hand side for the unsteady linear advection
198  * diffusion problem.
199  *
200  * @param inarray Given fields.
201  * @param outarray Calculated solution.
202  * @param time Time.
203  */
205  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
206  Array<OneD, Array<OneD, NekDouble> >&outarray,
207  const NekDouble time)
208  {
209  // Number of fields (variables of the problem)
210  int nVariables = inarray.num_elements();
211 
212  // Number of solution points
213  int nSolutionPts = GetNpoints();
214 
215  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
216 
217  for (int i = 0; i < nVariables; ++i)
218  {
219  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
220  }
221 
222  // RHS computation using the new advection base class
223  m_advObject->Advect(nVariables, m_fields, m_velocity,
224  inarray, outarray, time);
225 
226  // Negate the RHS
227  for (int i = 0; i < nVariables; ++i)
228  {
229  Vmath::Neg(nSolutionPts, outarray[i], 1);
230  }
231 
232  // No explicit diffusion for CG
234  {
235  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
236 
237  for (int i = 0; i < nVariables; ++i)
238  {
239  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
240  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
241  }
242  }
243 
244  }
245 
246  /**
247  * @brief Compute the projection for the unsteady advection
248  * diffusion problem.
249  *
250  * @param inarray Given fields.
251  * @param outarray Calculated solution.
252  * @param time Time.
253  */
255  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
256  Array<OneD, Array<OneD, NekDouble> > &outarray,
257  const NekDouble time)
258  {
259  int i;
260  int nvariables = inarray.num_elements();
261  SetBoundaryConditions(time);
262  switch(m_projectionType)
263  {
265  {
266  // Just copy over array
267  int npoints = GetNpoints();
268 
269  for(i = 0; i < nvariables; ++i)
270  {
271  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
272  }
273  break;
274  }
277  {
278  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
279 
280  for(i = 0; i < nvariables; ++i)
281  {
282  m_fields[i]->FwdTrans(inarray[i], coeffs);
283  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
284  }
285  break;
286  }
287  default:
288  {
289  ASSERTL0(false, "Unknown projection scheme");
290  break;
291  }
292  }
293  }
294 
295 
296  /* @brief Compute the diffusion term implicitly.
297  *
298  * @param inarray Given fields.
299  * @param outarray Calculated solution.
300  * @param time Time.
301  * @param lambda Diffusion coefficient.
302  */
304  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
305  Array<OneD, Array<OneD, NekDouble> >&outarray,
306  const NekDouble time,
307  const NekDouble lambda)
308  {
309  int nvariables = inarray.num_elements();
310  int nq = m_fields[0]->GetNpoints();
311 
313  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
314 
315  if(m_useSpecVanVisc)
316  {
319  }
320 
321  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
322  F[0] = Array<OneD, NekDouble> (nq*nvariables);
323 
324  for (int n = 1; n < nvariables; ++n)
325  {
326  F[n] = F[n-1] + nq;
327  }
328 
329  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
330  // inarray = input: \hat{rhs} -> output: \hat{Y}
331  // outarray = output: nabla^2 \hat{Y}
332  // where \hat = modal coeffs
333  for (int i = 0; i < nvariables; ++i)
334  {
335  // Multiply 1.0/timestep/lambda
336  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda],
337  inarray[i], 1, F[i], 1);
338  }
339 
340  //Setting boundary conditions
341  SetBoundaryConditions(time);
342 
343  for (int i = 0; i < nvariables; ++i)
344  {
345  // Solve a system of equations with Helmholtz solver
346  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
347  NullFlagList, factors);
348 
349  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
350  }
351  }
352 
353  /**
354  * @brief Return the flux vector for the advection part.
355  *
356  * @param physfield Fields.
357  * @param flux Resulting flux.
358  */
360  const Array<OneD, Array<OneD, NekDouble> > &physfield,
361  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &flux)
362  {
363  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
364  "Dimension of flux array and velocity array do not match");
365 
366  const int nq = m_fields[0]->GetNpoints();
367 
368  for (int i = 0; i < flux.num_elements(); ++i)
369  {
370  for (int j = 0; j < flux[0].num_elements(); ++j)
371  {
372  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
373  flux[i][j], 1);
374  }
375  }
376  }
377 
378  /**
379  * @brief Return the flux vector for the diffusion part.
380  *
381  * @param i Equation number.
382  * @param j Spatial direction.
383  * @param physfield Fields.
384  * @param derivatives First order derivatives.
385  * @param flux Resulting flux.
386  */
388  const int i,
389  const int j,
390  const Array<OneD, Array<OneD, NekDouble> > &physfield,
391  Array<OneD, Array<OneD, NekDouble> > &derivatives,
392  Array<OneD, Array<OneD, NekDouble> > &flux)
393  {
394  for (int k = 0; k < flux.num_elements(); ++k)
395  {
396  Vmath::Zero(GetNpoints(), flux[k], 1);
397  }
398  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
399  }
400 
403  {
405  }
406 }