Nektar++
FilterIntegral.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File FilterIntegral.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: Outputs integrals of fields during time-stepping.
32//
33///////////////////////////////////////////////////////////////////////////////
34
40
41#include <boost/algorithm/string.hpp>
42
43namespace Nektar::SolverUtils
44{
45std::string FilterIntegral::className =
48
49/**
50 * @brief Initialise the filter and parse the session file parameters.
51 */
54 const std::weak_ptr<EquationSystem> &pEquation, const ParamMap &pParams)
55 : Filter(pSession, pEquation)
56{
57 std::string outName;
58
59 // OutputFile
60 auto it = pParams.find("OutputFile");
61 if (it == pParams.end())
62 {
63 outName = m_session->GetSessionName();
64 }
65 else
66 {
67 ASSERTL0(it->second.length() > 0, "Empty parameter 'OutputFile'.");
68 outName = it->second;
69 }
70 outName += ".int";
71
72 // Composites (to calculate integrals on)
73 it = pParams.find("Composites");
74 ASSERTL0(it != pParams.end(), "Missing parameter 'Composites'.");
75 ASSERTL0(it->second.length() > 0, "Empty parameter 'Composites'.");
76
77 std::vector<std::string> splitComposite;
78 boost::split(m_splitCompString, it->second, boost::is_any_of(","),
79 boost::token_compress_on);
80
81 for (auto &comp : m_splitCompString)
82 {
83 boost::trim(comp);
84 size_t first = comp.find_first_of('[') + 1;
85 size_t last = comp.find_last_of(']') - 1;
86 auto tmpString = comp.substr(first, last - first + 1);
87
88 std::vector<unsigned int> tmpVec;
89 bool parseGood = ParseUtils::GenerateSeqVector(tmpString, tmpVec);
90
91 ASSERTL0(parseGood && !tmpVec.empty(),
92 "Unable to read composite regions index range for "
93 "FilterIntegral: " +
94 comp);
95
96 m_compVector.emplace_back(tmpVec);
97 }
98
99 // OutputPrecision
100 size_t precision;
101 it = pParams.find("OutputPrecision");
102 if (it == pParams.end())
103 {
104 precision = 7;
105 }
106 else
107 {
108 ASSERTL0(it->second.length() > 0, "Empty parameter 'OutputPrecision'.");
109 precision = std::stoi(it->second);
110 }
111
112 // Lock equation system pointer
113 auto equationSys = m_equ.lock();
114 ASSERTL0(equationSys, "Weak pointer expired");
115
116 m_numVariables = equationSys->GetNvariables();
117
118 m_comm = pSession->GetComm();
119 if (m_comm->GetRank() == 0)
120 {
121 m_outFile.open(outName);
122 ASSERTL0(m_outFile.good(), "Unable to open: '" + outName + "'");
123 m_outFile.setf(std::ios::scientific, std::ios::floatfield);
124 m_outFile.precision(precision);
125 m_outFile << "#Time";
126
127 for (auto &compName : m_splitCompString)
128 {
129 for (size_t j = 0; j < m_numVariables; ++j)
130 {
131 std::string varName = equationSys->GetVariable(j);
132 m_outFile << " " + compName + "_" + varName + "_integral";
133 }
134 }
135 m_outFile << std::endl;
136 }
137
138 // OutputFrequency
139 it = pParams.find("OutputFrequency");
140 if (it == pParams.end())
141 {
143 }
144 else
145 {
146 ASSERTL0(it->second.length() > 0, "Empty parameter 'OutputFrequency'.");
147 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
148 m_outputFrequency = round(equ.Evaluate());
149 }
150}
151
152/**
153 * @brief Parse composite list and geometric entities to integrate over.
154 */
157 const NekDouble &time)
158{
159
160 // Create map from element ID -> expansion list ID
161 auto expList = pFields[0]->GetExp();
162 auto meshGraph = pFields[0]->GetGraph();
163
164 for (size_t i = 0; i < expList->size(); ++i)
165 {
166 auto exp = (*expList)[i];
167 m_geomElmtIdToExpId[exp->GetGeom()->GetGlobalID()] = i;
168 }
169
170 // Create a map from geom ID -> trace expansion list ID
171 std::map<size_t, size_t> geomIdToTraceId;
172 auto trace = pFields[0]->GetTrace()->GetExp();
173 for (size_t i = 0; i < trace->size(); ++i)
174 {
175 auto exp = (*trace)[i];
176 geomIdToTraceId[exp->GetGeom()->GetGlobalID()] = i;
177 }
178
179 // Get comp list dimension from first composite & element
180 auto composites = pFields[0]->GetGraph()->GetComposites();
181 size_t meshDim = pFields[0]->GetGraph()->GetMeshDimension();
182
183 for (size_t i = 0; i < m_compVector.size(); ++i)
184 {
185 // Check composite is present in the rank
186 if (composites.find(m_compVector[i][0]) == composites.end())
187 {
188 continue;
189 }
190
191 std::vector<std::shared_ptr<SpatialDomains::Geometry>> geomVec =
192 composites[m_compVector[i][0]]->m_geomVec;
193 size_t dim =
194 composites[m_compVector[i][0]]->m_geomVec[0]->GetShapeDim();
195
196 // Vector of all geometry IDs contained within the composite list
197 std::vector<size_t> compGeomIds;
198 for (auto compNum : m_compVector[i])
199 {
200 ASSERTL0(composites.find(compNum) != composites.end(),
201 "In FilterIntegral defined composite C[" +
202 std::to_string(compNum) +
203 "] does not exist in the mesh.")
204
205 auto compGeom = composites[compNum]->m_geomVec;
206
207 for (auto &geom : compGeom)
208 {
209 compGeomIds.emplace_back(geom->GetGlobalID());
210 }
211
212 // Only check first element in each comp for dimension
213 ASSERTL0(
214 dim == compGeom[0]->GetShapeDim(),
215 "Differing geometry dimensions specified in FilterIntegral '" +
216 m_splitCompString[i] + "'.");
217 }
218
219 std::vector<std::pair<LocalRegions::ExpansionSharedPtr, int>>
220 tmpCompExp(compGeomIds.size());
221
222 // If dimension of composite == dimension of mesh then we only need the
223 // expansion of the element
224 if (dim == pFields[0]->GetShapeDimension())
225 {
226 for (size_t j = 0; j < compGeomIds.size(); ++j)
227 {
228 tmpCompExp[j] = std::make_pair(
229 (*expList)[m_geomElmtIdToExpId[compGeomIds[j]]], -1);
230 }
231 }
232 // however if the dimension is less we need the expansion of the element
233 // containing the global composite geometry and the face/edge local ID
234 // within that. 3D mesh -> 2D, 2D -> 1D.
235 else if (meshDim > 1 && dim == meshDim - 1)
236 {
237 for (size_t j = 0; j < compGeomIds.size(); ++j)
238 {
240 (*trace)[geomIdToTraceId[compGeomIds[j]]];
241
242 LocalRegions::ExpansionSharedPtr leftAdjElmtExp =
243 exp->GetLeftAdjacentElementExp();
244 int leftAdjElmtFace = exp->GetLeftAdjacentElementTrace();
245
246 tmpCompExp[j] = std::make_pair(leftAdjElmtExp, leftAdjElmtFace);
247 }
248 }
249 else
250 {
252 "FilterIntegral: Only composite dimensions equal to or "
253 "one lower than the mesh dimension are supported.")
254 }
255
256 m_compExpMap[i] = tmpCompExp;
257 }
258
259 v_Update(pFields, time);
260}
261
262/**
263 * @brief Calculate integral over requested composites.
264 */
267 const NekDouble &time)
268{
269 if (m_index++ % m_outputFrequency > 0)
270 {
271 return;
272 }
273
274 if (m_comm->GetRank() == 0)
275 {
276 m_outFile << time;
277 }
278
279 for (size_t j = 0; j < m_compVector.size(); ++j)
280 {
281 for (size_t i = 0; i < m_numVariables; ++i)
282 {
284 phys = pFields[i]->GetPhys();
285
286 NekDouble sum = 0.0;
287 NekDouble c = 0.0;
288
289 // Check if composite is on the rank
290 if (m_compExpMap.find(j) != m_compExpMap.end())
291 {
292 auto compExp = m_compExpMap[j];
293 size_t dim = compExp[0].first->GetGeom()->GetShapeDim();
294 size_t meshDim = pFields[i]->GetGraph()->GetMeshDimension();
295
296 // Evaluate integral using improved Kahan–Babuška summation
297 // algorithm to reduce numerical error from adding floating
298 // points
299 for (auto &expPair : compExp)
300 {
301 NekDouble input = 0;
302 auto exp = expPair.first;
303
304 if (meshDim == dim)
305 {
306 size_t offset = pFields[i]->GetPhys_Offset(
307 m_geomElmtIdToExpId[exp->GetGeom()->GetGlobalID()]);
308 input = exp->Integral(phys + offset);
309 }
310 else if (meshDim > 1 && dim == meshDim - 1)
311 {
312 Array<OneD, NekDouble> tracePhys;
313 exp->GetTracePhysVals(expPair.second, exp, phys,
314 tracePhys);
315 input =
316 pFields[i]
317 ->GetTrace()
318 ->GetExp(exp->GetGeom()->GetTid(expPair.second))
319 ->Integral(tracePhys);
320 }
321 else
322 {
324 "FilterIntegral: Only composite dimensions "
325 "equal to or one lower than the mesh "
326 "dimension are supported.")
327 }
328
329 NekDouble t = sum + input;
330 c += fabs(sum) >= fabs(input) ? (sum - t) + input
331 : (input - t) + sum;
332 sum = t;
333 }
334 }
335
336 // Sum integral values from all ranks
337 Array<OneD, NekDouble> sumArray(1, sum + c);
338 m_comm->AllReduce(sumArray, LibUtilities::ReduceSum);
339 if (m_comm->GetRank() == 0)
340 {
341 m_outFile << " " << sumArray[0];
342 }
343 }
344 }
345
346 if (m_comm->GetRank() == 0)
347 {
348 m_outFile << std::endl;
349 }
350}
351
352/**
353 * @brief Finalise the filter.
354 */
357 &pFields,
358 [[maybe_unused]] const NekDouble &time)
359{
360 if (m_comm->GetRank() == 0)
361 {
362 m_outFile.close();
363 }
364}
365
366/**
367 * @brief This is a time-dependent filter.
368 */
370{
371 return true;
372}
373
374} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
Definition: ParseUtils.cpp:104
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
std::map< int, std::vector< std::pair< LocalRegions::ExpansionSharedPtr, int > > > m_compExpMap
Map of composite ID to vector of expansions an face/edge local ID.
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) final
Calculate integral over requested composites.
std::map< size_t, size_t > m_geomElmtIdToExpId
Mapping from geometry ID to expansion ID.
size_t m_outputFrequency
Frequency to write to output data file in timesteps.
size_t m_numVariables
Number of fields to perform integral on.
std::ofstream m_outFile
Out file.
bool v_IsTimeDependent() final
Returns true as filter depends on time.
LibUtilities::CommSharedPtr m_comm
Global communicator.
static std::string className
Name of the class.
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) final
Finalise the filter.
std::vector< std::vector< unsigned int > > m_compVector
Vector of vector of composites IDs as integers.
std::vector< std::string > m_splitCompString
Vector of composite IDs as a single string.
SOLVER_UTILS_EXPORT FilterIntegral(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
Constructs the integral filter and parses filter options, opens file.
void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) final
Parse composite list and geometric entities to integrate over.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39
double NekDouble