Nektar++
ForcingAbsorption.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ForcingAbsorption.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2016 Kilian Lackhove
10// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11// Department of Aeronautics, Imperial College London (UK), and Scientific
12// Computing and Imaging Institute, University of Utah (USA).
13//
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Absorption layer forcing
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <boost/core/ignore_unused.hpp>
37
39
41
42using namespace std;
43
44namespace Nektar
45{
46namespace SolverUtils
47{
48
51 "Absorption", ForcingAbsorption::create, "Forcing Absorption");
52
55 const std::weak_ptr<EquationSystem> &pEquation)
56 : Forcing(pSession, pEquation), m_hasRefFlow(false), m_hasRefFlowTime(false)
57{
58}
59
62 const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
63{
64 m_NumVariable = pNumForcingFields;
65 int npts = pFields[0]->GetTotPoints();
66
67 CalcAbsorption(pFields, pForce);
68
70 for (int i = 0; i < m_NumVariable; ++i)
71 {
72 m_Forcing[i] = Array<OneD, NekDouble>(npts, 0.0);
73 }
74
75 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("REFFLOW");
76 if (funcNameElmt)
77 {
78 string funcName = funcNameElmt->GetText();
79 ASSERTL0(m_session->DefinesFunction(funcName),
80 "Function '" + funcName + "' not defined.");
82 for (int i = 0; i < m_NumVariable; ++i)
83 {
84 std::string s_FieldStr = m_session->GetVariable(i);
85 ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
86 "Variable '" + s_FieldStr + "' not defined.");
87 m_Refflow[i] = Array<OneD, NekDouble>(npts, 0.0);
88 GetFunction(pFields, m_session, funcName)
89 ->Evaluate(s_FieldStr, m_Refflow[i]);
90 }
91 m_hasRefFlow = true;
92 }
93
94 funcNameElmt = pForce->FirstChildElement("REFFLOWTIME");
95 if (funcNameElmt)
96 {
97 m_funcNameTime = funcNameElmt->GetText();
98 m_hasRefFlowTime = true;
99 m_hasRefFlow = true;
101 for (int i = 0; i < m_NumVariable; ++i)
102 {
103 std::string s_FieldStr = m_session->GetVariable(i);
104 ASSERTL0(m_session->DefinesFunction(m_funcNameTime, s_FieldStr),
105 "Variable '" + s_FieldStr + "' not defined.");
106 m_Refflow[i] = Array<OneD, NekDouble>(npts, 0.0);
107 }
108 }
109}
110
113 &pFields,
114 const TiXmlElement *pForce)
115{
116 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("COEFF");
117 ASSERTL0(funcNameElmt,
118 "Requires COEFF tag, specifying function "
119 "name which prescribes absorption layer coefficient.");
120 string funcName = funcNameElmt->GetText();
121 ASSERTL0(m_session->DefinesFunction(funcName),
122 "Function '" + funcName + "' not defined.");
123
124 int npts = pFields[0]->GetTotPoints();
125
127 for (int i = 0; i < m_NumVariable; ++i)
128 {
129 m_Absorption[i] = Array<OneD, NekDouble>(npts, 0.0);
130 }
131
132 funcNameElmt = pForce->FirstChildElement("BOUNDARYREGIONS");
133 if (funcNameElmt)
134 {
135 ASSERTL0(
136 ParseUtils::GenerateVector(funcNameElmt->GetText(), m_bRegions),
137 "Unable to process list of BOUNDARYREGIONS in Absorption "
138 "Forcing: " +
139 std::string(funcNameElmt->GetText()));
140
141 // alter m_bRegions so that it contains the boundaryRegions of this rank
142 std::vector<unsigned int> localBRegions;
144 pFields[0]->GetGraph());
146 bcs.GetBoundaryRegions();
147 SpatialDomains::BoundaryRegionCollection::iterator it1;
148 int n = 0;
149 for (it1 = regions.begin(); it1 != regions.end(); ++it1)
150 {
151 if (std::find(m_bRegions.begin(), m_bRegions.end(), it1->first) !=
152 m_bRegions.end())
153 {
154 localBRegions.push_back(n);
155 }
156 n++;
157 }
158 m_bRegions = localBRegions;
159
160 if (m_bRegions.size() == 0)
161 {
162 return;
163 }
164
165 std::vector<Array<OneD, const NekDouble>> points;
166
168 for (int i = 0; i < 3; i++)
169 {
170 x[i] = Array<OneD, NekDouble>(npts, 0.0);
171 }
172 pFields[0]->GetCoords(x[0], x[1], x[2]);
173 for (int i = 0; i < 3; i++)
174 {
175 points.push_back(x[i]);
176 }
177
178 Array<OneD, NekDouble> t(npts, 0.0);
179 points.push_back(t);
180
181 Array<OneD, NekDouble> r(npts, 0.0);
182 std::vector<unsigned int>::iterator it;
183 std::vector<BPointPair> inPoints;
185 for (it = m_bRegions.begin(); it != m_bRegions.end(); ++it)
186 {
187 int bpts = pFields[0]->GetBndCondExpansions()[*it]->GetNpoints();
188 for (int i = 0; i < 3; i++)
189 {
190 b[i] = Array<OneD, NekDouble>(bpts, 0.0);
191 }
192 pFields[0]->GetBndCondExpansions()[*it]->GetCoords(b[0], b[1],
193 b[2]);
194 for (int i = 0;
195 i < pFields[0]->GetBndCondExpansions()[*it]->GetNpoints(); ++i)
196 {
197 inPoints.push_back(
198 BPointPair(BPoint(b[0][i], b[1][i], b[2][i]), i));
199 }
200 }
202 m_rtree->insert(inPoints.begin(), inPoints.end());
203
204 for (int i = 0; i < npts; ++i)
205 {
206 std::vector<BPointPair> result;
207 BPoint sPoint(x[0][i], x[1][i], x[2][i]);
208 m_rtree->query(bgi::nearest(sPoint, 1), std::back_inserter(result));
209 r[i] = bg::distance(sPoint, result[0].first);
210 }
211 points.push_back(r);
212
213 std::string s_FieldStr;
214 for (int i = 0; i < m_NumVariable; ++i)
215 {
216 s_FieldStr = m_session->GetVariable(i);
217 ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
218 "Variable '" + s_FieldStr + "' not defined.");
219
221 m_session->GetFunction(funcName, s_FieldStr);
222 ASSERTL0(ffunc->GetVlist() == "x y z t r",
223 "EVARS of " + funcName + " must be 'r'");
224
225 ffunc->Evaluate(points, m_Absorption[i]);
226 }
227 }
228 else
229 {
230 for (int i = 0; i < m_NumVariable; ++i)
231 {
232 std::string s_FieldStr = m_session->GetVariable(i);
233 GetFunction(pFields, m_session, funcName)
234 ->Evaluate(s_FieldStr, m_Absorption[i]);
235 }
236 }
237}
238
241 const Array<OneD, Array<OneD, NekDouble>> &inarray,
242 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time)
243{
244 int nq = m_Forcing[0].size();
245 CalculateForcing(fields, inarray, time);
246 for (int i = 0; i < m_NumVariable; ++i)
247 {
248 Vmath::Vadd(nq, m_Forcing[i], 1, outarray[i], 1, outarray[i], 1);
249 }
250}
251
254 const Array<OneD, Array<OneD, NekDouble>> &inarray,
255 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time)
256{
257 int ncoeff = outarray[m_NumVariable - 1].size();
258 Array<OneD, NekDouble> tmp(ncoeff, 0.0);
259 CalculateForcing(fields, inarray, time);
260
261 for (int i = 0; i < m_NumVariable; ++i)
262 {
263 fields[i]->FwdTrans(m_Forcing[i], tmp);
264 Vmath::Vadd(ncoeff, tmp, 1, outarray[i], 1, outarray[i], 1);
265 }
266}
267
270 const Array<OneD, Array<OneD, NekDouble>> &inarray, const NekDouble &time)
271{
272 boost::ignore_unused(fields);
273 int nq = m_Forcing[0].size();
274
275 std::string s_FieldStr;
276 Array<OneD, NekDouble> TimeScale(1);
278
279 if (m_hasRefFlow)
280 {
281 for (int i = 0; i < m_NumVariable; i++)
282 {
283 RefflowScaled[i] = Array<OneD, NekDouble>(nq);
285 {
286 s_FieldStr = m_session->GetVariable(i);
287
288 std::string s_FieldStr = m_session->GetVariable(i);
290 ->Evaluate(s_FieldStr, m_Refflow[i], time);
291 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
292 }
293 else
294 {
295 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
296 }
297
298 Vmath::Vsub(nq, inarray[i], 1, RefflowScaled[i], 1, m_Forcing[i],
299 1);
300 Vmath::Vmul(nq, m_Absorption[i], 1, m_Forcing[i], 1, m_Forcing[i],
301 1);
302 }
303 }
304 else
305 {
306 for (int i = 0; i < m_NumVariable; i++)
307 {
308 Vmath::Vmul(nq, m_Absorption[i], 1, inarray[i], 1, m_Forcing[i], 1);
309 }
310 }
311}
312
313} // namespace SolverUtils
314} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:131
bg::model::point< NekDouble, 3, bg::cs::cartesian > BPoint
virtual SOLVER_UTILS_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
virtual SOLVER_UTILS_EXPORT void v_ApplyCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
ForcingAbsorption(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
static SOLVER_UTILS_EXPORT ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
virtual SOLVER_UTILS_EXPORT void v_Apply(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
Array< OneD, Array< OneD, NekDouble > > m_Refflow
std::pair< BPoint, unsigned int > BPointPair
void CalcAbsorption(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
void CalculateForcing(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const NekDouble &time)
static std::string className
Name of the class.
std::vector< unsigned int > m_bRegions
Array< OneD, Array< OneD, NekDouble > > m_Absorption
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:73
int m_NumVariable
Number of variables.
Definition: Forcing.h:123
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:121
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, std::string pName, bool pCache=false)
Get a SessionFunction by name.
Definition: Forcing.cpp:194
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:117
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:234
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:44
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:453
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
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
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414