41 #include <boost/algorithm/string.hpp>
42 #include <boost/core/ignore_unused.hpp>
57 const std::weak_ptr<EquationSystem> &pEquation,
const ParamMap &pParams)
58 :
Filter(pSession, pEquation)
63 auto it = pParams.find(
"OutputFile");
64 if (it == pParams.end())
70 ASSERTL0(it->second.length() > 0,
"Empty parameter 'OutputFile'.");
76 it = pParams.find(
"Composites");
77 ASSERTL0(it != pParams.end(),
"Missing parameter 'Composites'.");
78 ASSERTL0(it->second.length() > 0,
"Empty parameter 'Composites'.");
80 std::vector<std::string> splitComposite;
82 boost::token_compress_on);
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);
91 std::vector<unsigned int> tmpVec;
94 ASSERTL0(parseGood && !tmpVec.empty(),
95 "Unable to read composite regions index range for "
104 it = pParams.find(
"OutputPrecision");
105 if (it == pParams.end())
111 ASSERTL0(it->second.length() > 0,
"Empty parameter 'OutputPrecision'.");
112 precision = std::stoi(it->second);
116 auto equationSys =
m_equ.lock();
117 ASSERTL0(equationSys,
"Weak pointer expired");
121 m_comm = pSession->GetComm();
122 if (
m_comm->GetRank() == 0)
126 m_outFile.setf(std::ios::scientific, std::ios::floatfield);
134 std::string varName = equationSys->GetVariable(j);
135 m_outFile <<
" " + compName +
"_" + varName +
"_integral";
142 it = pParams.find(
"OutputFrequency");
143 if (it == pParams.end())
149 ASSERTL0(it->second.length() > 0,
"Empty parameter 'OutputFrequency'.");
164 auto expList = pFields[0]->GetExp();
165 auto meshGraph = pFields[0]->GetGraph();
167 for (
size_t i = 0; i < expList->size(); ++i)
169 auto exp = (*expList)[i];
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)
178 auto exp = (*trace)[i];
179 geomIdToTraceId[exp->GetGeom()->GetGlobalID()] = i;
183 auto composites = pFields[0]->GetGraph()->GetComposites();
184 size_t meshDim = pFields[0]->GetGraph()->GetMeshDimension();
189 if (composites.find(
m_compVector[i][0]) == composites.end())
194 std::vector<std::shared_ptr<SpatialDomains::Geometry>> geomVec =
197 composites[
m_compVector[i][0]]->m_geomVec[0]->GetShapeDim();
200 std::vector<size_t> compGeomIds;
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.")
208 auto compGeom = composites[compNum]->m_geomVec;
210 for (
auto &geom : compGeom)
212 compGeomIds.emplace_back(geom->GetGlobalID());
217 dim == compGeom[0]->GetShapeDim(),
218 "Differing geometry dimensions specified in FilterIntegral '" +
222 std::vector<std::pair<LocalRegions::ExpansionSharedPtr, int>>
223 tmpCompExp(compGeomIds.size());
227 if (dim == pFields[0]->GetShapeDimension())
229 for (
size_t j = 0; j < compGeomIds.size(); ++j)
231 tmpCompExp[j] = std::make_pair(
238 else if (meshDim > 1 && dim == meshDim - 1)
240 for (
size_t j = 0; j < compGeomIds.size(); ++j)
243 (*trace)[geomIdToTraceId[compGeomIds[j]]];
246 exp->GetLeftAdjacentElementExp();
247 int leftAdjElmtFace = exp->GetLeftAdjacentElementTrace();
249 tmpCompExp[j] = std::make_pair(leftAdjElmtExp, leftAdjElmtFace);
255 "FilterIntegral: Only composite dimensions equal to or "
256 "one lower than the mesh dimension are supported.")
277 if (
m_comm->GetRank() == 0)
287 phys = pFields[i]->GetPhys();
296 size_t dim = compExp[0].first->GetGeom()->GetShapeDim();
297 size_t meshDim = pFields[i]->GetGraph()->GetMeshDimension();
302 for (
auto &expPair : compExp)
305 auto exp = expPair.first;
309 size_t offset = pFields[i]->GetPhys_Offset(
311 input = exp->Integral(phys + offset);
313 else if (meshDim > 1 && dim == meshDim - 1)
316 exp->GetTracePhysVals(expPair.second, exp, phys,
321 ->GetExp(exp->GetGeom()->GetTid(expPair.second))
322 ->Integral(tracePhys);
327 "FilterIntegral: Only composite dimensions "
328 "equal to or one lower than the mesh "
329 "dimension are supported.")
333 c += fabs(sum) >= fabs(input) ? (sum - t) + input
342 if (
m_comm->GetRank() == 0)
349 if (
m_comm->GetRank() == 0)
362 boost::ignore_unused(pFields, time);
364 if (
m_comm->GetRank() == 0)
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
NekDouble Evaluate() const
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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.
LibUtilities::SessionReaderSharedPtr m_session
const std::weak_ptr< EquationSystem > m_equ
std::map< std::string, std::string > ParamMap
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
FilterFactory & GetFilterFactory()
The above copyright notice and this permission notice shall be included.