Nektar++
Diffusion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Diffusion.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: Abstract base class for diffusion.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 
37 namespace Nektar
38 {
39 namespace SolverUtils
40 {
42 {
43  static DiffusionFactory instance;
44  return instance;
45 }
46 
49 {
50  v_InitObject(pSession, pFields);
51 
52  // Div curl storage
53  int nPts = pFields[0]->GetTotPoints();
54  m_divVel = Array<OneD, NekDouble>(nPts, 0.0);
57 }
58 
60  const std::size_t nConvectiveFields,
62  const Array<OneD, Array<OneD, NekDouble>> &inarray,
63  Array<OneD, Array<OneD, NekDouble>> &outarray,
64  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
65  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
66 {
67  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
68 }
69 
70 /**
71  * @brief Similar with Diffusion::Diffuse(): calculate diffusion flux
72  * The difference is in the outarray:
73  * it is the coefficients of basis for DiffuseCoeffs()
74  * it is the physics on quadrature points for Diffuse()
75  */
77  const std::size_t nConvectiveFields,
79  const Array<OneD, Array<OneD, NekDouble>> &inarray,
80  Array<OneD, Array<OneD, NekDouble>> &outarray,
81  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
82  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
83 {
84  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
85 }
86 
88  const std::size_t nConvectiveFields,
90  const Array<OneD, Array<OneD, NekDouble>> &inarray,
91  Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
92  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
93  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
94 {
95  m_time = time;
96  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
97 }
98 
100  const std::size_t nConvectiveFields,
102  const Array<OneD, Array<OneD, NekDouble>> &inarray,
103  Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
104  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
105  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
106 {
107  m_time = time;
108  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
109 }
110 
112  const int nConvectiveFields, const int nDim, const int nPts,
113  const int nTracePts,
115  const Array<OneD, const int> &nonZeroIndex,
117  Array<OneD, Array<OneD, NekDouble>> &outarray)
118 {
119  boost::ignore_unused(nTracePts);
120 
121  int nCoeffs = outarray[nConvectiveFields - 1].size();
122  Array<OneD, NekDouble> tmpCoeff(nCoeffs, 0.0);
123  Array<OneD, Array<OneD, NekDouble>> tmpfield(nDim);
124 
125  for (int i = 0; i < nDim; i++)
126  {
127  tmpfield[i] = Array<OneD, NekDouble>(nPts, 0.0);
128  }
129  int nv = 0;
130 
131  for (int j = 0; j < nonZeroIndex.size(); j++)
132  {
133  nv = nonZeroIndex[j];
134  MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
135  for (int nd = 0; nd < nDim; nd++)
136  {
137  Vmath::Zero(nPts, tmpfield[nd], 1);
138 
139  tracelist->MultiplyByQuadratureMetric(Fwdflux[nd][nv],
140  Fwdflux[nd][nv]);
141  tracelist->MultiplyByQuadratureMetric(Bwdflux[nd][nv],
142  Bwdflux[nd][nv]);
143 
144  fields[nv]->AddTraceQuadPhysToOffDiag(
145  Fwdflux[nd][nv], Bwdflux[nd][nv], tmpfield[nd]);
146  fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
147  }
148  fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
149  Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
150  }
151 }
152 
154  const std::size_t nConvectiveFields,
156  const Array<OneD, Array<OneD, NekDouble>> &inarray,
157  Array<OneD, Array<OneD, NekDouble>> &outarray,
158  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
159  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
160 {
161  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, pFwd,
162  pBwd);
163  NEKERROR(ErrorUtil::efatal, "v_Diffuse not defined");
164 }
165 
167  const std::size_t nConvectiveFields,
169  const Array<OneD, Array<OneD, NekDouble>> &inarray,
170  Array<OneD, Array<OneD, NekDouble>> &outarray,
171  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
172  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
173 {
174  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, pFwd,
175  pBwd);
176  NEKERROR(ErrorUtil::efatal, "v_DiffuseCoeffs not defined");
177 }
178 
180  const std::size_t nConvectiveFields,
182  const Array<OneD, Array<OneD, NekDouble>> &inarray,
183  Array<OneD, Array<OneD, NekDouble>> &outarray,
184  const Array<OneD, Array<OneD, NekDouble>> &vFwd,
185  const Array<OneD, Array<OneD, NekDouble>> &vBwd,
186  TensorOfArray3D<NekDouble> &qfield, Array<OneD, int> &nonZeroIndex)
187 {
188  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, vFwd,
189  vBwd, qfield, nonZeroIndex);
190  NEKERROR(ErrorUtil::efatal, "v_DiffuseCoeffs not defined");
191 }
192 
193 // Diffusion Calculate the physical derivatives
196  const Array<OneD, Array<OneD, NekDouble>> &inarray,
198  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
199  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
200 {
201  v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
202 }
203 
205 {
206  NEKERROR(ErrorUtil::efatal, "v_GetTraceNormal not defined");
208 }
209 
211  const std::size_t nConvectiveFields, const size_t npnts,
212  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
213  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
216 {
217  boost::ignore_unused(nConvectiveFields, npnts, vFwd, vBwd, aver, jump);
218  NEKERROR(ErrorUtil::efatal, "v_ConsVarAveJump not defined");
219 }
220 
223  const Array<OneD, Array<OneD, NekDouble>> &inarray,
225  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
226  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
227 {
228  boost::ignore_unused(fields, inarray, qfields, pFwd, pBwd);
229  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
230 }
231 
234  const Array<OneD, Array<OneD, NekDouble>> &inarray,
236  Array<OneD, int> &nonZeroIndex)
237 {
238  boost::ignore_unused(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
239  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
240 }
241 
244  const Array<OneD, Array<OneD, NekDouble>> &inarray,
246  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
247  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
248  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
249  Array<OneD, int> &nonZeroIndex)
250 {
251  boost::ignore_unused(fields, inarray, qfields, VolumeFlux, TraceFlux, pFwd,
252  pBwd, nonZeroIndex);
253  NEKERROR(ErrorUtil::efatal, "Not defined function DiffuseTraceFLux.");
254 }
255 
256 /**
257  * @brief Compute primary variables
258  *
259  */
262  const Array<OneD, Array<OneD, NekDouble>> &inarray,
263  Array<OneD, Array<OneD, NekDouble>> &primVar)
264 {
265 
266  int nDim = fields[0]->GetCoordim(0);
267  int nPts = fields[0]->GetTotPoints();
268  for (int i = 0; i < nDim; ++i)
269  {
270  primVar[i] = Array<OneD, NekDouble>(nPts, 0.0);
271  Vmath::Vdiv(nPts, inarray[i + 1], 1, inarray[0], 1, primVar[i], 1);
272  }
273 }
274 
275 /**
276  * @brief Compute and store scalars needed for shock sensor, i.e.
277  * divergence and curl squared
278  * @param dimensions, points, primary variables derivative
279  *
280  */
283  const TensorOfArray3D<NekDouble> &pVarDer)
284 {
285  int nDim = fields[0]->GetCoordim(0);
286  int nPts = fields[0]->GetTotPoints();
287  // Compute and store scalars needed for shock sensor
288  Array<OneD, NekDouble> tmp3(nPts, 0.0);
289  Array<OneD, NekDouble> tmp4(nPts, 0.0);
290 
291  // Zero the variables for the current iteration step
292  Vmath::Zero(nPts, m_divVel, 1);
293  Vmath::Zero(nPts, m_divVelSquare, 1);
294  Vmath::Zero(nPts, m_curlVelSquare, 1);
295 
296  // div vel
297  for (int j = 0; j < nDim; ++j)
298  {
299  Vmath::Vadd(nPts, m_divVel, 1, pVarDer[j][j], 1, m_divVel, 1);
300  }
301  // (div vel)**2
302  Vmath::Vmul(nPts, m_divVel, 1, m_divVel, 1, m_divVelSquare, 1);
303 
304  // curl vel
305  if (nDim > 2)
306  {
307  // curl[0] 3/2
308  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][1], 1, tmp4, 1);
309  // curl[0]-2/3
310  Vmath::Vcopy(nPts, pVarDer[1][2], 1, tmp3, 1);
311  Vmath::Neg(nPts, tmp3, 1);
312  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
313  // square curl[0]
314  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
315  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
316 
317  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
318  // curl[1] 3/1
319  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][0], 1, tmp4, 1);
320  // curl[1]-1/3
321  Vmath::Vcopy(nPts, pVarDer[0][2], 1, tmp3, 1);
322  Vmath::Neg(nPts, tmp3, 1);
323  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
324  // square curl[1]
325  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
326  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
327 
328  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
329  // curl[2] 1/2
330  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
331  // curl[2]-2/1
332  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
333  Vmath::Neg(nPts, tmp3, 1);
334  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
335  // square curl[2]
336  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
337  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
338  }
339  else if (nDim > 1)
340  {
341  // curl[2] 1/2
342  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
343  // curl[2]-2/1
344  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
345  Vmath::Neg(nPts, tmp3, 1);
346  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
347  // square curl[2]
348  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
349  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
350  }
351  else
352  {
353  Vmath::Fill(nPts, 0.0, m_curlVelSquare, 1);
354  }
355 }
356 } // namespace SolverUtils
357 } // namespace Nektar
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
Provides a generic Factory class.
Definition: NekFactory.hpp:105
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:132
virtual SOLVER_UTILS_EXPORT void v_ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump)
Definition: Diffusion.cpp:210
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:221
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.cpp:59
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)
Similar with Diffusion::Diffuse(): calculate diffusion flux The difference is in the outarray: it is ...
Definition: Diffusion.cpp:76
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.cpp:194
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:232
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:166
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:131
SOLVER_UTILS_EXPORT void GetDivCurl(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
Compute divergence and curl squared.
Definition: Diffusion.cpp:281
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:133
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.cpp:47
virtual SOLVER_UTILS_EXPORT void v_GetPrimVar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
Compute primary derivatives.
Definition: Diffusion.cpp:260
SOLVER_UTILS_EXPORT void AddSymmFluxIntegralToOffDiag(const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &Fwdflux, TensorOfArray3D< NekDouble > &Bwdflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: Diffusion.cpp:111
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:242
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:204
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)
Definition: Diffusion.cpp:153
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.h:412
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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
double NekDouble
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:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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:359
void Vdiv(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:284
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255