41#include <boost/algorithm/string.hpp>
54 const std::shared_ptr<EquationSystem> &pEquation,
const ParamMap &pParams)
55 :
Filter(pSession, pEquation)
58 std::string ext =
".int";
62 auto it = pParams.find(
"Composites");
63 ASSERTL0(it != pParams.end(),
"Missing parameter 'Composites'.");
64 ASSERTL0(it->second.length() > 0,
"Empty parameter 'Composites'.");
66 std::vector<std::string> splitComposite;
68 boost::token_compress_on);
73 size_t first = comp.find_first_of(
'[') + 1;
74 size_t last = comp.find_last_of(
']') - 1;
75 auto tmpString = comp.substr(first, last - first + 1);
77 std::vector<unsigned int> tmpVec;
80 ASSERTL0(parseGood && !tmpVec.empty(),
81 "Unable to read composite regions index range for "
90 it = pParams.find(
"OutputPrecision");
91 if (it == pParams.end())
97 ASSERTL0(it->second.length() > 0,
"Empty parameter 'OutputPrecision'.");
98 precision = std::stoi(it->second);
102 auto equationSys =
m_equ.lock();
103 ASSERTL0(equationSys,
"Weak pointer expired");
107 m_comm = pSession->GetComm();
108 if (
m_comm->GetRank() == 0)
112 m_outFile.setf(std::ios::scientific, std::ios::floatfield);
120 std::string varName = equationSys->GetVariable(j);
121 m_outFile <<
" " + compName +
"_" + varName +
"_integral";
128 it = pParams.find(
"OutputFrequency");
129 if (it == pParams.end())
135 ASSERTL0(it->second.length() > 0,
"Empty parameter 'OutputFrequency'.");
150 auto expList = pFields[0]->GetExp();
151 auto meshGraph = pFields[0]->GetGraph();
153 for (
size_t i = 0; i < expList->size(); ++i)
155 auto exp = (*expList)[i];
160 std::map<size_t, size_t> geomIdToTraceId;
161 auto trace = pFields[0]->GetTrace()->GetExp();
162 for (
size_t i = 0; i < trace->size(); ++i)
164 auto exp = (*trace)[i];
165 geomIdToTraceId[exp->GetGeom()->GetGlobalID()] = i;
169 auto composites = pFields[0]->GetGraph()->GetComposites();
170 size_t meshDim = pFields[0]->GetGraph()->GetMeshDimension();
175 if (composites.find(
m_compVector[i][0]) == composites.end())
180 std::vector<std::shared_ptr<SpatialDomains::Geometry>> geomVec =
183 composites[
m_compVector[i][0]]->m_geomVec[0]->GetShapeDim();
186 std::vector<size_t> compGeomIds;
189 ASSERTL0(composites.find(compNum) != composites.end(),
190 "In FilterIntegral defined composite C[" +
191 std::to_string(compNum) +
192 "] does not exist in the mesh.")
194 auto compGeom = composites[compNum]->m_geomVec;
196 for (
auto &geom : compGeom)
198 compGeomIds.emplace_back(geom->GetGlobalID());
203 dim == compGeom[0]->GetShapeDim(),
204 "Differing geometry dimensions specified in FilterIntegral '" +
208 std::vector<std::pair<LocalRegions::ExpansionSharedPtr, int>>
209 tmpCompExp(compGeomIds.size());
213 if (dim == pFields[0]->GetShapeDimension())
215 for (
size_t j = 0; j < compGeomIds.size(); ++j)
217 tmpCompExp[j] = std::make_pair(
226 else if (meshDim == 3 && dim == 2)
228 for (
size_t j = 0; j < compGeomIds.size(); ++j)
231 (*trace)[geomIdToTraceId[compGeomIds[j]]];
233 std::dynamic_pointer_cast<LocalRegions::Expansion2D>(exp);
236 std::dynamic_pointer_cast<LocalRegions::Expansion>(
237 exp2D->GetLeftAdjacentElementExp());
238 int leftAdjElmtFace = exp2D->GetLeftAdjacentElementTrace();
240 tmpCompExp[j] = std::make_pair(leftAdjElmtExp, leftAdjElmtFace);
243 else if (meshDim == 2 && dim == 1)
245 for (
size_t j = 0; j < compGeomIds.size(); ++j)
248 (*trace)[geomIdToTraceId[compGeomIds[j]]];
250 std::dynamic_pointer_cast<LocalRegions::Expansion1D>(exp);
253 std::dynamic_pointer_cast<LocalRegions::Expansion>(
254 exp1D->GetLeftAdjacentElementExp());
255 int leftAdjElmtEdge = exp1D->GetLeftAdjacentElementTrace();
257 tmpCompExp[j] = std::make_pair(leftAdjElmtExp, leftAdjElmtEdge);
263 "FilterIntegral: Only composite dimensions equal to or "
264 "one lower than the mesh dimension are supported.")
287 if (
m_comm->GetRank() == 0)
297 phys = pFields[i]->GetPhys();
306 size_t dim = compExp[0].first->GetGeom()->GetShapeDim();
307 size_t meshDim = pFields[i]->GetGraph()->GetMeshDimension();
312 for (
auto &expPair : compExp)
315 auto exp = expPair.first;
319 size_t offset = pFields[i]->GetPhys_Offset(
321 input = exp->Integral(phys + offset);
323 else if (meshDim == 3 && dim == 2)
326 exp->GetTracePhysVals(expPair.second, exp, phys,
331 ->GetExp(exp->GetGeom()->GetTid(expPair.second))
332 ->Integral(facePhys);
334 else if (meshDim == 2 && dim == 1)
337 exp->GetTracePhysVals(expPair.second, exp, phys,
342 ->GetExp(exp->GetGeom()->GetTid(expPair.second))
343 ->Integral(edgePhys);
348 "FilterIntegral: Only composite dimensions "
349 "equal to or one lower than the mesh "
350 "dimension are supported.")
354 c += fabs(sum) >= fabs(input) ? (sum - t) + input
363 if (
m_comm->GetRank() == 0)
370 if (
m_comm->GetRank() == 0)
384 if (
m_comm->GetRank() == 0)
#define ASSERTL0(condition, msg)
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.
SOLVER_UTILS_EXPORT std::string SetupOutput(const std::string ext, const ParamMap &pParams)
LibUtilities::SessionReaderSharedPtr m_session
const std::weak_ptr< EquationSystem > m_equ
std::map< std::string, std::string > ParamMap
bool m_updateOnInitialise
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.
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.
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
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
This is a time-dependent filter.
std::vector< std::vector< unsigned int > > m_compVector
Vector of vector of composites IDs as integers.
SOLVER_UTILS_EXPORT FilterIntegral(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
Constructs the integral filter and parses filter options, opens file.
std::vector< std::string > m_splitCompString
Vector of composite IDs as a single string.
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
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
FilterFactory & GetFilterFactory()