Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
Nektar::FilterCellHistoryPoints Class Reference

#include <FilterCellHistoryPoints.h>

Inheritance diagram for Nektar::FilterCellHistoryPoints:
[legend]

Public Member Functions

 FilterCellHistoryPoints (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 
 ~FilterCellHistoryPoints () override
 
void SetCellModel (CellModelSharedPtr &pCellModel)
 
- Public Member Functions inherited from Nektar::SolverUtils::FilterHistoryPoints
SOLVER_UTILS_EXPORT FilterHistoryPoints (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterHistoryPoints () override
 
- Public Member Functions inherited from Nektar::SolverUtils::Filter
SOLVER_UTILS_EXPORT Filter (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
 
virtual SOLVER_UTILS_EXPORT ~Filter ()
 
SOLVER_UTILS_EXPORT void Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT void Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT void Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT bool IsTimeDependent ()
 

Static Public Member Functions

static SolverUtils::FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::SolverUtils::FilterHistoryPoints
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. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::FilterHistoryPoints
static std::string className
 Name of the class. More...
 

Protected Member Functions

void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
virtual void v_WriteData (const int &rank, const Array< OneD, NekDouble > &data, const int &numFields, const NekDouble &time)
 
- Protected Member Functions inherited from Nektar::SolverUtils::FilterHistoryPoints
SOLVER_UTILS_EXPORT void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT bool v_IsTimeDependent () override
 
bool GetPoint (Array< OneD, NekDouble > gloCoord, int I)
 
SOLVER_UTILS_EXPORT void v_WriteData (const int &rank, const Array< OneD, NekDouble > &data, const int &numFields, const NekDouble &time)
 
virtual void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual bool v_IsTimeDependent ()=0
 

Protected Attributes

CellModelSharedPtr m_cell
 
- Protected Attributes inherited from Nektar::SolverUtils::FilterHistoryPoints
Array< OneD, Array< OneD, const NekDouble > > m_historyPoints
 
size_t m_historyPointsSize = 0
 
unsigned int m_index = 0
 
unsigned int m_outputFrequency
 
int m_outputPlane
 plane to take history point from if using a homogeneous1D expansion More...
 
std::vector< int > m_planeIDs
 
bool m_isHomogeneous1D
 
bool m_waveSpace
 
std::string m_outputFile
 
std::ofstream m_outputStream
 
std::stringstream m_historyPointStream
 
std::list< std::tuple< Array< OneD, const NekDouble >, Array< OneD, const NekDouble >, int, int > > m_historyList
 
std::map< LibUtilities::PtsType, Array< OneD, NekDouble > > m_pointDatMap
 
std::map< LibUtilities::PtsType, Array< OneD, int > > m_pointNumMap
 
unsigned int m_outputIndex = 0
 
bool m_outputOneFile
 
bool m_adaptive
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Friends

class MemoryManager< FilterCellHistoryPoints >
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Filter
typedef std::map< std::string, std::string > ParamMap
 

Detailed Description

Definition at line 44 of file FilterCellHistoryPoints.h.

Constructor & Destructor Documentation

◆ FilterCellHistoryPoints()

Nektar::FilterCellHistoryPoints::FilterCellHistoryPoints ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< SolverUtils::EquationSystem > &  pEquation,
const ParamMap pParams 
)

Definition at line 52 of file FilterCellHistoryPoints.cpp.

56 : FilterHistoryPoints(pSession, pEquation, pParams)
57{
58}
SOLVER_UTILS_EXPORT FilterHistoryPoints(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)

◆ ~FilterCellHistoryPoints()

Nektar::FilterCellHistoryPoints::~FilterCellHistoryPoints ( )
override

Definition at line 63 of file FilterCellHistoryPoints.cpp.

64{
65}

Member Function Documentation

◆ create()

static SolverUtils::FilterSharedPtr Nektar::FilterCellHistoryPoints::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< SolverUtils::EquationSystem > &  pEquation,
const ParamMap pParams 
)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file FilterCellHistoryPoints.h.

54 {
57 pSession, pEquation, pParams);
58 return p;
59 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Filter > FilterSharedPtr
A shared pointer to a Driver object.
Definition: Filter.h:51

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ SetCellModel()

void Nektar::FilterCellHistoryPoints::SetCellModel ( CellModelSharedPtr pCellModel)
inline

Definition at line 70 of file FilterCellHistoryPoints.h.

71 {
72 m_cell = pCellModel;
73 }

References m_cell.

◆ v_Update()

void Nektar::FilterCellHistoryPoints::v_Update ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterHistoryPoints.

Definition at line 70 of file FilterCellHistoryPoints.cpp.

73{
74 // Only output every m_outputFrequency.
75 if ((m_index++) % m_outputFrequency)
76 {
77 return;
78 }
79
80 const size_t numPoints = m_historyPoints.size();
81 const size_t numFields = m_cell->GetNumCellVariables();
82 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
83 Array<OneD, NekDouble> data(numPoints * numFields, 0.0);
84 Array<OneD, NekDouble> gloCoord(3, 0.0);
85 Array<OneD, NekDouble> physvals;
86 Array<OneD, NekDouble> locCoord;
87
88 // Pull out data values field by field
89 for (size_t j = 0; j < numFields; ++j)
90 {
92 {
93 // Number of points per plane
94 const int nppp = pFields[0]->GetPlane(0)->GetTotPoints();
95
96 for (auto &x : m_historyList)
97 {
98 locCoord = std::get<1>(x);
99 const int indx = std::get<2>(x);
100 const int expId = std::get<3>(x);
101
102 physvals = m_cell->GetCellSolution(j) + m_outputPlane * nppp +
103 pFields[j]->GetPhys_Offset(expId);
104
105 // interpolate point
106 data[indx * numFields + j] =
107 pFields[0]->GetExp(expId)->StdPhysEvaluate(locCoord,
108 physvals);
109 }
110 }
111 else
112 {
113 for (auto &x : m_historyList)
114 {
115 locCoord = std::get<1>(x);
116 const int indx = std::get<2>(x);
117 const int expId = std::get<3>(x);
118
119 physvals = m_cell->GetCellSolution(j) +
120 pFields[0]->GetPhys_Offset(expId);
121
122 // interpolate point
123 data[indx * numFields + j] =
124 pFields[0]->GetExp(expId)->StdPhysEvaluate(locCoord,
125 physvals);
126 }
127 }
128 }
129
130 // Exchange history data
131 // This could be improved to reduce communication but works for now
132 vComm->AllReduce(data, LibUtilities::ReduceSum);
133 v_WriteData(vComm->GetRank(), data, numFields, time);
134}
virtual void v_WriteData(const int &rank, const Array< OneD, NekDouble > &data, const int &numFields, const NekDouble &time)
Array< OneD, Array< OneD, const NekDouble > > m_historyPoints
int m_outputPlane
plane to take history point from if using a homogeneous1D expansion
std::list< std::tuple< Array< OneD, const NekDouble >, Array< OneD, const NekDouble >, int, int > > m_historyList
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55

References m_cell, Nektar::SolverUtils::FilterHistoryPoints::m_historyList, Nektar::SolverUtils::FilterHistoryPoints::m_historyPoints, Nektar::SolverUtils::FilterHistoryPoints::m_index, Nektar::SolverUtils::FilterHistoryPoints::m_isHomogeneous1D, Nektar::SolverUtils::FilterHistoryPoints::m_outputFrequency, Nektar::SolverUtils::FilterHistoryPoints::m_outputPlane, Nektar::LibUtilities::ReduceSum, and v_WriteData().

◆ v_WriteData()

void Nektar::FilterCellHistoryPoints::v_WriteData ( const int &  rank,
const Array< OneD, NekDouble > &  data,
const int &  numFields,
const NekDouble time 
)
protectedvirtual

Definition at line 136 of file FilterCellHistoryPoints.cpp.

140{
141 // Only the root process writes out history data
142 if (rank > 0)
143 {
144 return;
145 }
146
147 Array<OneD, NekDouble> gloCoord(3, 0.0);
148 if (!m_outputOneFile || m_index == 1)
149 {
150 std::stringstream vOutputFilename;
151 if (m_outputOneFile)
152 {
153 vOutputFilename << m_outputFile << ".his";
154 }
155 else
156 {
157 vOutputFilename << m_outputFile << "_" << m_outputIndex << ".his";
158 }
160 if (m_adaptive)
161 {
162 m_outputStream.open(vOutputFilename.str().c_str(), ofstream::app);
163 }
164 else
165 {
166 m_outputStream.open(vOutputFilename.str().c_str());
167 }
168 // m_outputStream << "# History data for variables (:";
169
170 // for (int i = 0; i < numFields; ++i)
171 //{
172 // m_outputStream << m_session->GetVariable(i) <<",";
173 // }
174
176 {
177 m_outputStream << ") at points:" << endl;
178 }
179 else
180 {
181 m_outputStream << ") at points:" << endl;
182 }
183
184 for (int i = 0; i < m_historyPoints.size(); ++i)
185 {
186 gloCoord = m_historyPoints[i];
187
188 m_outputStream << "# " << boost::format("%6.0f") % i;
189 m_outputStream << " " << boost::format("%15.9e") % gloCoord[0];
190 m_outputStream << " " << boost::format("%15.9e") % gloCoord[1];
191 m_outputStream << " " << boost::format("%15.9e") % gloCoord[2];
192 m_outputStream << endl;
193 }
194
196 {
197 if (m_waveSpace)
198 {
199 m_outputStream << "# (in Wavespace)" << endl;
200 }
201 }
202 }
203
204 // Write data values point by point
205 for (int k = 0; k < m_historyPoints.size(); ++k)
206 {
207 m_outputStream << boost::format("%15.9e") % time;
208 for (int j = 0; j < numFields; ++j)
209 {
210 m_outputStream << " "
211 << boost::format("%15.9e") % data[k * numFields + j];
212 }
213 m_outputStream << endl;
214 }
215
216 if (!m_outputOneFile)
217 {
218 m_outputStream.close();
219 }
220}

References CellMLToNektar.pycml::format, Nektar::SolverUtils::FilterHistoryPoints::m_adaptive, Nektar::SolverUtils::FilterHistoryPoints::m_historyPoints, Nektar::SolverUtils::FilterHistoryPoints::m_index, Nektar::SolverUtils::FilterHistoryPoints::m_isHomogeneous1D, Nektar::SolverUtils::FilterHistoryPoints::m_outputFile, Nektar::SolverUtils::FilterHistoryPoints::m_outputIndex, Nektar::SolverUtils::FilterHistoryPoints::m_outputOneFile, Nektar::SolverUtils::FilterHistoryPoints::m_outputStream, and Nektar::SolverUtils::FilterHistoryPoints::m_waveSpace.

Referenced by v_Update().

Friends And Related Function Documentation

◆ MemoryManager< FilterCellHistoryPoints >

friend class MemoryManager< FilterCellHistoryPoints >
friend

Definition at line 1 of file FilterCellHistoryPoints.h.

Member Data Documentation

◆ className

std::string Nektar::FilterCellHistoryPoints::className
static
Initial value:
=
"CellHistoryPoints", FilterCellHistoryPoints::create)
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 62 of file FilterCellHistoryPoints.h.

◆ m_cell

CellModelSharedPtr Nektar::FilterCellHistoryPoints::m_cell
protected

Definition at line 83 of file FilterCellHistoryPoints.h.

Referenced by SetCellModel(), and v_Update().