Nektar++
Diffusion.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Diffusion.h
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: Abstract base class for diffusion.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERUTILS_DIFFUSION
36#define NEKTAR_SOLVERUTILS_DIFFUSION
37
38#include <functional>
39#include <string>
40
48
49namespace Nektar::SolverUtils
50{
51
52class Diffusion;
53
54/// A shared pointer to an EquationSystem object
55typedef std::shared_ptr<Diffusion> DiffusionSharedPtr;
56
57/// Datatype of the NekFactory used to instantiate classes derived
58/// from the Diffusion class.
62
63typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
67
68typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
72
73typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
77
78/**
79 * Parameter list meaning:
80 * 1st: field conservative variables
81 * 2th: Devrivatives of field conservative varialbes
82 * 3rd: the current time for time-dependent boundary
83 * 4th: Fwd of field conservative variables optional
84 * 5th: Fwd of Devrivatives(2nd) optional
85 *
86 * a null pointer need to be passed for optional parameters
87 */
88typedef std::function<void(
89 const Array<OneD, const Array<OneD, NekDouble>> &,
91 const Array<OneD, const Array<OneD, NekDouble>> &,
92 const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &)>
94
95/**
96 * Parameter list meaning:
97 * 1st: nvariables
98 * 2nd: nspaceDimension
99 * 3rd: field conservative variables
100 * 4th: Devrivatives of field conservative varialbes
101 * 5th: nonzero flux index array, optional
102 * 6th: normal vectors optional
103 *
104 * a null pointer need to be passed for optional parameters
105 */
106typedef std::function<void(
107 const int, const Array<OneD, Array<OneD, NekDouble>> &,
111
112/**
113 * Parameter list meaning:
114 * 1st: nvariables
115 * 2nd: nspaceDimension
116 * 3rd: trace conservative variables for Diffusion Flux Jacobian
117 * 4th: trace conservative variables( usually the jump of trace value)
118 * 5th: trace symmetric flux
119 * 6th: nonzero flux index array, optional
120 *
121 * a null pointer need to be passed for optional parameters
122 */
123typedef std::function<void(
124 const int, const Array<OneD, Array<OneD, NekDouble>> &,
128
129/**
130 * Parameter list meaning:
131 * 1rd: trace conservative variables
132 */
133typedef std::function<void(Array<OneD, Array<OneD, NekDouble>> &)>
135
137{
138public:
142 {
143 v_InitObject(pSession, pFields);
144
145 // Div curl storage
146 int nPts = pFields[0]->GetTotPoints();
147 m_divVel = Array<OneD, NekDouble>(nPts, 0.0);
150 }
151
153 const std::size_t nConvectiveFields,
155 const Array<OneD, Array<OneD, NekDouble>> &inarray,
157 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
159 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
161 {
162 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
163 }
164
166 const std::size_t nConvectiveFields,
168 const Array<OneD, Array<OneD, NekDouble>> &inarray,
169 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
170 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
172 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
174 {
175 m_time = time;
176 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
177 }
178
180 const std::size_t nConvectiveFields,
182 const Array<OneD, Array<OneD, NekDouble>> &inarray,
184 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
186 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
188 {
189 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
190 pBwd);
191 }
192
194 const std::size_t nConvectiveFields,
196 const Array<OneD, Array<OneD, NekDouble>> &inarray,
197 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
198 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
200 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
202 {
203 m_time = time;
204 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
205 pBwd);
206 }
207
208 // Diffusion Calculate the physical derivatives
211 const Array<OneD, Array<OneD, NekDouble>> &inarray,
213 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
215 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
217 {
218 v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
219 }
220
221 /// Diffusion Volume FLux
224 const Array<OneD, Array<OneD, NekDouble>> &inarray,
226 TensorOfArray3D<NekDouble> &VolumeFlux,
227 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
228 {
229 v_DiffuseVolumeFlux(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
230 }
231
232 /// Diffusion term Trace Flux
235 const Array<OneD, Array<OneD, NekDouble>> &inarray,
237 TensorOfArray3D<NekDouble> &VolumeFlux,
238 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
239 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
241 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
243 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
244 {
245 v_DiffuseTraceFlux(fields, inarray, qfields, VolumeFlux, TraceFlux,
246 pFwd, pBwd, nonZeroIndex);
247 }
248
251 {
252 v_SetHomoDerivs(deriv);
253 }
254
256 {
257 return v_GetFluxTensor();
258 }
259
260 /// Get trace normal
263 {
264 return v_GetTraceNormal();
265 }
266
267 template <typename FuncPointerT, typename ObjectPointerT>
268 void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
269 {
270 m_fluxVector = std::bind(func, obj, std::placeholders::_1,
271 std::placeholders::_2, std::placeholders::_3);
272 }
273
275 {
276 m_fluxVector = fluxVector;
277 }
278
279 template <typename FuncPointerT, typename ObjectPointerT>
280 void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
281 {
283 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
284 std::placeholders::_3);
285 }
286
288 {
289 m_fluxVectorNS = fluxVector;
290 }
291
292 template <typename FuncPointerT, typename ObjectPointerT>
293 void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
294 {
296 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
297 std::placeholders::_3);
298 }
299
301 {
302 m_fluxPenaltyNS = flux;
303 }
304
305 template <typename FuncPointerT, typename ObjectPointerT>
306 void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
307 {
309 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
310 std::placeholders::_3, std::placeholders::_4,
311 std::placeholders::_5, std::placeholders::_6);
312 }
313
315 {
317 }
318
319 template <typename FuncPointerT, typename ObjectPointerT>
320 void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
321 {
323 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
324 std::placeholders::_3, std::placeholders::_4,
325 std::placeholders::_5, std::placeholders::_6);
326 }
327
329 {
331 }
332
333 template <typename FuncPointerT, typename ObjectPointerT>
334 void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
335 {
336 m_SpecialBndTreat = std::bind(func, obj, std::placeholders::_1);
337 }
338
339 template <typename FuncPointerT, typename ObjectPointerT>
340 void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
341 {
343 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
344 std::placeholders::_3, std::placeholders::_4,
345 std::placeholders::_5, std::placeholders::_6);
346 }
347
349 Array<OneD, Array<OneD, NekDouble>> &gridVelocityTrace)
350 {
351 m_gridVelocityTrace = gridVelocityTrace;
352 }
353
354protected:
355 /// Params for Ducros sensor
359
369
370 SOLVER_UTILS_EXPORT virtual ~Diffusion() = default;
371
375
377 const std::size_t nConvectiveFields,
379 const Array<OneD, Array<OneD, NekDouble>> &inarray,
381 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
382 const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
383
385 const std::size_t nConvectiveFields,
387 const Array<OneD, Array<OneD, NekDouble>> &inarray,
389 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
390 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
391
392 /// Diffusion Flux, calculate the physical derivatives
395 const Array<OneD, Array<OneD, NekDouble>> &inarray,
397 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
398 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
399
400 /// Diffusion Volume Flux
403 const Array<OneD, Array<OneD, NekDouble>> &inarray,
405 TensorOfArray3D<NekDouble> &VolumeFlux, Array<OneD, int> &nonZeroIndex);
406
407 /// Diffusion term Trace Flux
410 const Array<OneD, Array<OneD, NekDouble>> &inarray,
412 TensorOfArray3D<NekDouble> &VolumeFlux,
413 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
414 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
415 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
416 Array<OneD, int> &nonZeroIndex);
417
418 virtual void v_SetHomoDerivs(
419 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &deriv)
420 {
421 }
422
424 {
426 return tmp;
427 }
428
429 SOLVER_UTILS_EXPORT virtual const Array<OneD,
432};
433
434} // namespace Nektar::SolverUtils
435
436#endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:357
void SetFluxVectorNS(DiffusionFluxVecCBNS fluxVector)
Definition: Diffusion.h:287
void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:334
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:209
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal()
Get trace normal.
Definition: Diffusion.h:262
void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:293
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:364
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs(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, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Definition: Diffusion.cpp:45
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Definition: Diffusion.h:367
virtual SOLVER_UTILS_EXPORT void v_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, const Array< OneD, Array< OneD, NekDouble > > &pBwd)=0
void SetFluxPenaltyNS(DiffusionFluxPenaltyNS flux)
Definition: Diffusion.h:300
SOLVER_UTILS_EXPORT void SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:249
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:356
virtual SOLVER_UTILS_EXPORT 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.
Definition: Diffusion.cpp:72
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)
Definition: Diffusion.h:152
SOLVER_UTILS_EXPORT void SetGridVelocityTrace(Array< OneD, Array< OneD, NekDouble > > &gridVelocityTrace)
Definition: Diffusion.h:348
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, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.h:165
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:358
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:222
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:361
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.h:139
void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:280
void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:320
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:268
void SetDiffusionFluxCons(DiffusionFluxCons flux)
Definition: Diffusion.h:314
SOLVER_UTILS_EXPORT TensorOfArray3D< NekDouble > & GetFluxTensor()
Definition: Diffusion.h:255
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:233
virtual SOLVER_UTILS_EXPORT ~Diffusion()=default
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:363
void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:340
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:365
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:418
SOLVER_UTILS_EXPORT void DiffuseCoeffs(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.h:193
virtual SOLVER_UTILS_EXPORT 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.
Definition: Diffusion.cpp:62
void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:306
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)=0
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:366
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:362
virtual SOLVER_UTILS_EXPORT 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.
Definition: Diffusion.cpp:82
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:423
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:56
void SetDiffusionFluxConsTrace(DiffusionFluxCons flux)
Definition: Diffusion.h:328
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:360
SOLVER_UTILS_EXPORT void DiffuseCoeffs(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)
Definition: Diffusion.h:179
SOLVER_UTILS_EXPORT void SetFluxVector(DiffusionFluxVecCB fluxVector)
Definition: Diffusion.h:274
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:55
std::function< void(Array< OneD, Array< OneD, NekDouble > > &)> SpecialBndTreat
Definition: Diffusion.h:134
std::function< void(const int, const Array< OneD, Array< OneD, NekDouble > > &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &, Array< OneD, int > &, const Array< OneD, Array< OneD, NekDouble > > &)> DiffusionFluxCons
Definition: Diffusion.h:110
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const Array< OneD, Array< OneD, NekDouble > > &, Array< OneD, Array< OneD, NekDouble > > &)> DiffusionFluxPenaltyNS
Definition: Diffusion.h:76
LibUtilities::NekFactory< std::string, Diffusion, std::string > DiffusionFactory
Datatype of the NekFactory used to instantiate classes derived from the Diffusion class.
Definition: Diffusion.h:60
std::function< void(const Array< OneD, const Array< OneD, NekDouble > > &, const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > &, NekDouble, const Array< OneD, const Array< OneD, NekDouble > > &, const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > &)> FunctorDerivBndCond
Definition: Diffusion.h:93
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCB
Definition: Diffusion.h:66
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCBNS
Definition: Diffusion.h:71
std::function< void(const int, const Array< OneD, Array< OneD, NekDouble > > &, const Array< OneD, Array< OneD, NekDouble > > &, TensorOfArray3D< NekDouble > &, Array< OneD, int > &, const Array< OneD, Array< OneD, NekDouble > > &)> DiffusionSymmFluxCons
Definition: Diffusion.h:127
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble