Nektar++
Advection.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Advection.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 advection.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
38
39namespace Nektar
40{
41namespace SolverUtils
42{
43
44/**
45 * @returns The advection factory.
46 */
48{
49 static AdvectionFactory instance;
50 return instance;
51}
52
53/**
54 * @param pSession Session configuration data.
55 * @param pFields Array of ExpList objects.
56 */
59{
60 v_InitObject(pSession, pFields);
61}
62
63/**
64 * @param nConvectiveFields Number of velocity components.
65 * @param pFields Expansion lists for scalar fields.
66 * @param pAdvVel Advection velocity.
67 * @param pInarray Scalar data to advect.
68 * @param pOutarray Advected scalar data.
69 * @param pTime Simulation time.
70 */
72 const int nConvectiveFields,
74 const Array<OneD, Array<OneD, NekDouble>> &pAdvVel,
75 const Array<OneD, Array<OneD, NekDouble>> &pInarray,
76 Array<OneD, Array<OneD, NekDouble>> &pOutarray, const NekDouble &pTime,
77 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
78 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
79{
80 v_Advect(nConvectiveFields, pFields, pAdvVel, pInarray, pOutarray, pTime,
81 pFwd, pBwd);
82}
83
84template <typename DataType, typename TypeNekBlkMatSharedPtr>
86 const int nConvectiveFields, const int nSpaceDim,
88 const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacCons,
90 const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacGrad,
91 const Array<OneD, Array<OneD, DataType>> &TracePntJacGradSign)
92{
93 MultiRegions::ExpListSharedPtr tracelist = pFields[0]->GetTrace();
94 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
95 tracelist->GetExp();
96 int ntotTrac = (*traceExp).size();
97 int nTracePnts = tracelist->GetTotPoints();
98 int nTracPnt, nTracCoef;
99
100 DNekMatSharedPtr tmp2Add;
101 Array<OneD, DataType> MatData0;
102 Array<OneD, NekDouble> MatData1;
103 Array<OneD, NekDouble> MatData2;
104
105 Array<OneD, DNekMatSharedPtr> TraceJacFwd(ntotTrac);
106 Array<OneD, DNekMatSharedPtr> TraceJacBwd(ntotTrac);
107
108 for (int nelmt = 0; nelmt < ntotTrac; nelmt++)
109 {
110 nTracCoef = (*traceExp)[nelmt]->GetNcoeffs();
111 nTracPnt = (*traceExp)[nelmt]->GetTotPoints();
112 TraceJacFwd[nelmt] =
113 MemoryManager<DNekMat>::AllocateSharedPtr(nTracCoef, nTracPnt, 0.0);
114 TraceJacBwd[nelmt] =
115 MemoryManager<DNekMat>::AllocateSharedPtr(nTracCoef, nTracPnt, 0.0);
116 }
117
118 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
119 pFields[0]->GetExp();
120 int nTotElmt = (*fieldExp).size();
121 int nElmtPnt, nElmtCoef;
122
123 Array<OneD, DNekMatSharedPtr> ElmtJacQuad(nTotElmt);
124 Array<OneD, DNekMatSharedPtr> ElmtJacCoef(nTotElmt);
125
126 Array<OneD, NekDouble> SymmMatData;
127
128 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
129 {
130 nElmtCoef = (*fieldExp)[nelmt]->GetNcoeffs();
131 nElmtPnt = (*fieldExp)[nelmt]->GetTotPoints();
132
133 ElmtJacQuad[nelmt] =
134 MemoryManager<DNekMat>::AllocateSharedPtr(nElmtCoef, nElmtPnt, 0.0);
136 nElmtCoef, nElmtCoef, 0.0);
137 }
138
139 bool TracePntJacGradflag = true;
140
141 Array<OneD, Array<OneD, DataType>> TraceJacConsSign(2);
142 for (int i = 0; i < 2; i++)
143 {
144 TraceJacConsSign[i] = Array<OneD, DataType>(nTracePnts, 1.0);
145 }
146
147 if (0 == TracePntJacGrad.size())
148 {
149 TracePntJacGradflag = false;
150 }
151
152 for (int m = 0; m < nConvectiveFields; m++)
153 {
154 for (int n = 0; n < nConvectiveFields; n++)
155 {
156 // ElmtJacCons to set 0
157 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
158 {
159 (*ElmtJacCoef[nelmt]) = 0.0;
160 (*ElmtJacQuad[nelmt]) = 0.0;
161 }
162
163 if (TracePntJacGradflag)
164 {
165 // only BJac is stored here because BJac = - FJac
166 for (int ndir = 0; ndir < nSpaceDim; ndir++)
167 {
168 // ElmtJacGrad to set 0
169 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
170 {
171 (*ElmtJacQuad[nelmt]) = 0.0;
172 }
173 int ngrad = n * nSpaceDim + ndir;
174
175 CalcJacobTraceInteg(pFields, m, ngrad, TracePntJacGrad,
176 TracePntJacGradSign, TraceJacFwd,
177 TraceJacBwd);
178 pFields[0]->AddTraceJacToElmtJac(TraceJacFwd, TraceJacBwd,
179 ElmtJacQuad);
180 // cost can be lower down by doing 3 directions together
181 pFields[0]->AddRightIPTPhysDerivBase(ndir, ElmtJacQuad,
182 ElmtJacCoef);
183 }
184 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
185 {
186 nElmtCoef = (*fieldExp)[nelmt]->GetNcoeffs();
187 if (SymmMatData.size() < nElmtCoef)
188 {
189 SymmMatData =
190 Array<OneD, NekDouble>(nElmtCoef * nElmtCoef);
191 }
192 MatData1 = ElmtJacCoef[nelmt]->GetPtr();
193 for (int i = 0; i < nElmtCoef; i++)
194 {
195 Vmath::Vcopy(nElmtCoef, &MatData1[i * nElmtCoef], 1,
196 &SymmMatData[i], nElmtCoef);
197 }
198 Vmath::Vadd(nElmtCoef * nElmtCoef, SymmMatData, 1, MatData1,
199 1, MatData1, 1);
200 }
201 }
202
203 CalcJacobTraceInteg(pFields, m, n, TracePntJacCons,
204 TraceJacConsSign, TraceJacFwd, TraceJacBwd);
205
206 pFields[0]->AddTraceJacToElmtJac(TraceJacFwd, TraceJacBwd,
207 ElmtJacQuad);
208 // Memory can be save to explore reusing ElmtJacQuad as output
209 pFields[0]->AddRightIPTBaseMatrix(ElmtJacQuad, ElmtJacCoef);
210
211 // add ElmtJacCons to gmtxarray[m][n]
212 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
213 {
214 tmp2Add = ElmtJacCoef[nelmt];
215 MatData0 = gmtxarray[m][n]->GetBlock(nelmt, nelmt)->GetPtr();
216 MatData1 = tmp2Add->GetPtr();
217 for (int i = 0; i < MatData0.size(); i++)
218 {
219 MatData0[i] += DataType(MatData1[i]);
220 }
221 }
222 }
223 }
224}
225
227 const int nConvectiveFields, const int nSpaceDim,
229 const Array<OneD, DNekBlkMatSharedPtr> &TracePntJacCons,
231 const Array<OneD, DNekBlkMatSharedPtr> &TracePntJacGrad,
232 const Array<OneD, Array<OneD, NekDouble>> &TracePntJacGradSign);
234 const int nConvectiveFields, const int nSpaceDim,
236 const Array<OneD, SNekBlkMatSharedPtr> &TracePntJacCons,
238 const Array<OneD, SNekBlkMatSharedPtr> &TracePntJacGrad,
239 const Array<OneD, Array<OneD, NekSingle>> &TracePntJacGradSign);
240
241/**
242 *
243 */
244template <typename DataType, typename TypeNekBlkMatSharedPtr>
246 const Array<OneD, MultiRegions::ExpListSharedPtr> &pFields, const int m,
247 const int n, const Array<OneD, const TypeNekBlkMatSharedPtr> &PntJac,
248 const Array<OneD, const Array<OneD, DataType>> &PntJacSign,
251{
252 MultiRegions::ExpListSharedPtr tracelist = pFields[0]->GetTrace();
253 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
254 tracelist->GetExp();
255 int ntotTrac = (*traceExp).size();
256 int nTracPnt, noffset, pntoffset;
257
258 Array<OneD, int> tracepnts(ntotTrac);
259 Array<OneD, Array<OneD, NekDouble>> JacFwd(ntotTrac);
260 Array<OneD, Array<OneD, NekDouble>> JacBwd(ntotTrac);
261
262 for (int nelmt = 0; nelmt < ntotTrac; nelmt++)
263 {
264 nTracPnt = (*traceExp)[nelmt]->GetTotPoints();
265 tracepnts[nelmt] = nTracPnt;
266
267 JacFwd[nelmt] = Array<OneD, NekDouble>(nTracPnt, 0.0);
268 JacBwd[nelmt] = Array<OneD, NekDouble>(nTracPnt, 0.0);
269 }
270
271 DataType ftmp, btmp;
272 for (int nelmt = 0; nelmt < ntotTrac; nelmt++)
273 {
274 nTracPnt = tracepnts[nelmt];
275 noffset = tracelist->GetPhys_Offset(nelmt);
276 for (int npnt = 0; npnt < nTracPnt; npnt++)
277 {
278 pntoffset = noffset + npnt;
279 ftmp = (*PntJac[0]->GetBlock(pntoffset, pntoffset))(m, n);
280 JacFwd[nelmt][npnt] = NekDouble(PntJacSign[0][pntoffset] * ftmp);
281
282 btmp = (*PntJac[1]->GetBlock(pntoffset, pntoffset))(m, n);
283 JacBwd[nelmt][npnt] = NekDouble(PntJacSign[1][pntoffset] * btmp);
284 }
285 }
286 tracelist->GetDiagMatIpwrtBase(JacFwd, TraceJacFwd);
287 tracelist->GetDiagMatIpwrtBase(JacBwd, TraceJacBwd);
288}
289
290/**
291 * This function should be overridden in derived classes to initialise the
292 * specific advection data members. However, this base class function should
293 * be called as the first statement of the overridden function to ensure the
294 * base class is correctly initialised in order.
295 *
296 * @param pSession Session information.
297 * @param pFields Expansion lists for scalar fields.
298 */
302{
303 m_spaceDim = pFields[0]->GetCoordim(0);
304
305 if (pSession->DefinesSolverInfo("HOMOGENEOUS"))
306 {
307 std::string HomoStr = pSession->GetSolverInfo("HOMOGENEOUS");
308 if (HomoStr == "HOMOGENEOUS1D" || HomoStr == "Homogeneous1D" ||
309 HomoStr == "1D" || HomoStr == "Homo1D" ||
310 HomoStr == "HOMOGENEOUS2D" || HomoStr == "Homogeneous2D" ||
311 HomoStr == "2D" || HomoStr == "Homo2D")
312 {
313 m_spaceDim++;
314 }
315 else
316 {
318 "Only 1D homogeneous dimension supported.");
319 }
320 }
321}
322
323/**
324 *
325 */
327 const Array<OneD, Array<OneD, NekDouble>> &inarray,
329{
330 boost::ignore_unused(inarray, fields);
332 "A baseflow is not appropriate for this advection type.");
333}
334
335} // namespace SolverUtils
336} // 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
#define SOLVER_UTILS_EXPORT
size_type size() const
Returns the array's size.
Provides a generic Factory class.
Definition: NekFactory.hpp:105
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:174
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:299
void CalcJacobTraceInteg(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int m, const int n, const Array< OneD, const TypeNekBlkMatSharedPtr > &PntJac, const Array< OneD, const Array< OneD, DataType > > &PntJacSign, Array< OneD, DNekMatSharedPtr > &TraceJacFwd, Array< OneD, DNekMatSharedPtr > &TraceJacBwd)
Definition: Advection.cpp:245
virtual SOLVER_UTILS_EXPORT void v_Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)=0
Advects a vector field.
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Interface function to initialise the advection object.
Definition: Advection.cpp:57
virtual SOLVER_UTILS_EXPORT void v_SetBaseFlow(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Overrides the base flow used during linearised advection.
Definition: Advection.cpp:326
SOLVER_UTILS_EXPORT void AddTraceJacToMat(const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType > > &TracePntJacGradSign)
SOLVER_UTILS_EXPORT void Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
Interface function to advect the vector field.
Definition: Advection.cpp:71
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191