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