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
42namespace Nektar
43{
44namespace 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 {
72 }
73 pFields[0]->GetTrace()->GetNormals(m_traceNormals);
74}
75
77 const std::size_t nConvectiveFields,
79 const Array<OneD, Array<OneD, NekDouble>> &inarray,
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,
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};
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}
213
216 const Array<OneD, Array<OneD, NekDouble>> &inarray,
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,
281 Array<OneD, NekDouble> &penaltyflux)
282{
283 boost::ignore_unused(ufield, Bwd);
284 // Number of boundary regions
285 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
286 std::size_t cnt = 0;
287 for (std::size_t i = 0; i < nBndRegions; ++i)
288 {
289 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
291 {
292 continue;
293 }
294
295 // Number of boundary expansion related to that region
296 std::size_t nBndEdges =
297 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
298
299 // Weakly impose boundary conditions by modifying flux values
300 for (std::size_t e = 0; e < nBndEdges; ++e)
301 {
302 std::size_t nBndEdgePts = fields[var]
303 ->GetBndCondExpansions()[i]
304 ->GetExp(e)
305 ->GetTotPoints();
306
307 std::size_t id1 =
308 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
309
310 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
311 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
312
313 // AV boundary conditions
314 if (boost::iequals(
315 fields[var]->GetBndConditions()[i]->GetUserDefined(),
316 "Wall") ||
317 boost::iequals(
318 fields[var]->GetBndConditions()[i]->GetUserDefined(),
319 "Symmetry") ||
320 boost::iequals(
321 fields[var]->GetBndConditions()[i]->GetUserDefined(),
322 "WallViscous") ||
323 boost::iequals(
324 fields[var]->GetBndConditions()[i]->GetUserDefined(),
325 "WallAdiabatic"))
326 {
327 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
328 }
329 // For Dirichlet boundary condition: uflux = g_D
330 else if (fields[var]
331 ->GetBndConditions()[i]
332 ->GetBoundaryConditionType() ==
334 {
336 nBndEdgePts,
337 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
338 1, &penaltyflux[id2], 1);
339 }
340 // For Neumann boundary condition: uflux = u+
341 else if ((fields[var]->GetBndConditions()[i])
342 ->GetBoundaryConditionType() ==
344 {
345 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
346 }
347 }
348 }
349}
350
351/**
352 * @brief Build the numerical flux for the 2nd order derivatives
353 * todo: add variable coeff and h dependence to penalty term
354 */
357 const Array<OneD, Array<OneD, NekDouble>> &ufield,
360{
361 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
362 std::size_t nvariables = fields.size();
363 std::size_t nDim = qfield.size();
364
365 Array<OneD, NekDouble> Fwd{nTracePts};
366 Array<OneD, NekDouble> Bwd{nTracePts};
367 Array<OneD, NekDouble> qFwd{nTracePts};
368 Array<OneD, NekDouble> qBwd{nTracePts};
369 Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
370 Array<OneD, NekDouble> uterm{nTracePts};
371
372 // Evaulate upwind flux:
373 // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
374 for (std::size_t i = 0; i < nvariables; ++i)
375 {
376 // Generate Stability term = - C11 ( u- - u+ )
377 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
378 Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
379 Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
380
381 qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
382 for (std::size_t j = 0; j < nDim; ++j)
383 {
384 // Compute Fwd and Bwd value of ufield of jth direction
385 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
386
387 // Downwind
388 Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
389
390 Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
391 qfluxtemp, 1);
392
393 // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
394 Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
395
396 // Imposing weak boundary condition with flux
397 if (fields[0]->GetBndCondExpansions().size())
398 {
399 ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
400 qfluxtemp);
401 }
402
403 // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
404 // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
405 // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
406 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
407 }
408 }
409}
410
411/**
412 * Diffusion: Imposing weak boundary condition for q with flux
413 * uflux = g_D on Dirichlet boundary condition
414 * uflux = u_Fwd on Neumann boundary condition
415 */
418 const std::size_t var, const std::size_t dir,
419 const Array<OneD, const NekDouble> &qfield,
422 Array<OneD, NekDouble> &penaltyflux)
423{
424 boost::ignore_unused(qfield, qBwd);
425
426 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
427 std::size_t cnt = 0;
428
429 for (std::size_t i = 0; i < nBndRegions; ++i)
430 {
431 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
433 {
434 continue;
435 }
436 std::size_t nBndEdges =
437 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
438
439 // Weakly impose boundary conditions by modifying flux values
440 for (std::size_t e = 0; e < nBndEdges; ++e)
441 {
442 std::size_t nBndEdgePts = fields[var]
443 ->GetBndCondExpansions()[i]
444 ->GetExp(e)
445 ->GetTotPoints();
446
447 std::size_t id1 =
448 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
449
450 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
451 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
452
453 // AV boundary conditions
454 if (boost::iequals(
455 fields[var]->GetBndConditions()[i]->GetUserDefined(),
456 "Wall") ||
457 boost::iequals(
458 fields[var]->GetBndConditions()[i]->GetUserDefined(),
459 "Symmetry") ||
460 boost::iequals(
461 fields[var]->GetBndConditions()[i]->GetUserDefined(),
462 "WallViscous") ||
463 boost::iequals(
464 fields[var]->GetBndConditions()[i]->GetUserDefined(),
465 "WallAdiabatic"))
466 {
467 Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
468 }
469 // For Dirichlet boundary condition:
470 // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
471 else if (fields[var]
472 ->GetBndConditions()[i]
473 ->GetBoundaryConditionType() ==
475 {
476 Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
477 &qFwd[id2], 1, &penaltyflux[id2], 1);
478 }
479 // For Neumann boundary condition: qflux = g_N
480 else if ((fields[var]->GetBndConditions()[i])
481 ->GetBoundaryConditionType() ==
483 {
485 nBndEdgePts, &m_traceNormals[dir][id2], 1,
486 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
487 1, &penaltyflux[id2], 1);
488 }
489 }
490 }
491}
492
493} // namespace SolverUtils
494} // 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 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.h:205
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:218
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:229
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:350
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
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)
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.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:110
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)
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 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_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
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:111
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) 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...
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:108
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.
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:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:354
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:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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:414