Nektar++
DiffusionLDG.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: DiffusionLDG.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: LDG diffusion class.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iomanip>
36 #include <iostream>
37 
38 #include <boost/algorithm/string/predicate.hpp>
39 
41 
42 namespace Nektar
43 {
44 namespace SolverUtils
45 {
47  "LDG", DiffusionLDG::create, "Local Discontinuous Galkerin");
48 
50 {
51 }
52 
56 {
57  m_session = pSession;
58 
59  m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
60 
61  // Set up penalty term for LDG
62  m_session->LoadParameter("LDGc11", m_C11, 1.0);
63 
64  // Setting up the normals
65  std::size_t nDim = pFields[0]->GetCoordim(0);
66  std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
67 
69  for (std::size_t i = 0; i < nDim; ++i)
70  {
71  m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts};
72  }
73  pFields[0]->GetTrace()->GetNormals(m_traceNormals);
74 }
75 
77  const std::size_t nConvectiveFields,
79  const Array<OneD, Array<OneD, NekDouble>> &inarray,
80  Array<OneD, Array<OneD, NekDouble>> &outarray,
81  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
82  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
83 {
84  std::size_t nCoeffs = fields[0]->GetNcoeffs();
85 
86  Array<OneD, Array<OneD, NekDouble>> tmp{nConvectiveFields};
87  for (std::size_t i = 0; i < nConvectiveFields; ++i)
88  {
89  tmp[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
90  }
91 
92  DiffusionLDG::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, tmp, pFwd,
93  pBwd);
94 
95  for (std::size_t i = 0; i < nConvectiveFields; ++i)
96  {
97  fields[i]->BwdTrans(tmp[i], outarray[i]);
98  }
99 }
100 
102  const std::size_t nConvectiveFields,
104  const Array<OneD, Array<OneD, NekDouble>> &inarray,
105  Array<OneD, Array<OneD, NekDouble>> &outarray,
106  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
107  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
108 {
109  std::size_t nDim = fields[0]->GetCoordim(0);
110  std::size_t nPts = fields[0]->GetTotPoints();
111  std::size_t nCoeffs = fields[0]->GetNcoeffs();
112  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
113 
114  Array<OneD, NekDouble> tmp{nCoeffs};
115 
116  TensorOfArray3D<NekDouble> qfield{nDim};
117  for (std::size_t j = 0; j < nDim; ++j)
118  {
119  qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
120  for (std::size_t i = 0; i < nConvectiveFields; ++i)
121  {
122  qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
123  }
124  }
125 
126  Array<OneD, Array<OneD, NekDouble>> traceflux{nConvectiveFields};
127  for (std::size_t i = 0; i < nConvectiveFields; ++i)
128  {
129  traceflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
130  }
131 
132  DiffuseCalcDerivative(fields, inarray, qfield, pFwd, pBwd);
133 
134  // Initialize viscous tensor
136  for (std::size_t j = 0; j < nDim; ++j)
137  {
138  viscTensor[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
139  for (std::size_t i = 0; i < nConvectiveFields; ++i)
140  {
141  viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
142  }
143  }
144  DiffuseVolumeFlux(fields, inarray, qfield, viscTensor);
145  DiffuseTraceFlux(fields, inarray, qfield, viscTensor, traceflux, pFwd,
146  pBwd);
147 
149 
150  for (std::size_t i = 0; i < nConvectiveFields; ++i)
151  {
152  for (std::size_t j = 0; j < nDim; ++j)
153  {
154  qdbase[j] = viscTensor[j][i];
155  }
156  fields[i]->IProductWRTDerivBase(qdbase, tmp);
157 
158  Vmath::Neg(nCoeffs, tmp, 1);
159  fields[i]->AddTraceIntegral(traceflux[i], tmp);
160  fields[i]->SetPhysState(false);
161  fields[i]->MultiplyByElmtInvMass(tmp, outarray[i]);
162  }
163 }
164 
167  const Array<OneD, Array<OneD, NekDouble>> &inarray,
169  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
170  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
171 {
172  std::size_t nConvectiveFields = fields.size();
173  std::size_t nDim = fields[0]->GetCoordim(0);
174  std::size_t nCoeffs = fields[0]->GetNcoeffs();
175  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
176 
177  Array<OneD, NekDouble> tmp{nCoeffs};
178  TensorOfArray3D<NekDouble> flux{nDim};
179  for (std::size_t j = 0; j < nDim; ++j)
180  {
181  flux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
182  for (std::size_t i = 0; i < nConvectiveFields; ++i)
183  {
184  flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
185  }
186  }
187 
188  NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
189 
190  for (std::size_t j = 0; j < nDim; ++j)
191  {
192  for (std::size_t i = 0; i < nConvectiveFields; ++i)
193  {
194  fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
195  Vmath::Neg(nCoeffs, tmp, 1);
196  fields[i]->AddTraceIntegral(flux[j][i], tmp);
197  fields[i]->SetPhysState(false);
198  fields[i]->MultiplyByElmtInvMass(tmp, tmp);
199  fields[i]->BwdTrans(tmp, qfield[j][i]);
200  }
201  }
202 }
203 
206  const Array<OneD, Array<OneD, NekDouble>> &inarray,
208  Array<OneD, int> &nonZeroIndex)
209 {
210  boost::ignore_unused(fields, nonZeroIndex);
211  m_fluxVector(inarray, qfield, viscTensor);
212 }
215  const Array<OneD, Array<OneD, NekDouble>> &inarray,
217  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
218  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
219  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
220  Array<OneD, int> &nonZeroIndex)
221 {
222  boost::ignore_unused(qfield, pFwd, pBwd, nonZeroIndex);
223  NumFluxforVector(fields, inarray, viscTensor, TraceFlux);
224 }
225 
228  const Array<OneD, Array<OneD, NekDouble>> &ufield,
230  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
231  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
232 {
233  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
234  std::size_t nvariables = fields.size();
235  std::size_t nDim = fields[0]->GetCoordim(0);
236 
237  Array<OneD, NekDouble> Fwd{nTracePts};
238  Array<OneD, NekDouble> Bwd{nTracePts};
239  Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
240 
241  // Get the sign of (v \cdot n), v = an arbitrary vector
242  // Evaluate upwind flux:
243  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
244  for (std::size_t i = 0; i < nvariables; ++i)
245  {
246  // Compute Fwd and Bwd value of ufield of i direction
247  if (pFwd == NullNekDoubleArrayOfArray ||
249  {
250  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
251  }
252  else
253  {
254  Fwd = pFwd[i];
255  Bwd = pBwd[i];
256  }
257 
258  // Upwind
259  Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
260 
261  // Imposing weak boundary condition with flux
262  if (fields[0]->GetBndCondExpansions().size())
263  {
264  ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
265  }
266 
267  for (std::size_t j = 0; j < nDim; ++j)
268  {
269  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
270  uflux[j][i], 1);
271  }
272  }
273 }
274 
277  const std::size_t var, const Array<OneD, const NekDouble> &ufield,
278  const Array<OneD, const NekDouble> &Fwd,
279  const Array<OneD, const NekDouble> &Bwd,
280  Array<OneD, NekDouble> &penaltyflux)
281 {
282  boost::ignore_unused(ufield, Bwd);
283  // Number of boundary regions
284  std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
285  std::size_t cnt = 0;
286  for (std::size_t i = 0; i < nBndRegions; ++i)
287  {
288  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
290  {
291  continue;
292  }
293 
294  // Number of boundary expansion related to that region
295  std::size_t nBndEdges =
296  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
297 
298  // Weakly impose boundary conditions by modifying flux values
299  for (std::size_t e = 0; e < nBndEdges; ++e)
300  {
301  std::size_t nBndEdgePts = fields[var]
302  ->GetBndCondExpansions()[i]
303  ->GetExp(e)
304  ->GetTotPoints();
305 
306  std::size_t id1 =
307  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
308 
309  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
310  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
311 
312  // AV boundary conditions
313  if (boost::iequals(
314  fields[var]->GetBndConditions()[i]->GetUserDefined(),
315  "Wall") ||
316  boost::iequals(
317  fields[var]->GetBndConditions()[i]->GetUserDefined(),
318  "Symmetry") ||
319  boost::iequals(
320  fields[var]->GetBndConditions()[i]->GetUserDefined(),
321  "WallViscous") ||
322  boost::iequals(
323  fields[var]->GetBndConditions()[i]->GetUserDefined(),
324  "WallAdiabatic"))
325  {
326  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
327  }
328  // For Dirichlet boundary condition: uflux = g_D
329  else if (fields[var]
330  ->GetBndConditions()[i]
331  ->GetBoundaryConditionType() ==
333  {
334  Vmath::Vcopy(
335  nBndEdgePts,
336  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
337  1, &penaltyflux[id2], 1);
338  }
339  // For Neumann boundary condition: uflux = u+
340  else if ((fields[var]->GetBndConditions()[i])
341  ->GetBoundaryConditionType() ==
343  {
344  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
345  }
346  }
347  }
348 }
349 
350 /**
351  * @brief Build the numerical flux for the 2nd order derivatives
352  * todo: add variable coeff and h dependence to penalty term
353  */
356  const Array<OneD, Array<OneD, NekDouble>> &ufield,
359 {
360  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
361  std::size_t nvariables = fields.size();
362  std::size_t nDim = qfield.size();
363 
364  Array<OneD, NekDouble> Fwd{nTracePts};
365  Array<OneD, NekDouble> Bwd{nTracePts};
366  Array<OneD, NekDouble> qFwd{nTracePts};
367  Array<OneD, NekDouble> qBwd{nTracePts};
368  Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
369  Array<OneD, NekDouble> uterm{nTracePts};
370 
371  // Evaulate upwind flux:
372  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
373  for (std::size_t i = 0; i < nvariables; ++i)
374  {
375  // Generate Stability term = - C11 ( u- - u+ )
376  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
377  Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
378  Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
379 
380  qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
381  for (std::size_t j = 0; j < nDim; ++j)
382  {
383  // Compute Fwd and Bwd value of ufield of jth direction
384  fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
385 
386  // Downwind
387  Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
388 
389  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
390  qfluxtemp, 1);
391 
392  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
393  Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
394 
395  // Imposing weak boundary condition with flux
396  if (fields[0]->GetBndCondExpansions().size())
397  {
398  ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
399  qfluxtemp);
400  }
401 
402  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
403  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
404  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
405  Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
406  }
407  }
408 }
409 
410 /**
411  * Diffusion: Imposing weak boundary condition for q with flux
412  * uflux = g_D on Dirichlet boundary condition
413  * uflux = u_Fwd on Neumann boundary condition
414  */
417  const std::size_t var, const std::size_t dir,
418  const Array<OneD, const NekDouble> &qfield,
419  const Array<OneD, const NekDouble> &qFwd,
420  const Array<OneD, const NekDouble> &qBwd,
421  Array<OneD, NekDouble> &penaltyflux)
422 {
423  boost::ignore_unused(qfield, qBwd);
424 
425  std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
426  std::size_t cnt = 0;
427 
428  for (std::size_t i = 0; i < nBndRegions; ++i)
429  {
430  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
432  {
433  continue;
434  }
435  std::size_t nBndEdges =
436  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
437 
438  // Weakly impose boundary conditions by modifying flux values
439  for (std::size_t e = 0; e < nBndEdges; ++e)
440  {
441  std::size_t nBndEdgePts = fields[var]
442  ->GetBndCondExpansions()[i]
443  ->GetExp(e)
444  ->GetTotPoints();
445 
446  std::size_t id1 =
447  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
448 
449  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
450  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
451 
452  // AV boundary conditions
453  if (boost::iequals(
454  fields[var]->GetBndConditions()[i]->GetUserDefined(),
455  "Wall") ||
456  boost::iequals(
457  fields[var]->GetBndConditions()[i]->GetUserDefined(),
458  "Symmetry") ||
459  boost::iequals(
460  fields[var]->GetBndConditions()[i]->GetUserDefined(),
461  "WallViscous") ||
462  boost::iequals(
463  fields[var]->GetBndConditions()[i]->GetUserDefined(),
464  "WallAdiabatic"))
465  {
466  Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
467  }
468  // For Dirichlet boundary condition:
469  // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
470  else if (fields[var]
471  ->GetBndConditions()[i]
472  ->GetBoundaryConditionType() ==
474  {
475  Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
476  &qFwd[id2], 1, &penaltyflux[id2], 1);
477  }
478  // For Neumann boundary condition: qflux = g_N
479  else if ((fields[var]->GetBndConditions()[i])
480  ->GetBoundaryConditionType() ==
482  {
483  Vmath::Vmul(
484  nBndEdgePts, &m_traceNormals[dir][id2], 1,
485  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
486  1, &penaltyflux[id2], 1);
487  }
488  }
489  }
490 }
491 
492 } // namespace SolverUtils
493 } // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
Definition: Diffusion.h:206
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.cpp:194
SOLVER_UTILS_EXPORT void DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
Definition: Diffusion.h:217
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:402
virtual void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void ApplyScalarBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionLDG.h:49
virtual void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
void ApplyVectorBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:66
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
virtual void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63
virtual void v_Diffuse(const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
virtual void v_DiffuseCoeffs(const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble >> &qflux)
Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to p...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
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:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419