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
47
48namespace Nektar::SolverUtils
49{
50
51class Diffusion;
52
53/// A shared pointer to an EquationSystem object
54typedef std::shared_ptr<Diffusion> DiffusionSharedPtr;
55
56/// Datatype of the NekFactory used to instantiate classes derived
57/// from the Diffusion class.
61
62typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
66
67typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
71
72typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
76
77/**
78 * Parameter list meaning:
79 * 1st: field conservative variables
80 * 2th: Devrivatives of field conservative varialbes
81 * 3rd: the current time for time-dependent boundary
82 * 4th: Fwd of field conservative variables optional
83 * 5th: Fwd of Devrivatives(2nd) optional
84 *
85 * a null pointer need to be passed for optional parameters
86 */
87typedef std::function<void(
88 const Array<OneD, const Array<OneD, NekDouble>> &,
90 const Array<OneD, const Array<OneD, NekDouble>> &,
91 const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &)>
93
94/**
95 * Parameter list meaning:
96 * 1st: nvariables
97 * 2nd: nspaceDimension
98 * 3rd: field conservative variables
99 * 4th: Devrivatives of field conservative varialbes
100 * 5th: nonzero flux index array, optional
101 * 6th: normal vectors optional
102 *
103 * a null pointer need to be passed for optional parameters
104 */
105typedef std::function<void(
106 const int, const Array<OneD, Array<OneD, NekDouble>> &,
110
111/**
112 * Parameter list meaning:
113 * 1st: nvariables
114 * 2nd: nspaceDimension
115 * 3rd: trace conservative variables for Diffusion Flux Jacobian
116 * 4th: trace conservative variables( usually the jump of trace value)
117 * 5th: trace symmetric flux
118 * 6th: nonzero flux index array, optional
119 *
120 * a null pointer need to be passed for optional parameters
121 */
122typedef std::function<void(
123 const int, const Array<OneD, Array<OneD, NekDouble>> &,
127
128/**
129 * Parameter list meaning:
130 * 1rd: trace conservative variables
131 */
132typedef std::function<void(Array<OneD, Array<OneD, NekDouble>> &)>
134
136{
137public:
139
143
145 const std::size_t nConvectiveFields,
147 const Array<OneD, Array<OneD, NekDouble>> &inarray,
149 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
151 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
153 {
154 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
155 }
156
158 const std::size_t nConvectiveFields,
160 const Array<OneD, Array<OneD, NekDouble>> &inarray,
161 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
162 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
164 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
166 {
167 m_time = time;
168 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
169 }
170
172 const std::size_t nConvectiveFields,
174 const Array<OneD, Array<OneD, NekDouble>> &inarray,
176 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
178 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
180 {
181 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
182 pBwd);
183 }
184
186 const std::size_t nConvectiveFields,
188 const Array<OneD, Array<OneD, NekDouble>> &inarray,
189 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
190 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
192 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
194 {
195 m_time = time;
196 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
197 pBwd);
198 }
199
200 // Diffusion Calculate the physical derivatives
203 const Array<OneD, Array<OneD, NekDouble>> &inarray,
205 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
207 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
209 {
210 v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
211 }
212
213 /// Diffusion Volume FLux
216 const Array<OneD, Array<OneD, NekDouble>> &inarray,
218 TensorOfArray3D<NekDouble> &VolumeFlux,
219 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
220 {
221 v_DiffuseVolumeFlux(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
222 }
223
224 /// Diffusion term Trace Flux
227 const Array<OneD, Array<OneD, NekDouble>> &inarray,
229 TensorOfArray3D<NekDouble> &VolumeFlux,
230 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
231 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
233 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
235 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
236 {
237 v_DiffuseTraceFlux(fields, inarray, qfields, VolumeFlux, TraceFlux,
238 pFwd, pBwd, nonZeroIndex);
239 }
240
243 {
244 v_SetHomoDerivs(deriv);
245 }
246
248 {
249 return v_GetFluxTensor();
250 }
251
252 /// Get trace normal
255 {
256 return v_GetTraceNormal();
257 }
258
259 template <typename FuncPointerT, typename ObjectPointerT>
260 void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
261 {
262 m_fluxVector = std::bind(func, obj, std::placeholders::_1,
263 std::placeholders::_2, std::placeholders::_3);
264 }
265
267 {
268 m_fluxVector = fluxVector;
269 }
270
271 template <typename FuncPointerT, typename ObjectPointerT>
272 void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
273 {
275 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
276 std::placeholders::_3);
277 }
278
280 {
281 m_fluxVectorNS = fluxVector;
282 }
283
284 template <typename FuncPointerT, typename ObjectPointerT>
285 void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
286 {
288 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
289 std::placeholders::_3);
290 }
291
293 {
294 m_fluxPenaltyNS = flux;
295 }
296
297 template <typename FuncPointerT, typename ObjectPointerT>
298 void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
299 {
301 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
302 std::placeholders::_3, std::placeholders::_4,
303 std::placeholders::_5, std::placeholders::_6);
304 }
305
307 {
309 }
310
311 template <typename FuncPointerT, typename ObjectPointerT>
312 void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
313 {
315 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
316 std::placeholders::_3, std::placeholders::_4,
317 std::placeholders::_5, std::placeholders::_6);
318 }
319
321 {
323 }
324
325 template <typename FuncPointerT, typename ObjectPointerT>
326 void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
327 {
328 m_SpecialBndTreat = std::bind(func, obj, std::placeholders::_1);
329 }
330
331 template <typename FuncPointerT, typename ObjectPointerT>
332 void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
333 {
335 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
336 std::placeholders::_3, std::placeholders::_4,
337 std::placeholders::_5, std::placeholders::_6);
338 }
339
340protected:
341 /// Params for Ducros sensor
345
353
355
359
361 const std::size_t nConvectiveFields,
363 const Array<OneD, Array<OneD, NekDouble>> &inarray,
365 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
366 const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
367
369 const std::size_t nConvectiveFields,
371 const Array<OneD, Array<OneD, NekDouble>> &inarray,
373 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
374 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
375
376 /// Diffusion Flux, calculate the physical derivatives
379 const Array<OneD, Array<OneD, NekDouble>> &inarray,
381 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
382 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
383
384 /// Diffusion Volume Flux
387 const Array<OneD, Array<OneD, NekDouble>> &inarray,
389 TensorOfArray3D<NekDouble> &VolumeFlux, Array<OneD, int> &nonZeroIndex);
390
391 /// Diffusion term Trace Flux
394 const Array<OneD, Array<OneD, NekDouble>> &inarray,
396 TensorOfArray3D<NekDouble> &VolumeFlux,
397 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
398 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
399 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
400 Array<OneD, int> &nonZeroIndex);
401
402 virtual void v_SetHomoDerivs(
403 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &deriv)
404 {
405 }
406
408 {
410 return tmp;
411 }
412
415};
416
417} // namespace Nektar::SolverUtils
418
419#endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:104
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:343
void SetFluxVectorNS(DiffusionFluxVecCBNS fluxVector)
Definition: Diffusion.h:279
void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:326
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:201
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal()
Get trace normal.
Definition: Diffusion.h:254
void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:285
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:350
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:57
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:292
SOLVER_UTILS_EXPORT void SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:241
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:342
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:84
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:144
virtual SOLVER_UTILS_EXPORT ~Diffusion()
Definition: Diffusion.h:138
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:157
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:344
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:214
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:347
void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:272
void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:312
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:260
void SetDiffusionFluxCons(DiffusionFluxCons flux)
Definition: Diffusion.h:306
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.cpp:45
SOLVER_UTILS_EXPORT TensorOfArray3D< NekDouble > & GetFluxTensor()
Definition: Diffusion.h:247
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:225
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:349
void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:332
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:351
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:402
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:185
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:74
void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:298
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)=0
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:352
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:348
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:94
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:407
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:68
void SetDiffusionFluxConsTrace(DiffusionFluxCons flux)
Definition: Diffusion.h:320
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:346
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:171
SOLVER_UTILS_EXPORT void SetFluxVector(DiffusionFluxVecCB fluxVector)
Definition: Diffusion.h:266
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:54
std::function< void(Array< OneD, Array< OneD, NekDouble > > &)> SpecialBndTreat
Definition: Diffusion.h:133
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:109
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:75
LibUtilities::NekFactory< std::string, Diffusion, std::string > DiffusionFactory
Datatype of the NekFactory used to instantiate classes derived from the Diffusion class.
Definition: Diffusion.h:59
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:92
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCB
Definition: Diffusion.h:65
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCBNS
Definition: Diffusion.h:70
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:126
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble