Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
41#include <boost/core/ignore_unused.hpp>
42
49
50namespace Nektar
51{
52namespace SolverUtils
53{
54
55class Diffusion;
56
57/// A shared pointer to an EquationSystem object
58typedef std::shared_ptr<Diffusion> DiffusionSharedPtr;
59
60/// Datatype of the NekFactory used to instantiate classes derived
61/// from the Diffusion class.
65
66typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
70
71typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
75
76typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
80
81/**
82 * Parameter list meaning:
83 * 1st: field conservative variables
84 * 2th: Devrivatives of field conservative varialbes
85 * 3rd: the current time for time-dependent boundary
86 * 4th: Fwd of field conservative variables optional
87 * 5th: Fwd of Devrivatives(2nd) optional
88 *
89 * a null pointer need to be passed for optional parameters
90 */
91typedef std::function<void(
92 const Array<OneD, const Array<OneD, NekDouble>> &,
94 const Array<OneD, const Array<OneD, NekDouble>> &,
95 const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &)>
97
98/**
99 * Parameter list meaning:
100 * 1st: nvariables
101 * 2nd: nspaceDimension
102 * 3rd: field conservative variables
103 * 4th: Devrivatives of field conservative varialbes
104 * 5th: nonzero flux index array, optional
105 * 6th: normal vectors optional
106 *
107 * a null pointer need to be passed for optional parameters
108 */
109typedef std::function<void(
110 const int, const Array<OneD, Array<OneD, NekDouble>> &,
114
115/**
116 * Parameter list meaning:
117 * 1st: nvariables
118 * 2nd: nspaceDimension
119 * 3rd: trace conservative variables for Diffusion Flux Jacobian
120 * 4th: trace conservative variables( usually the jump of trace value)
121 * 5th: trace symmetric flux
122 * 6th: nonzero flux index array, optional
123 *
124 * a null pointer need to be passed for optional parameters
125 */
126typedef std::function<void(
127 const int, const Array<OneD, Array<OneD, NekDouble>> &,
131
132/**
133 * Parameter list meaning:
134 * 1rd: trace conservative variables
135 */
136typedef std::function<void(Array<OneD, Array<OneD, NekDouble>> &)>
138
140{
141public:
143
147
149 const std::size_t nConvectiveFields,
151 const Array<OneD, Array<OneD, NekDouble>> &inarray,
153 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
155 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
157 {
158 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
159 }
160
162 const std::size_t nConvectiveFields,
164 const Array<OneD, Array<OneD, NekDouble>> &inarray,
165 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
166 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
168 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
170 {
171 m_time = time;
172 v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
173 }
174
176 const std::size_t nConvectiveFields,
178 const Array<OneD, Array<OneD, NekDouble>> &inarray,
180 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
182 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
184 {
185 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
186 pBwd);
187 }
188
190 const std::size_t nConvectiveFields,
192 const Array<OneD, Array<OneD, NekDouble>> &inarray,
193 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
194 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
196 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
198 {
199 m_time = time;
200 v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
201 pBwd);
202 }
203
204 // Diffusion Calculate the physical derivatives
207 const Array<OneD, Array<OneD, NekDouble>> &inarray,
209 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
211 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
213 {
214 v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
215 }
216
217 /// Diffusion Volume FLux
220 const Array<OneD, Array<OneD, NekDouble>> &inarray,
222 TensorOfArray3D<NekDouble> &VolumeFlux,
223 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
224 {
225 v_DiffuseVolumeFlux(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
226 }
227
228 /// Diffusion term Trace Flux
231 const Array<OneD, Array<OneD, NekDouble>> &inarray,
233 TensorOfArray3D<NekDouble> &VolumeFlux,
234 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
235 const Array<OneD, Array<OneD, NekDouble>> &pFwd =
237 const Array<OneD, Array<OneD, NekDouble>> &pBwd =
239 Array<OneD, int> &nonZeroIndex = NullInt1DArray)
240 {
241 v_DiffuseTraceFlux(fields, inarray, qfields, VolumeFlux, TraceFlux,
242 pFwd, pBwd, nonZeroIndex);
243 }
244
247 {
248 v_SetHomoDerivs(deriv);
249 }
250
252 {
253 return v_GetFluxTensor();
254 }
255
256 /// Get trace normal
259 {
260 return v_GetTraceNormal();
261 }
262
263 template <typename FuncPointerT, typename ObjectPointerT>
264 void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
265 {
266 m_fluxVector = std::bind(func, obj, std::placeholders::_1,
267 std::placeholders::_2, std::placeholders::_3);
268 }
269
271 {
272 m_fluxVector = fluxVector;
273 }
274
275 template <typename FuncPointerT, typename ObjectPointerT>
276 void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
277 {
279 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
280 std::placeholders::_3);
281 }
282
284 {
285 m_fluxVectorNS = fluxVector;
286 }
287
288 template <typename FuncPointerT, typename ObjectPointerT>
289 void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
290 {
292 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
293 std::placeholders::_3);
294 }
295
297 {
298 m_fluxPenaltyNS = flux;
299 }
300
301 template <typename FuncPointerT, typename ObjectPointerT>
302 void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
303 {
305 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
306 std::placeholders::_3, std::placeholders::_4,
307 std::placeholders::_5, std::placeholders::_6);
308 }
309
311 {
313 }
314
315 template <typename FuncPointerT, typename ObjectPointerT>
316 void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
317 {
319 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
320 std::placeholders::_3, std::placeholders::_4,
321 std::placeholders::_5, std::placeholders::_6);
322 }
323
325 {
327 }
328
329 template <typename FuncPointerT, typename ObjectPointerT>
330 void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
331 {
332 m_SpecialBndTreat = std::bind(func, obj, std::placeholders::_1);
333 }
334
335 template <typename FuncPointerT, typename ObjectPointerT>
336 void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
337 {
339 std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
340 std::placeholders::_3, std::placeholders::_4,
341 std::placeholders::_5, std::placeholders::_6);
342 }
343
344protected:
345 /// Params for Ducros sensor
349
357
359
363
365 const std::size_t nConvectiveFields,
367 const Array<OneD, Array<OneD, NekDouble>> &inarray,
369 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
370 const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
371
373 const std::size_t nConvectiveFields,
375 const Array<OneD, Array<OneD, NekDouble>> &inarray,
377 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
378 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
379
380 /// Diffusion Flux, calculate the physical derivatives
383 const Array<OneD, Array<OneD, NekDouble>> &inarray,
385 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
386 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
387
388 /// Diffusion Volume Flux
391 const Array<OneD, Array<OneD, NekDouble>> &inarray,
393 TensorOfArray3D<NekDouble> &VolumeFlux, Array<OneD, int> &nonZeroIndex);
394
395 /// Diffusion term Trace Flux
398 const Array<OneD, Array<OneD, NekDouble>> &inarray,
400 TensorOfArray3D<NekDouble> &VolumeFlux,
401 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
402 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
403 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
404 Array<OneD, int> &nonZeroIndex);
405
407 {
408 boost::ignore_unused(deriv);
409 }
410
412 {
414 return tmp;
415 }
416
419};
420
421} // namespace SolverUtils
422} // namespace Nektar
423
424#endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:347
void SetFluxVectorNS(DiffusionFluxVecCBNS fluxVector)
Definition: Diffusion.h:283
void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:330
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 const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal()
Get trace normal.
Definition: Diffusion.h:258
void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:289
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:354
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:59
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:296
SOLVER_UTILS_EXPORT void SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:245
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:346
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:89
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:148
virtual SOLVER_UTILS_EXPORT ~Diffusion()
Definition: Diffusion.h:142
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:161
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:348
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
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:351
void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:276
void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:316
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:264
void SetDiffusionFluxCons(DiffusionFluxCons flux)
Definition: Diffusion.h:310
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.cpp:47
SOLVER_UTILS_EXPORT TensorOfArray3D< NekDouble > & GetFluxTensor()
Definition: Diffusion.h:251
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
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:353
void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:336
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:355
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:406
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:189
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:78
void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:302
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)=0
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:356
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:352
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:99
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:411
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:72
void SetDiffusionFluxConsTrace(DiffusionFluxCons flux)
Definition: Diffusion.h:324
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:350
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:175
SOLVER_UTILS_EXPORT void SetFluxVector(DiffusionFluxVecCB fluxVector)
Definition: Diffusion.h:270
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:58
std::function< void(Array< OneD, Array< OneD, NekDouble > > &)> SpecialBndTreat
Definition: Diffusion.h:137
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:113
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const Array< OneD, Array< OneD, NekDouble > > &, Array< OneD, Array< OneD, NekDouble > > &)> DiffusionFluxPenaltyNS
Definition: Diffusion.h:79
LibUtilities::NekFactory< std::string, Diffusion, std::string > DiffusionFactory
Datatype of the NekFactory used to instantiate classes derived from the Diffusion class.
Definition: Diffusion.h:63
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:96
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCB
Definition: Diffusion.h:69
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCBNS
Definition: Diffusion.h:74
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:130
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble