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