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