Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
46 std::string DiffusionLDG::type =
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,
93  pFwd, 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, pBwd);
146 
147  Array<OneD, Array<OneD, NekDouble> > qdbase{nDim};
148 
149  for (std::size_t i = 0; i < nConvectiveFields; ++i)
150  {
151  for (std::size_t j = 0; j < nDim; ++j)
152  {
153  qdbase[j] = viscTensor[j][i];
154  }
155  fields[i]->IProductWRTDerivBase(qdbase, tmp);
156 
157  Vmath::Neg (nCoeffs, tmp, 1);
158  fields[i]->AddTraceIntegral (traceflux[i], tmp);
159  fields[i]->SetPhysState (false);
160  fields[i]->MultiplyByElmtInvMass(tmp, outarray[i]);
161  }
162 }
163 
166  const Array<OneD, Array<OneD, NekDouble> > &inarray,
168  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
169  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
170 {
171  std::size_t nConvectiveFields = fields.size();
172  std::size_t nDim = fields[0]->GetCoordim(0);
173  std::size_t nCoeffs = fields[0]->GetNcoeffs();
174  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
175 
176  Array<OneD, NekDouble> tmp{nCoeffs};
177  TensorOfArray3D<NekDouble> flux {nDim};
178  for (std::size_t j = 0; j < nDim; ++j)
179  {
180  flux[j] = Array<OneD, Array<OneD, NekDouble> > {nConvectiveFields};
181  for (std::size_t i = 0; i < nConvectiveFields; ++i)
182  {
183  flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
184  }
185  }
186 
187  NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
188 
189  for (std::size_t j = 0; j < nDim; ++j)
190  {
191  for (std::size_t i = 0; i < nConvectiveFields; ++i)
192  {
193  fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
194  Vmath::Neg (nCoeffs, tmp, 1);
195  fields[i]->AddTraceIntegral (flux[j][i], tmp);
196  fields[i]->SetPhysState (false);
197  fields[i]->MultiplyByElmtInvMass(tmp, tmp);
198  fields[i]->BwdTrans (tmp, qfield[j][i]);
199  }
200  }
201 }
202 
205  const Array<OneD, Array<OneD, NekDouble>> &inarray,
207  TensorOfArray3D<NekDouble> &viscTensor,
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  TensorOfArray3D<NekDouble> &viscTensor,
218  Array<OneD, Array<OneD, NekDouble> > &TraceFlux,
219  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
220  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
221  Array< OneD, int > &nonZeroIndex)
222 {
223  boost::ignore_unused(qfield, pFwd, pBwd, nonZeroIndex);
224  NumFluxforVector(fields, inarray, viscTensor, TraceFlux);
225 }
226 
229  const Array<OneD, Array<OneD, NekDouble>> &ufield,
231  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
232  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
233 {
234  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
235  std::size_t nvariables = fields.size();
236  std::size_t nDim = fields[0]->GetCoordim(0);
237 
238  Array<OneD, NekDouble> Fwd{nTracePts};
239  Array<OneD, NekDouble> Bwd{nTracePts};
240  Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
241 
242  // Get the sign of (v \cdot n), v = an arbitrary vector
243  // Evaluate upwind flux:
244  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
245  for (std::size_t i = 0; i < nvariables; ++i)
246  {
247  // Compute Fwd and Bwd value of ufield of i direction
248  if (pFwd == NullNekDoubleArrayOfArray ||
250  {
251  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
252  }
253  else
254  {
255  Fwd = pFwd[i];
256  Bwd = pBwd[i];
257  }
258 
259  // Upwind
260  Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
261 
262  // Imposing weak boundary condition with flux
263  if (fields[0]->GetBndCondExpansions().size())
264  {
265  ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
266  }
267 
268  for (std::size_t j = 0; j < nDim; ++j)
269  {
270  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
271  uflux[j][i], 1);
272  }
273  }
274 }
275 
278  const std::size_t var, const Array<OneD, const NekDouble> &ufield,
279  const Array<OneD, const NekDouble> &Fwd,
280  const Array<OneD, const NekDouble> &Bwd,
281  Array<OneD, NekDouble> &penaltyflux)
282 {
283  boost::ignore_unused(ufield, Bwd);
284  // Number of boundary regions
285  std::size_t nBndRegions =
286  fields[var]->GetBndCondExpansions().size();
287  std::size_t cnt = 0;
288  for (std::size_t i = 0; i < nBndRegions; ++i)
289  {
290  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
292  {
293  continue;
294  }
295 
296  // Number of boundary expansion related to that region
297  std::size_t nBndEdges =
298  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
299 
300  // Weakly impose boundary conditions by modifying flux values
301  for (std::size_t e = 0; e < nBndEdges; ++e)
302  {
303  std::size_t nBndEdgePts = fields[var]
304  ->GetBndCondExpansions()[i]
305  ->GetExp(e)
306  ->GetTotPoints();
307 
308  std::size_t id1 =
309  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
310 
311  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
312  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
313  cnt++));
314 
315  // AV boundary conditions
316  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
317  GetUserDefined(),"Wall") ||
318  boost::iequals(fields[var]->GetBndConditions()[i]->
319  GetUserDefined(),"Symmetry") ||
320  boost::iequals(fields[var]->GetBndConditions()[i]->
321  GetUserDefined(),"WallViscous") ||
322  boost::iequals(fields[var]->GetBndConditions()[i]->
323  GetUserDefined(),"WallAdiabatic"))
324  {
325  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
326  }
327  // For Dirichlet boundary condition: uflux = g_D
328  else if (fields[var]
329  ->GetBndConditions()[i]
330  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
331  {
332  Vmath::Vcopy(
333  nBndEdgePts,
334  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
335  1, &penaltyflux[id2], 1);
336  }
337  // For Neumann boundary condition: uflux = u+
338  else if ((fields[var]->GetBndConditions()[i])
339  ->GetBoundaryConditionType() ==
341  {
342  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
343  }
344  }
345  }
346 }
347 
348 /**
349  * @brief Build the numerical flux for the 2nd order derivatives
350  * todo: add variable coeff and h dependence to penalty term
351  */
354  const Array<OneD, Array<OneD, NekDouble>> &ufield,
357 {
358  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
359  std::size_t nvariables = fields.size();
360  std::size_t nDim = qfield.size();
361 
362  Array<OneD, NekDouble> Fwd{nTracePts};
363  Array<OneD, NekDouble> Bwd{nTracePts};
364  Array<OneD, NekDouble> qFwd{nTracePts};
365  Array<OneD, NekDouble> qBwd{nTracePts};
366  Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
367  Array<OneD, NekDouble> uterm{nTracePts};
368 
369  // Evaulate upwind flux:
370  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
371  for (std::size_t i = 0; i < nvariables; ++i)
372  {
373  // Generate Stability term = - C11 ( u- - u+ )
374  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
375  Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
376  Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
377 
378  qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
379  for (std::size_t j = 0; j < nDim; ++j)
380  {
381  // Compute Fwd and Bwd value of ufield of jth direction
382  fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
383 
384  // Downwind
385  Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
386 
387  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
388  qfluxtemp, 1);
389 
390  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
391  Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
392 
393  // Imposing weak boundary condition with flux
394  if (fields[0]->GetBndCondExpansions().size())
395  {
396  ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
397  qfluxtemp);
398  }
399 
400  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
401  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
402  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
403  Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
404  }
405  }
406 }
407 
408 /**
409  * Diffusion: Imposing weak boundary condition for q with flux
410  * uflux = g_D on Dirichlet boundary condition
411  * uflux = u_Fwd on Neumann boundary condition
412  */
415  const std::size_t var, const std::size_t dir,
416  const Array<OneD, const NekDouble> &qfield,
417  const Array<OneD, const NekDouble> &qFwd,
418  const Array<OneD, const NekDouble> &qBwd,
419  Array<OneD, NekDouble> &penaltyflux)
420 {
421  boost::ignore_unused(qfield, qBwd);
422 
423  std::size_t nBndRegions =
424  fields[var]->GetBndCondExpansions().size();
425  std::size_t cnt = 0;
426 
427  for (std::size_t i = 0; i < nBndRegions; ++i)
428  {
429  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
431  {
432  continue;
433  }
434  std::size_t nBndEdges =
435  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
436 
437  // Weakly impose boundary conditions by modifying flux values
438  for (std::size_t e = 0; e < nBndEdges; ++e)
439  {
440  std::size_t nBndEdgePts = fields[var]
441  ->GetBndCondExpansions()[i]
442  ->GetExp(e)
443  ->GetTotPoints();
444 
445  std::size_t id1 =
446  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
447 
448  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
449  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
450  cnt++));
451 
452  // AV boundary conditions
453  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
454  GetUserDefined(),"Wall") ||
455  boost::iequals(fields[var]->GetBndConditions()[i]->
456  GetUserDefined(),"Symmetry") ||
457  boost::iequals(fields[var]->GetBndConditions()[i]->
458  GetUserDefined(),"WallViscous") ||
459  boost::iequals(fields[var]->GetBndConditions()[i]->
460  GetUserDefined(),"WallAdiabatic"))
461  {
462  Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
463  }
464  // For Dirichlet boundary condition:
465  // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
466  else if (fields[var]
467  ->GetBndConditions()[i]
468  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
469  {
470  Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
471  &qFwd[id2], 1, &penaltyflux[id2], 1);
472  }
473  // For Neumann boundary condition: qflux = g_N
474  else if ((fields[var]->GetBndConditions()[i])
475  ->GetBoundaryConditionType() ==
477  {
478  Vmath::Vmul(
479  nBndEdgePts, &m_traceNormals[dir][id2], 1,
480  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
481  1, &penaltyflux[id2], 1);
482  }
483  }
484  }
485 }
486 
487 } // namespace SolverUtils
488 } // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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:249
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:236
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:224
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:466
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)
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)
Diffusion Volume Flux.
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
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)
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)
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:66
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)
Diffusion term Trace Flux.
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...
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63
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)
Diffusion Flux, calculate the physical derivatives.
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)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
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:1
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:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:322
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:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
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:372