Nektar++
AdvectionWeakDG.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: AdvectionWeakDG.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: Weak DG advection class.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <iomanip>
39 #include <iostream>
40 
42 
43 namespace Nektar
44 {
45 namespace SolverUtils
46 {
47 std::string AdvectionWeakDG::type =
49  "WeakDG", AdvectionWeakDG::create, "Weak DG");
50 
52 {
53 }
54 
55 /**
56  * @brief Initialise AdvectionWeakDG objects and store them before
57  * starting the time-stepping.
58  *
59  * @param pSession Pointer to session reader.
60  * @param pFields Pointer to fields.
61  */
65 {
66  Advection::v_InitObject(pSession, pFields);
67 }
68 
69 /**
70  * @brief Compute the advection term at each time-step using the
71  * Discontinuous Galerkin approach (DG).
72  *
73  * @param nConvectiveFields Number of fields.
74  * @param fields Pointer to fields.
75  * @param advVel Advection velocities.
76  * @param inarray Solution at the previous time-step.
77  * @param outarray Advection term to be passed at the
78  * time integration class.
79  */
81  const int nConvectiveFields,
83  const Array<OneD, Array<OneD, NekDouble>> &advVel,
84  const Array<OneD, Array<OneD, NekDouble>> &inarray,
85  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
86  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
87  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
88 {
89  LibUtilities::Timer timer1;
90  timer1.Start();
91 
92  boost::ignore_unused(advVel, time);
93  size_t nCoeffs = fields[0]->GetNcoeffs();
94 
95  Array<OneD, Array<OneD, NekDouble>> tmp{size_t(nConvectiveFields)};
96  for (int i = 0; i < nConvectiveFields; ++i)
97  {
98  tmp[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
99  }
100  LibUtilities::Timer timer;
101  timer.Start();
102  AdvectionWeakDG::v_AdvectCoeffs(nConvectiveFields, fields, advVel, inarray,
103  tmp, time, pFwd, pBwd);
104  timer.Stop();
105  timer.AccumulateRegion("AdvWeakDG:v_AdvectCoeffs", 2);
106 
107  // why was this broken in many loops over convective fields?
108  // this is terrible for locality
109 
110  timer.Start();
111  for (int i = 0; i < nConvectiveFields; ++i)
112  {
113  fields[i]->BwdTrans(tmp[i], outarray[i]);
114  }
115  timer.Stop();
116  timer.AccumulateRegion("AdvWeakDG:_BwdTrans", 10);
117  timer1.Stop();
118  timer1.AccumulateRegion("AdvWeakDG:All", 10);
119 }
120 
122  const int nConvectiveFields,
124  const Array<OneD, Array<OneD, NekDouble>> &advVel,
125  const Array<OneD, Array<OneD, NekDouble>> &inarray,
126  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
127  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
128  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
129 {
130  size_t nPointsTot = fields[0]->GetTotPoints();
131  size_t nCoeffs = fields[0]->GetNcoeffs();
132  size_t nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
133 
135  nConvectiveFields);
136  // Allocate storage for flux vector F(u).
137  for (int i = 0; i < nConvectiveFields; ++i)
138  {
139  fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>{size_t(m_spaceDim)};
140  for (int j = 0; j < m_spaceDim; ++j)
141  {
142  fluxvector[i][j] = Array<OneD, NekDouble>(nPointsTot, 0.0);
143  }
144  }
145 
146  LibUtilities::Timer timer;
147  timer.Start();
148  v_AdvectVolumeFlux(nConvectiveFields, fields, advVel, inarray, fluxvector,
149  time);
150  timer.Stop();
151  timer.AccumulateRegion("AdvWeakDG:_fluxVector", 3);
152 
153  timer.Start();
154  // Get the advection part (without numerical flux)
155  for (int i = 0; i < nConvectiveFields; ++i)
156  {
157  Vmath::Fill(outarray[i].size(), 0.0, outarray[i], 1);
158  fields[i]->IProductWRTDerivBase(fluxvector[i], outarray[i]);
159  }
160  timer.Stop();
161  timer.AccumulateRegion("AdvWeakDG:_IProductWRTDerivBase", 10);
162 
163  Array<OneD, Array<OneD, NekDouble>> numflux{size_t(nConvectiveFields)};
164 
165  for (int i = 0; i < nConvectiveFields; ++i)
166  {
167  numflux[i] = Array<OneD, NekDouble>{nTracePointsTot, 0.0};
168  }
169 
170  v_AdvectTraceFlux(nConvectiveFields, fields, advVel, inarray, numflux, time,
171  pFwd, pBwd);
172 
173  // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
174  for (int i = 0; i < nConvectiveFields; ++i)
175  {
176  Vmath::Neg(nCoeffs, outarray[i], 1);
177 
178  timer.Start();
179  fields[i]->AddTraceIntegral(numflux[i], outarray[i]);
180  timer.Stop();
181  timer.AccumulateRegion("AdvWeakDG:_AddTraceIntegral", 10);
182 
183  timer.Start();
184  fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
185  timer.Stop();
186  timer.AccumulateRegion("AdvWeakDG:_MultiplyByElmtInvMass", 10);
187  }
188 }
189 
191  const int nConvectiveFields,
193  const Array<OneD, Array<OneD, NekDouble>> &advVel,
194  const Array<OneD, Array<OneD, NekDouble>> &inarray,
195  Array<OneD, Array<OneD, NekDouble>> &TraceFlux, const NekDouble &time,
196  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
197  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
198 {
199  boost::ignore_unused(advVel, time);
200  int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
201 
202  ASSERTL1(m_riemann, "Riemann solver must be provided for AdvectionWeakDG.");
203 
204  // Store forwards/backwards space along trace space
205  Array<OneD, Array<OneD, NekDouble>> Fwd(nConvectiveFields);
206  Array<OneD, Array<OneD, NekDouble>> Bwd(nConvectiveFields);
207 
209  {
210  for (int i = 0; i < nConvectiveFields; ++i)
211  {
212  Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
213  Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
214  fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
215  }
216  }
217  else
218  {
219  for (int i = 0; i < nConvectiveFields; ++i)
220  {
221  Fwd[i] = pFwd[i];
222  Bwd[i] = pBwd[i];
223  }
224  }
225 
226  LibUtilities::Timer timer;
227  timer.Start();
228  m_riemann->Solve(m_spaceDim, Fwd, Bwd, TraceFlux);
229  timer.Stop();
230  timer.AccumulateRegion("AdvWeakDG:_Riemann", 10);
231 }
232 
235  const int &nConvectiveFields,
236  const TensorOfArray5D<NekDouble> &ElmtJacArray,
238 {
239  MultiRegions::ExpListSharedPtr explist = pFields[0];
240  std::shared_ptr<LocalRegions::ExpansionVector> pexp = explist->GetExp();
241  int nTotElmt = (*pexp).size();
242  int nElmtPnt, nElmtCoef;
243 
244  SNekMatSharedPtr tmpGmtx;
245  DNekMatSharedPtr ElmtMat;
246 
247  Array<OneD, NekSingle> GMat_data;
248  Array<OneD, NekDouble> Elmt_data;
249  Array<OneD, NekSingle> Elmt_dataSingle;
250 
251  Array<OneD, DNekMatSharedPtr> mtxPerVar(nTotElmt);
252  Array<OneD, int> elmtpnts(nTotElmt);
253  Array<OneD, int> elmtcoef(nTotElmt);
254  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
255  {
256  nElmtCoef = (*pexp)[nelmt]->GetNcoeffs();
257  nElmtPnt = (*pexp)[nelmt]->GetTotPoints();
258  elmtpnts[nelmt] = nElmtPnt;
259  elmtcoef[nelmt] = nElmtCoef;
260  mtxPerVar[nelmt] =
261  MemoryManager<DNekMat>::AllocateSharedPtr(nElmtCoef, nElmtPnt);
262  }
263 
264  Array<OneD, DNekMatSharedPtr> mtxPerVarCoeff(nTotElmt);
265  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
266  {
267  nElmtCoef = elmtcoef[nelmt];
268  mtxPerVarCoeff[nelmt] =
269  MemoryManager<DNekMat>::AllocateSharedPtr(nElmtCoef, nElmtCoef);
270  }
271 
272  for (int m = 0; m < nConvectiveFields; m++)
273  {
274  for (int n = 0; n < nConvectiveFields; n++)
275  {
276  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
277  {
278  (*mtxPerVarCoeff[nelmt]) = 0.0;
279  (*mtxPerVar[nelmt]) = 0.0;
280  }
281 
282  explist->GetMatIpwrtDeriveBase(ElmtJacArray[m][n], mtxPerVar);
283  // Memory can be saved by reusing mtxPerVar
284  explist->AddRightIPTBaseMatrix(mtxPerVar, mtxPerVarCoeff);
285 
286  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
287  {
288  nElmtCoef = elmtcoef[nelmt];
289  nElmtPnt = elmtpnts[nelmt];
290  int ntotDofs = nElmtCoef * nElmtCoef;
291 
292  if (Elmt_dataSingle.size() < ntotDofs)
293  {
294  Elmt_dataSingle = Array<OneD, NekSingle>(ntotDofs);
295  }
296 
297  tmpGmtx = gmtxarray[m][n]->GetBlock(nelmt, nelmt);
298  ElmtMat = mtxPerVarCoeff[nelmt];
299 
300  GMat_data = tmpGmtx->GetPtr();
301  Elmt_data = ElmtMat->GetPtr();
302 
303  for (int i = 0; i < ntotDofs; i++)
304  {
305  Elmt_dataSingle[i] = NekSingle(Elmt_data[i]);
306  }
307 
308  Vmath::Vadd(ntotDofs, GMat_data, 1, Elmt_dataSingle, 1,
309  GMat_data, 1);
310  }
311  }
312  }
313 }
314 
315 } // end of namespace SolverUtils
316 } // end of namespace Nektar
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:73
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:215
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:353
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:213
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Initialise AdvectionWeakDG objects and store them before starting the time-stepping.
virtual void v_AddVolumJacToMat(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray) override
virtual void v_AdvectTraceFlux(const int nConvective, 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 >> &TraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray) override
Advects Trace Flux.
virtual void v_AdvectVolumeFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &VolumeFlux, const NekDouble &time) override
Advects Volume Flux.
static AdvectionSharedPtr create(std::string advType)
virtual void v_Advect(const int nConvective, 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) override
Compute the advection term at each time-step using the Discontinuous Galerkin approach (DG).
virtual void v_AdvectCoeffs(const int nConvective, 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) override
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< SNekMat > SNekMatSharedPtr
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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 Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45