Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::DiffusionLDG Class Reference

#include <DiffusionLDG.h>

Inheritance diagram for Nektar::SolverUtils::DiffusionLDG:
[legend]

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionLDG ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
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=NullNekDoubleArrayofArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayofArray)
 
void NumFluxforScalar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayofArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayofArray)
 
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)
 
void NumFluxforVector (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, 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 penalty term. More...
 
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)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual Array< OneD, Array< OneD, Array< OneD, NekDouble > > > & v_GetFluxTensor ()
 

Protected Attributes

std::string m_shockCaptureType
 
NekDouble m_C11
 Coefficient of penalty term. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
LibUtilities::SessionReaderSharedPtr m_session
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, 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=NullNekDoubleArrayofArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayofArray)
 
SOLVER_UTILS_EXPORT void FluxVec (Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &fluxvector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxVector (DiffusionFluxVecCB fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVectorNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxVectorNS (DiffusionFluxVecCBNS fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxPenaltyNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxPenaltyNS (DiffusionFluxPenaltyNS flux)
 
void SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual Array< OneD, Array< OneD, Array< OneD, NekDouble > > > & GetFluxTensor ()
 

Detailed Description

Definition at line 46 of file DiffusionLDG.h.

Constructor & Destructor Documentation

◆ DiffusionLDG()

Nektar::SolverUtils::DiffusionLDG::DiffusionLDG ( )
protected

Definition at line 49 of file DiffusionLDG.cpp.

Referenced by create().

50 {
51 }

Member Function Documentation

◆ ApplyScalarBCs()

void Nektar::SolverUtils::DiffusionLDG::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 
)
protected

Definition at line 207 of file DiffusionLDG.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, and Vmath::Vcopy().

Referenced by NumFluxforScalar().

213 {
214  boost::ignore_unused(ufield, Bwd);
215  // Number of boundary regions
216  std::size_t nBndRegions =
217  fields[var]->GetBndCondExpansions().num_elements();
218  std::size_t cnt = 0;
219  for (std::size_t i = 0; i < nBndRegions; ++i)
220  {
221  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
223  {
224  continue;
225  }
226 
227  // Number of boundary expansion related to that region
228  std::size_t nBndEdges =
229  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
230 
231  // Weakly impose boundary conditions by modifying flux values
232  for (std::size_t e = 0; e < nBndEdges; ++e)
233  {
234  std::size_t nBndEdgePts = fields[var]
235  ->GetBndCondExpansions()[i]
236  ->GetExp(e)
237  ->GetTotPoints();
238 
239  std::size_t id1 =
240  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
241 
242  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
243  fields[0]->GetTraceMap()->GetBndCondTraceToGlobalTraceMap(
244  cnt++));
245 
246  // AV boundary conditions
247  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
248  GetUserDefined(),"Wall") ||
249  boost::iequals(fields[var]->GetBndConditions()[i]->
250  GetUserDefined(),"Symmetry") ||
251  boost::iequals(fields[var]->GetBndConditions()[i]->
252  GetUserDefined(),"WallViscous") ||
253  boost::iequals(fields[var]->GetBndConditions()[i]->
254  GetUserDefined(),"WallAdiabatic"))
255  {
256  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
257  }
258  // For Dirichlet boundary condition: uflux = g_D
259  else if (fields[var]
260  ->GetBndConditions()[i]
261  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
262  {
263  Vmath::Vcopy(
264  nBndEdgePts,
265  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
266  1, &penaltyflux[id2], 1);
267  }
268  // For Neumann boundary condition: uflux = u+
269  else if ((fields[var]->GetBndConditions()[i])
270  ->GetBoundaryConditionType() ==
272  {
273  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
274  }
275  }
276  }
277 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ ApplyVectorBCs()

void Nektar::SolverUtils::DiffusionLDG::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 
)
protected

Diffusion: Imposing weak boundary condition for q with flux uflux = g_D on Dirichlet boundary condition uflux = u_Fwd on Neumann boundary condition

Definition at line 344 of file DiffusionLDG.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, m_traceNormals, Vmath::Vmul(), and Vmath::Zero().

Referenced by NumFluxforVector().

351 {
352  boost::ignore_unused(qfield, qBwd);
353 
354  std::size_t nBndRegions =
355  fields[var]->GetBndCondExpansions().num_elements();
356  std::size_t cnt = 0;
357 
358  for (std::size_t i = 0; i < nBndRegions; ++i)
359  {
360  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
362  {
363  continue;
364  }
365  std::size_t nBndEdges =
366  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
367 
368  // Weakly impose boundary conditions by modifying flux values
369  for (std::size_t e = 0; e < nBndEdges; ++e)
370  {
371  std::size_t nBndEdgePts = fields[var]
372  ->GetBndCondExpansions()[i]
373  ->GetExp(e)
374  ->GetTotPoints();
375 
376  std::size_t id1 =
377  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
378 
379  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
380  fields[0]->GetTraceMap()->GetBndCondTraceToGlobalTraceMap(
381  cnt++));
382 
383  // AV boundary conditions
384  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
385  GetUserDefined(),"Wall") ||
386  boost::iequals(fields[var]->GetBndConditions()[i]->
387  GetUserDefined(),"Symmetry") ||
388  boost::iequals(fields[var]->GetBndConditions()[i]->
389  GetUserDefined(),"WallViscous") ||
390  boost::iequals(fields[var]->GetBndConditions()[i]->
391  GetUserDefined(),"WallAdiabatic"))
392  {
393  Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
394  }
395  // For Dirichlet boundary condition:
396  // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
397  else if (fields[var]
398  ->GetBndConditions()[i]
399  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
400  {
401  Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
402  &qFwd[id2], 1, &penaltyflux[id2], 1);
403  }
404  // For Neumann boundary condition: qflux = g_N
405  else if ((fields[var]->GetBndConditions()[i])
406  ->GetBoundaryConditionType() ==
408  {
409  Vmath::Vmul(
410  nBndEdgePts, &m_traceNormals[dir][id2], 1,
411  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
412  1, &penaltyflux[id2], 1);
413  }
414  }
415  }
416 }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
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:186

◆ create()

static DiffusionSharedPtr Nektar::SolverUtils::DiffusionLDG::create ( std::string  diffType)
inlinestatic

Definition at line 49 of file DiffusionLDG.h.

References DiffusionLDG().

50  {
51  boost::ignore_unused(diffType);
52  return DiffusionSharedPtr(new DiffusionLDG());
53  }
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:183

◆ NumFluxforScalar()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforScalar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  uflux,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd = NullNekDoubleArrayofArray,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd = NullNekDoubleArrayofArray 
)
protected

Definition at line 158 of file DiffusionLDG.cpp.

References ApplyScalarBCs(), m_traceNormals, Nektar::NullNekDoubleArrayofArray, Vmath::Vcopy(), and Vmath::Vmul().

Referenced by v_Diffuse().

164 {
165  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
166  std::size_t nvariables = fields.num_elements();
167  std::size_t nDim = fields[0]->GetCoordim(0);
168 
169  Array<OneD, NekDouble> Fwd{nTracePts};
170  Array<OneD, NekDouble> Bwd{nTracePts};
171  Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
172 
173  // Get the sign of (v \cdot n), v = an arbitrary vector
174  // Evaluate upwind flux:
175  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
176  for (std::size_t i = 0; i < nvariables; ++i)
177  {
178  // Compute Fwd and Bwd value of ufield of i direction
179  if (pFwd == NullNekDoubleArrayofArray ||
181  {
182  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
183  }
184  else
185  {
186  Fwd = pFwd[i];
187  Bwd = pBwd[i];
188  }
189 
190  // Upwind
191  Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
192 
193  // Imposing weak boundary condition with flux
194  if (fields[0]->GetBndCondExpansions().num_elements())
195  {
196  ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
197  }
198 
199  for (std::size_t j = 0; j < nDim; ++j)
200  {
201  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
202  uflux[j][i], 1);
203  }
204  }
205 }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
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)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayofArray
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:186

◆ NumFluxforVector()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforVector ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
protected

Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to penalty term.

Definition at line 283 of file DiffusionLDG.cpp.

References ApplyVectorBCs(), m_C11, m_traceNormals, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_Diffuse().

288 {
289  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
290  std::size_t nvariables = fields.num_elements();
291  std::size_t nDim = qfield.num_elements();
292 
293  Array<OneD, NekDouble> Fwd{nTracePts};
294  Array<OneD, NekDouble> Bwd{nTracePts};
295  Array<OneD, NekDouble> qFwd{nTracePts};
296  Array<OneD, NekDouble> qBwd{nTracePts};
297  Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
298  Array<OneD, NekDouble> uterm{nTracePts};
299 
300  // Evaulate upwind flux:
301  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
302  for (std::size_t i = 0; i < nvariables; ++i)
303  {
304  // Generate Stability term = - C11 ( u- - u+ )
305  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
306  Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
307  Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
308 
309  qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
310  for (std::size_t j = 0; j < nDim; ++j)
311  {
312  // Compute Fwd and Bwd value of ufield of jth direction
313  fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
314 
315  // Downwind
316  Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
317 
318  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
319  qfluxtemp, 1);
320 
321  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
322  Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
323 
324  // Imposing weak boundary condition with flux
325  if (fields[0]->GetBndCondExpansions().num_elements())
326  {
327  ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
328  qfluxtemp);
329  }
330 
331  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
332  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
333  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
334  Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
335  }
336  }
337 }
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 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:216
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:346
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:302
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:186

◆ v_Diffuse()

void Nektar::SolverUtils::DiffusionLDG::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 = NullNekDoubleArrayofArray,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd = NullNekDoubleArrayofArray 
)
protectedvirtual

Implements Nektar::SolverUtils::Diffusion.

Definition at line 76 of file DiffusionLDG.cpp.

References Nektar::SolverUtils::Diffusion::m_fluxVector, Vmath::Neg(), NumFluxforScalar(), and NumFluxforVector().

83 {
84  std::size_t nDim = fields[0]->GetCoordim(0);
85  std::size_t nPts = fields[0]->GetTotPoints();
86  std::size_t nCoeffs = fields[0]->GetNcoeffs();
87  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
88 
89  Array<OneD, NekDouble> tmp{nCoeffs};
90  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> flux{nDim};
91  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> qfield{nDim};
92 
93  for (std::size_t j = 0; j < nDim; ++j)
94  {
95  qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
96  flux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
97 
98  for (std::size_t i = 0; i < nConvectiveFields; ++i)
99  {
100  qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
101  flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
102  }
103  }
104 
105  // Compute q_{\eta} and q_{\xi}
106  // Obtain numerical fluxes
107 
108  NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
109 
110  for (std::size_t j = 0; j < nDim; ++j)
111  {
112  for (std::size_t i = 0; i < nConvectiveFields; ++i)
113  {
114  fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
115  Vmath::Neg(nCoeffs, tmp, 1);
116  fields[i]->AddTraceIntegral(flux[j][i], tmp);
117  fields[i]->SetPhysState(false);
118  fields[i]->MultiplyByElmtInvMass(tmp, tmp);
119  fields[i]->BwdTrans(tmp, qfield[j][i]);
120  }
121  }
122 
123  // Initialize viscous tensor
124  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> viscTensor{nDim};
125  for (std::size_t j = 0; j < nDim; ++j)
126  {
127  viscTensor[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
128  for (std::size_t i = 0; i < nConvectiveFields; ++i)
129  {
130  viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
131  }
132  }
133  // Get viscous tensor
134  m_fluxVector(inarray, qfield, viscTensor);
135 
136  // Compute u from q_{\eta} and q_{\xi}
137  // Obtain numerical fluxes
138  NumFluxforVector(fields, inarray, viscTensor, flux[0]);
139 
140  Array<OneD, Array<OneD, NekDouble>> qdbase{nDim};
141  for (std::size_t i = 0; i < nConvectiveFields; ++i)
142  {
143  for (std::size_t j = 0; j < nDim; ++j)
144  {
145  qdbase[j] = viscTensor[j][i];
146  }
147  fields[i]->IProductWRTDerivBase(qdbase, tmp);
148 
149  // Evaulate <\phi, \hat{F}\cdot n> - outarray[i]
150  Vmath::Neg(nCoeffs, tmp, 1);
151  fields[i]->AddTraceIntegral(flux[0][i], tmp);
152  fields[i]->SetPhysState(false);
153  fields[i]->MultiplyByElmtInvMass(tmp, tmp);
154  fields[i]->BwdTrans(tmp, outarray[i]);
155  }
156 }
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayofArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayofArray)
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, 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...
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:150

◆ v_InitObject()

void Nektar::SolverUtils::DiffusionLDG::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 53 of file DiffusionLDG.cpp.

References m_C11, m_session, m_shockCaptureType, and m_traceNormals.

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 
68  m_traceNormals = Array<OneD, Array<OneD, NekDouble>>{nDim};
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 }
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:66
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63

Member Data Documentation

◆ m_C11

NekDouble Nektar::SolverUtils::DiffusionLDG::m_C11
protected

Coefficient of penalty term.

Definition at line 63 of file DiffusionLDG.h.

Referenced by NumFluxforVector(), and v_InitObject().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionLDG::m_session
protected

Definition at line 66 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_shockCaptureType

std::string Nektar::SolverUtils::DiffusionLDG::m_shockCaptureType
protected

Definition at line 60 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_traceNormals

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDG::m_traceNormals
protected

Definition at line 65 of file DiffusionLDG.h.

Referenced by ApplyVectorBCs(), NumFluxforScalar(), NumFluxforVector(), and v_InitObject().

◆ type

std::string Nektar::SolverUtils::DiffusionLDG::type
static
Initial value:

Definition at line 55 of file DiffusionLDG.h.