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
43namespace Nektar
44{
45namespace SolverUtils
46{
47std::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{
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 }
101 timer.Start();
102 AdvectCoeffs(nConvectiveFields, fields, advVel, inarray, tmp, time, pFwd,
103 pBwd);
104 timer.Stop();
105 timer.AccumulateRegion("AdvectCoeffs", 2);
106
107 timer.Start();
108 for (int i = 0; i < nConvectiveFields; ++i)
109 {
110 fields[i]->BwdTrans(tmp[i], outarray[i]);
111 }
112 timer.Stop();
113 timer.AccumulateRegion("AdvWeakDG:_BwdTrans", 10);
114 timer1.Stop();
115 timer1.AccumulateRegion("AdvWeakDG:All", 10);
116}
117
118/**
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
147 timer.Start();
148 AdvectVolumeFlux(inarray, fluxvector);
149 timer.Stop();
150 timer.AccumulateRegion("AdvWeakDG:_fluxVector", 3);
151
152 timer.Start();
153 // Get the advection part (without numerical flux)
154 for (int i = 0; i < nConvectiveFields; ++i)
155 {
156 Vmath::Fill(outarray[i].size(), 0.0, outarray[i], 1);
157 fields[i]->IProductWRTDerivBase(fluxvector[i], outarray[i]);
158 }
159 timer.Stop();
160 timer.AccumulateRegion("AdvWeakDG:_IProductWRTDerivBase", 10);
161
162 Array<OneD, Array<OneD, NekDouble>> numflux{size_t(nConvectiveFields)};
163
164 for (int i = 0; i < nConvectiveFields; ++i)
165 {
166 numflux[i] = Array<OneD, NekDouble>{nTracePointsTot, 0.0};
167 }
168
169 AdvectTraceFlux(nConvectiveFields, fields, advVel, inarray, numflux, time,
170 pFwd, pBwd);
171
172 // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
173 for (int i = 0; i < nConvectiveFields; ++i)
174 {
175 Vmath::Neg(nCoeffs, outarray[i], 1);
176
177 timer.Start();
178 fields[i]->AddTraceIntegral(numflux[i], outarray[i]);
179 timer.Stop();
180 timer.AccumulateRegion("AdvWeakDG:_AddTraceIntegral", 10);
181
182 timer.Start();
183 fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
184 timer.Stop();
185 timer.AccumulateRegion("AdvWeakDG:_MultiplyByElmtInvMass", 10);
186 }
187}
188
189/**
190 *
191 */
193 const int nConvectiveFields,
195 const Array<OneD, Array<OneD, NekDouble>> &advVel,
196 const Array<OneD, Array<OneD, NekDouble>> &inarray,
197 Array<OneD, Array<OneD, NekDouble>> &TraceFlux, const NekDouble &time,
198 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
199 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
200{
201 boost::ignore_unused(advVel, time);
202 int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
203
204 ASSERTL1(m_riemann, "Riemann solver must be provided for AdvectionWeakDG.");
205
206 // Store forwards/backwards space along trace space
207 Array<OneD, Array<OneD, NekDouble>> Fwd(nConvectiveFields);
208 Array<OneD, Array<OneD, NekDouble>> Bwd(nConvectiveFields);
209
211 {
212 for (int i = 0; i < nConvectiveFields; ++i)
213 {
214 Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
215 Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
216 fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
217 }
218 }
219 else
220 {
221 for (int i = 0; i < nConvectiveFields; ++i)
222 {
223 Fwd[i] = pFwd[i];
224 Bwd[i] = pBwd[i];
225 }
226 }
227
229 timer.Start();
230 m_riemann->Solve(m_spaceDim, Fwd, Bwd, TraceFlux);
231 timer.Stop();
232 timer.AccumulateRegion("AdvWeakDG:_Riemann", 10);
233}
234
235} // end of namespace SolverUtils
236} // 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:72
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
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:172
SOLVER_UTILS_EXPORT void AdvectVolumeFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &VolumeFlux)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Initialise AdvectionWeakDG objects and store them before starting the time-stepping.
SOLVER_UTILS_EXPORT void 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)
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).
SOLVER_UTILS_EXPORT void 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)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43