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

#include <FilterHistoryPoints.h>

Inheritance diagram for Nektar::SolverUtils::FilterHistoryPoints:
[legend]

Public Member Functions

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 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...
 

Protected Member Functions

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

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< FilterHistoryPoints >
 

Additional Inherited Members

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

Detailed Description

Definition at line 45 of file FilterHistoryPoints.h.

Constructor & Destructor Documentation

◆ FilterHistoryPoints()

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

Definition at line 53 of file FilterHistoryPoints.cpp.

56 : Filter(pSession, pEquation)
57{
58 // OutputFile
59 auto it = pParams.find("OutputFile");
60 if (it == pParams.end())
61 {
62 m_outputFile = m_session->GetSessionName();
63 }
64 else
65 {
66 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
67 m_outputFile = it->second;
68 }
69 if (m_outputFile.length() >= 4 &&
70 m_outputFile.substr(m_outputFile.length() - 4) == ".his")
71 {
72 m_outputFile = m_outputFile.substr(0, m_outputFile.length() - 4);
73 }
74
75 // OutputFrequency
76 it = pParams.find("OutputFrequency");
77 if (it == pParams.end())
78 {
80 }
81 else
82 {
83 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
84 m_outputFrequency = round(equ.Evaluate());
85 }
86
87 // output all data into one file
88 it = pParams.find("OutputOneFile");
89 if (it == pParams.end())
90 {
91 m_outputOneFile = true;
92 }
93 else
94 {
95 std::string sOption = it->second.c_str();
96 m_outputOneFile = (boost::iequals(sOption, "true")) ||
97 (boost::iequals(sOption, "yes"));
98 }
99
100 // OutputPlane
101 m_session->MatchSolverInfo("Homogeneous", "1D", m_isHomogeneous1D, false);
103 {
104 it = pParams.find("OutputPlane");
105 if (it == pParams.end())
106 {
107 m_outputPlane = -1;
108 }
109 else
110 {
111 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
112 m_outputPlane = round(equ.Evaluate());
113 }
114
115 it = pParams.find("WaveSpace");
116 if (it == pParams.end())
117 {
118 m_waveSpace = false;
119 }
120 else
121 {
122 std::string sOption = it->second.c_str();
123 m_waveSpace = (boost::iequals(sOption, "true")) ||
124 (boost::iequals(sOption, "yes"));
125 }
126 }
127
128 // Points
129 if (pParams.end() != (it = pParams.find("Points")))
130 {
131 m_pointNumMap[LibUtilities::ePtsFile] = Array<OneD, int>(1, 0);
132 m_historyPointStream.str(it->second);
133 }
134 else if (pParams.end() != (it = pParams.find("line")))
135 {
136 vector<NekDouble> values;
137 ASSERTL0(ParseUtils::GenerateVector(it->second, values),
138 "Failed to interpret line string");
139
140 ASSERTL0(values.size() > 2, "line string should contain 2*Dim+1 values "
141 "N,x0,y0,z0,x1,y1,z1");
142
143 double tmp;
144 ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
145 ASSERTL0(values[0] > 1, "N is not a valid number");
146
147 int dim = (values.size() - 1) / 2;
148 int npts = values[0];
149
150 Array<OneD, int> num(1, npts);
151 Array<OneD, NekDouble> delta(6, 0.);
152 for (int i = 0; i < dim; ++i)
153 {
154 delta[i + 0] = values[i + 1];
155 delta[i + 3] = (values[dim + i + 1] - values[i + 1]) / (npts - 1);
156 }
159 }
160 else if (pParams.end() != (it = pParams.find("plane")))
161 {
162 vector<NekDouble> values;
163 ASSERTL0(ParseUtils::GenerateVector(it->second, values),
164 "Failed to interpret plane string");
165
166 ASSERTL0(values.size() > 9,
167 "plane string should contain 4 Dim+2 values "
168 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
169
170 double tmp;
171 ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
172 ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
173
174 ASSERTL0(values[0] > 1, "N1 is not a valid number");
175 ASSERTL0(values[1] > 1, "N2 is not a valid number");
176
177 int dim = (values.size() - 2) / 4;
178
179 Array<OneD, int> npts(3);
180 npts[0] = values[0];
181 npts[1] = values[0] * values[1];
182 npts[2] = values[1];
183
184 Array<OneD, NekDouble> delta(12, 0.);
185 for (int i = 0; i < dim; ++i)
186 {
187 delta[i + 0] = values[2 + i];
188 delta[i + 3] = values[2 + 3 * dim + i];
189 delta[i + 6] = (values[2 + 1 * dim + i] - values[2 + 0 * dim + i]) /
190 (values[0] - 1);
191 delta[i + 9] = (values[2 + 2 * dim + i] - values[2 + 3 * dim + i]) /
192 (values[0] - 1);
193 }
196 }
197 else if (pParams.end() != (it = pParams.find("box")))
198 {
199 vector<NekDouble> values;
200 ASSERTL0(ParseUtils::GenerateVector(it->second, values),
201 "Failed to interpret box string");
202
203 ASSERTL0(values.size() == 9, "box string should contain 9 values "
204 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
205
206 int dim = 3;
207 Array<OneD, int> npts(3);
208 npts[0] = values[0];
209 npts[1] = values[0] * values[1];
210 npts[2] = values[0] * values[1] * values[2];
211
212 Array<OneD, NekDouble> delta(6, 0.);
213 for (int i = 0; i < dim; ++i)
214 {
215 delta[i + 0] = values[3 + 2 * i];
216 delta[i + 3] =
217 (values[4 + 2 * i] - values[3 + 2 * i]) / (values[i] - 1);
218 }
221 }
222 else
223 {
224 ASSERTL0(false, "Missing parameter 'Points'.");
225 }
226}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
std::map< LibUtilities::PtsType, Array< OneD, int > > m_pointNumMap
int m_outputPlane
plane to take history point from if using a homogeneous1D expansion
std::map< LibUtilities::PtsType, Array< OneD, NekDouble > > m_pointDatMap
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
SOLVER_UTILS_EXPORT Filter(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Definition: Filter.cpp:45

References ASSERTL0, Nektar::LibUtilities::ePtsBox, Nektar::LibUtilities::ePtsFile, Nektar::LibUtilities::ePtsLine, Nektar::LibUtilities::ePtsPlane, Nektar::LibUtilities::Equation::Evaluate(), Nektar::ParseUtils::GenerateVector(), m_historyPointStream, m_isHomogeneous1D, m_outputFile, m_outputFrequency, m_outputOneFile, m_outputPlane, m_pointDatMap, m_pointNumMap, Nektar::SolverUtils::Filter::m_session, and m_waveSpace.

◆ ~FilterHistoryPoints()

Nektar::SolverUtils::FilterHistoryPoints::~FilterHistoryPoints ( )
override

Definition at line 231 of file FilterHistoryPoints.cpp.

232{
233}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 51 of file FilterHistoryPoints.h.

55 {
58 pSession, pEquation, pParams);
59 return p;
60 }
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.

◆ GetPoint()

bool Nektar::SolverUtils::FilterHistoryPoints::GetPoint ( Array< OneD, NekDouble gloCoord,
int  I 
)
protected

Definition at line 235 of file FilterHistoryPoints.cpp.

236{
238 {
239 m_historyPointStream >> gloCoord[0] >> gloCoord[1] >> gloCoord[2];
240 return !m_historyPointStream.fail();
241 }
243 {
245 {
246 return false;
247 }
248 Array<OneD, NekDouble> values = m_pointDatMap[LibUtilities::ePtsLine];
249 Array<OneD, NekDouble> delta =
251 for (int n = 0; n < 3; ++n)
252 {
253 gloCoord[n] = values[n] + I * delta[n];
254 }
255 return true;
256 }
258 {
260 {
261 return false;
262 }
263 Array<OneD, NekDouble> values = m_pointDatMap[LibUtilities::ePtsPlane];
264 Array<OneD, NekDouble> delta =
266 Array<OneD, int> i(2);
268 i[0] = I - i[1] * m_pointNumMap[LibUtilities::ePtsPlane][0];
269 double n1 = -1. + (NekDouble)m_pointNumMap[LibUtilities::ePtsPlane][2];
270 for (int n = 0; n < 3; ++n)
271 {
272 gloCoord[n] = (values[n] + i[0] * delta[n]) * (1.0 - i[1] / n1) +
273 (values[3 + n] + i[0] * delta[3 + n]) * (i[1] / n1);
274 }
275 return true;
276 }
277 else if (m_pointNumMap.count(LibUtilities::ePtsBox))
278 {
280 {
281 return false;
282 }
283 Array<OneD, NekDouble> values = m_pointDatMap[LibUtilities::ePtsBox];
284 Array<OneD, NekDouble> delta = m_pointDatMap[LibUtilities::ePtsBox] + 3;
285 Array<OneD, int> i(3);
287 I -= i[2] * m_pointNumMap[LibUtilities::ePtsBox][1];
289 i[0] = I - i[1] * m_pointNumMap[LibUtilities::ePtsBox][0];
290 for (int n = 0; n < 3; ++n)
291 {
292 gloCoord[n] = values[n] + i[n] * delta[n];
293 }
294 return true;
295 }
296 return false;
297}
double NekDouble

References Nektar::LibUtilities::ePtsBox, Nektar::LibUtilities::ePtsFile, Nektar::LibUtilities::ePtsLine, Nektar::LibUtilities::ePtsPlane, m_historyPointStream, m_pointDatMap, and m_pointNumMap.

Referenced by v_Initialise().

◆ v_Finalise()

void Nektar::SolverUtils::FilterHistoryPoints::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 685 of file FilterHistoryPoints.cpp.

688{
689 if (pFields[0]->GetComm()->GetRank() == 0 && m_outputOneFile)
690 {
691 m_outputStream.close();
692 }
693}

References m_outputOneFile, and m_outputStream.

◆ v_Initialise()

void Nektar::SolverUtils::FilterHistoryPoints::v_Initialise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 302 of file FilterHistoryPoints.cpp.

305{
306 ASSERTL0(!m_historyPointStream.fail(), "No history points in stream.");
307
308 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
309
310 // All processes read all the history points
311 Array<OneD, NekDouble> gloCoord(3, 0.0);
312 std::vector<Array<OneD, const NekDouble>> GloCoords;
313 const int coordim = pFields[0]->GetGraph()->GetSpaceDimension();
314
315 int i = -1;
316 while (GetPoint(gloCoord, ++i))
317 {
318 // Overwrite gloCoord[2] for 3DH1D using m_outputPlane
319 // if it is defined
321 {
322 int nplanes = pFields[0]->GetHomogeneousBasis()->GetZ().size();
323 NekDouble lhom = pFields[0]->GetHomoLen();
324 int plane;
325 if (m_outputPlane != -1)
326 {
327 plane = m_outputPlane;
328 }
329 else
330 {
331 // Pick plane near the point
332 plane = round((gloCoord[coordim] * nplanes) / lhom);
333 }
334 if (m_waveSpace)
335 {
336 m_planeIDs.push_back(plane);
337 }
338 NekDouble Z = (pFields[0]->GetHomogeneousBasis()->GetZ())[plane];
339 Z = (Z + 1) * lhom / 2;
340 if (fabs(gloCoord[coordim] - Z) >
342 vComm->GetRank() == 0)
343 {
344 cout << "Resetting History point from z = " << gloCoord[coordim]
345 << " to z = " << Z << endl;
346 }
347 gloCoord[coordim] = Z;
348 }
349
350 // Save a copy of the global coordinate
351 GloCoords.push_back(
352 Array<OneD, const NekDouble>(gloCoord.size(), gloCoord.data()));
353 }
354
355 // Check which history points this process is responsible for
356 const int vRank = vComm->GetRank();
357 const int vColumnRank = vComm->GetColumnComm()->GetRank();
358 const int vRowRank = vComm->GetRowComm()->GetRank();
359 m_historyPointsSize = GloCoords.size();
360 Array<OneD, int> procList(m_historyPointsSize, -1);
361 Array<OneD, int> idList(m_historyPointsSize, -1);
362
363 Array<OneD, NekDouble> locCoord(pFields[0]->GetShapeDimension(), 0.);
364 std::vector<Array<OneD, const NekDouble>> LocCoords;
365
366 // For each history point, find the element that contains it. If no element
367 // is available on this process, idList[i] = -1 and procList[i] = -1
368 for (i = 0; i < m_historyPointsSize; ++i)
369 {
370 // Create a copy of the global coordinate
371 gloCoord = GloCoords[i];
372
373 // Determine the expansion and local coordinates
375 {
376 idList[i] = pFields[0]->GetPlane(0)->GetExpIndex(
377 gloCoord, locCoord, NekConstants::kGeomFactorsTol);
378 }
379 else
380 {
381 idList[i] = pFields[0]->GetExpIndex(gloCoord, locCoord,
383 }
384
385 for (int j = 0; j < locCoord.size(); ++j)
386 {
387 locCoord[j] = std::max(locCoord[j], -1.0);
388 locCoord[j] = std::min(locCoord[j], 1.0);
389 }
390
391 // If this process owns an element that contains the history point, we
392 // flag this by updating the corresponding element in procList
393 if (idList[i] != -1 && vColumnRank == 0)
394 {
395 procList[i] = vRank;
396 }
397
398 // Save a copy of the local coordinates
399 LocCoords.push_back(
400 Array<OneD, const NekDouble>(locCoord.size(), locCoord.data()));
401 }
402
403 // Clear list before populating it
404 m_historyList.clear();
405
406 // Share process ranks among all processes. If multiple processes found
407 // that they own a history point, procList will only contain the rank of
408 // the process with the largest rank
409 vComm->AllReduce(procList, LibUtilities::ReduceMax);
410 for (i = 0; i < m_historyPointsSize; ++i)
411 {
412 // If process owns a history point, it saves a reference to its global
413 // and local coordinates. Note that, if multiple processes have found
414 // that they own a history point, which can happen if the point is
415 // located on a partition boundary, only the process with the highest
416 // rank will store the point.
417 if (procList[i] == vRowRank)
418 {
419 // Note that the reference counting in the Array class prevents the
420 // data stored in GloCoords and LocCoords from being deallocated
421 // when they go out of scope, since m_historyList will increase the
422 // reference counter by 1 when it gets a reference to the data.
423 m_historyList.push_back(
424 std::make_tuple(GloCoords[i], LocCoords[i], i, idList[i]));
425 }
426 }
427
428 // Root process must have access to all coordinates for writing the file
429 if (vComm->TreatAsRankZero())
430 {
432 Array<OneD, Array<OneD, const NekDouble>>(m_historyPointsSize);
433
434 for (i = 0; i < m_historyPointsSize; ++i)
435 {
436 // Create a reference to the underlying array. Note that the data
437 // stored in GloCoords won't be deallocated because the Array class
438 // has a reference counter.
439 m_historyPoints[i] = GloCoords[i];
440 }
441 }
442
443 // Check that each history point is allocated to a process
444 if (vComm->TreatAsRankZero())
445 {
446 for (i = 0; i < m_historyPointsSize; ++i)
447 {
448 gloCoord = m_historyPoints[i];
449
450 // Write an error if no process owns history point
451 ASSERTL0(procList[i] != -1,
452 "History point " +
453 boost::lexical_cast<std::string>(gloCoord[0]) + ", " +
454 boost::lexical_cast<std::string>(gloCoord[1]) + ", " +
455 boost::lexical_cast<std::string>(gloCoord[2]) +
456 " cannot be found in the mesh.");
457 }
458
459 m_session->MatchSolverInfo("Driver", "Adaptive", m_adaptive, false);
460 }
461 v_Update(pFields, time);
462}
Array< OneD, Array< OneD, const NekDouble > > m_historyPoints
bool GetPoint(Array< OneD, NekDouble > gloCoord, int I)
SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
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
static const NekDouble kVertexTheSameDouble
static const NekDouble kGeomFactorsTol

References ASSERTL0, GetPoint(), Nektar::NekConstants::kGeomFactorsTol, Nektar::NekConstants::kVertexTheSameDouble, m_adaptive, m_historyList, m_historyPoints, m_historyPointsSize, m_historyPointStream, m_isHomogeneous1D, m_outputPlane, m_planeIDs, Nektar::SolverUtils::Filter::m_session, m_waveSpace, Nektar::LibUtilities::ReduceMax, and v_Update().

◆ v_IsTimeDependent()

bool Nektar::SolverUtils::FilterHistoryPoints::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 698 of file FilterHistoryPoints.cpp.

699{
700 return true;
701}

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Reimplemented in Nektar::FilterCellHistoryPoints.

Definition at line 467 of file FilterHistoryPoints.cpp.

470{
471 // Only output every m_outputFrequency.
472 if ((m_index++) % m_outputFrequency)
473 {
474 return;
475 }
476
477 const int numFields = pFields.size();
478 const int coordim = pFields[0]->GetGraph()->GetSpaceDimension();
479 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
480 Array<OneD, NekDouble> data(m_historyPointsSize * numFields, 0.0);
481 Array<OneD, NekDouble> physvals;
482 Array<OneD, NekDouble> locCoord;
483
484 // Pull out data values field by field
485 Array<OneD, NekDouble> gloCoord(3, 0.0);
486 for (int j = 0; j < numFields; ++j)
487 {
489 {
490 Array<OneD, const unsigned int> planes = pFields[j]->GetZIDs();
491 int nPlanes = pFields[j]->GetHomogeneousBasis()->GetZ().size();
492 NekDouble lHom = pFields[j]->GetHomoLen();
493 for (auto &x : m_historyList)
494 {
495 ASSERTL0(pFields[j]->GetWaveSpace(),
496 "HistoryPoints in Homogeneous1D require that solution "
497 "is in wavespace");
498 gloCoord = std::get<0>(x);
499 locCoord = std::get<1>(x);
500 const int indx = std::get<2>(x);
501 const int expId = std::get<3>(x);
502
503 NekDouble value = 0.0;
504 NekDouble BetaT =
505 2. * M_PI * fmod(gloCoord[coordim], lHom) / lHom;
506 for (size_t n = 0; n < planes.size(); ++n)
507 {
508 if (m_waveSpace && planes[n] != m_planeIDs[indx])
509 {
510 continue;
511 }
512 physvals = pFields[j]->GetPlane(n)->UpdatePhys() +
513 pFields[j]->GetPhys_Offset(expId);
514
515 // transform elemental data if required.
516 if (pFields[j]->GetPhysState() == false)
517 {
518 pFields[j]->GetPlane(n)->GetExp(expId)->BwdTrans(
519 pFields[j]->GetPlane(n)->GetCoeffs() +
520 pFields[j]->GetCoeff_Offset(expId),
521 physvals);
522 }
523 // Interpolate data
524 NekDouble coeff =
525 pFields[j]->GetPlane(n)->GetExp(expId)->StdPhysEvaluate(
526 locCoord, physvals);
527
528 if (m_waveSpace)
529 {
530 value = coeff;
531 }
532 else
533 {
534 if (planes[n] == 0)
535 {
536 value += coeff;
537 }
538 else if (planes[n] == 1)
539 {
540 value += cos(0.5 * nPlanes * BetaT) * coeff;
541 }
542 else if (planes[n] % 2 == 0)
543 {
544 NekDouble phase = (planes[n] >> 1) * BetaT;
545 value += cos(phase) * coeff;
546 }
547 else
548 {
549 NekDouble phase = (planes[n] >> 1) * BetaT;
550 value += -sin(phase) * coeff;
551 }
552 }
553 }
554 // store data
555 data[indx * numFields + j] = value;
556 }
557 }
558 else
559 {
560 for (auto &x : m_historyList)
561 {
562 locCoord = std::get<1>(x);
563 const int indx = std::get<2>(x);
564 const int expId = std::get<3>(x);
565
566 physvals = pFields[j]->UpdatePhys() +
567 pFields[j]->GetPhys_Offset(expId);
568
569 // transform elemental data if required.
570 if (pFields[j]->GetPhysState() == false)
571 {
572 pFields[j]->GetExp(expId)->BwdTrans(
573 pFields[j]->GetCoeffs() +
574 pFields[j]->GetCoeff_Offset(expId),
575 physvals);
576 }
577
578 // interpolate point
579 data[indx * numFields + j] =
580 pFields[j]->GetExp(expId)->StdPhysEvaluate(locCoord,
581 physvals);
582 }
583 }
584 }
585
586 // Exchange history data
587 // This could be improved to reduce communication but works for now
588 vComm->AllReduce(data, LibUtilities::ReduceSum);
589
590 // TODO: Why not only call this routine if we are rank 0?
591 v_WriteData(vComm->GetRank(), data, numFields, time);
592}
SOLVER_UTILS_EXPORT void v_WriteData(const int &rank, const Array< OneD, NekDouble > &data, const int &numFields, const NekDouble &time)

References ASSERTL0, m_historyList, m_historyPointsSize, m_index, m_isHomogeneous1D, m_outputFrequency, m_planeIDs, m_waveSpace, Nektar::LibUtilities::ReduceSum, and v_WriteData().

Referenced by v_Initialise().

◆ v_WriteData()

void Nektar::SolverUtils::FilterHistoryPoints::v_WriteData ( const int &  rank,
const Array< OneD, NekDouble > &  data,
const int &  numFields,
const NekDouble time 
)
protected

Definition at line 594 of file FilterHistoryPoints.cpp.

598{
599 // Only the root process writes out history data
600 if (rank == 0)
601 {
602 Array<OneD, NekDouble> gloCoord(3, 0.0);
603 if (!m_outputOneFile || m_index == 1)
604 {
605 std::stringstream vOutputFilename;
606 if (m_outputOneFile)
607 {
608 vOutputFilename << m_outputFile << ".his";
609 }
610 else
611 {
612 vOutputFilename << m_outputFile << "_" << m_outputIndex
613 << ".his";
614 }
616 if (m_adaptive)
617 {
618 m_outputStream.open(vOutputFilename.str().c_str(),
619 ofstream::app);
620 }
621 else
622 {
623 m_outputStream.open(vOutputFilename.str().c_str());
624 }
625 m_outputStream << "# History data for variables (:";
626
627 for (int i = 0; i < numFields; ++i)
628 {
629 m_outputStream << m_session->GetVariable(i) << ",";
630 }
631
633 {
634 m_outputStream << ") at points:" << endl;
635 }
636 else
637 {
638 m_outputStream << ") at points:" << endl;
639 }
640
641 for (int i = 0; i < m_historyPoints.size(); ++i)
642 {
643 gloCoord = m_historyPoints[i];
644
645 m_outputStream << "# " << boost::format("%6.0f") % i;
646 m_outputStream << " " << boost::format("%25.19e") % gloCoord[0];
647 m_outputStream << " " << boost::format("%25.19e") % gloCoord[1];
648 m_outputStream << " " << boost::format("%25.19e") % gloCoord[2];
649 m_outputStream << endl;
650 }
651
653 {
654 if (m_waveSpace)
655 {
656 m_outputStream << "# (in Wavespace)" << endl;
657 }
658 }
659 }
660
661 // Write data values point by point
662 for (int k = 0; k < m_historyPoints.size(); ++k)
663 {
664 m_outputStream << boost::format("%25.19e") % time;
665 for (int j = 0; j < numFields; ++j)
666 {
668 << " "
669 << boost::format("%25.19e") % data[k * numFields + j];
670 }
671 m_outputStream << endl;
672 }
673
674 if (!m_outputOneFile)
675 {
676 m_outputStream.close();
677 ;
678 }
679 }
680}

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

Referenced by v_Update().

Friends And Related Function Documentation

◆ MemoryManager< FilterHistoryPoints >

friend class MemoryManager< FilterHistoryPoints >
friend

Definition at line 1 of file FilterHistoryPoints.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterHistoryPoints::className
static
Initial value:
=
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
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.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 63 of file FilterHistoryPoints.h.

◆ m_adaptive

bool Nektar::SolverUtils::FilterHistoryPoints::m_adaptive
protected

◆ m_historyList

std::list<std::tuple<Array<OneD, const NekDouble>, Array<OneD, const NekDouble>, int, int> > Nektar::SolverUtils::FilterHistoryPoints::m_historyList
protected

◆ m_historyPoints

Array<OneD, Array<OneD, const NekDouble> > Nektar::SolverUtils::FilterHistoryPoints::m_historyPoints
protected
Initial value:
=
Array<OneD, Array<OneD, const NekDouble>>(0)

Definition at line 88 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), Nektar::FilterCellHistoryPoints::v_Update(), v_WriteData(), and Nektar::FilterCellHistoryPoints::v_WriteData().

◆ m_historyPointsSize

size_t Nektar::SolverUtils::FilterHistoryPoints::m_historyPointsSize = 0
protected

Definition at line 90 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), and v_Update().

◆ m_historyPointStream

std::stringstream Nektar::SolverUtils::FilterHistoryPoints::m_historyPointStream
protected

Definition at line 100 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), GetPoint(), and v_Initialise().

◆ m_index

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_index = 0
protected

◆ m_isHomogeneous1D

bool Nektar::SolverUtils::FilterHistoryPoints::m_isHomogeneous1D
protected

◆ m_outputFile

std::string Nektar::SolverUtils::FilterHistoryPoints::m_outputFile
protected

◆ m_outputFrequency

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_outputFrequency
protected

◆ m_outputIndex

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_outputIndex = 0
protected

◆ m_outputOneFile

bool Nektar::SolverUtils::FilterHistoryPoints::m_outputOneFile
protected

◆ m_outputPlane

int Nektar::SolverUtils::FilterHistoryPoints::m_outputPlane
protected

plane to take history point from if using a homogeneous1D expansion

Definition at line 94 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), v_Initialise(), and Nektar::FilterCellHistoryPoints::v_Update().

◆ m_outputStream

std::ofstream Nektar::SolverUtils::FilterHistoryPoints::m_outputStream
protected

◆ m_planeIDs

std::vector<int> Nektar::SolverUtils::FilterHistoryPoints::m_planeIDs
protected

Definition at line 95 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), and v_Update().

◆ m_pointDatMap

std::map<LibUtilities::PtsType, Array<OneD, NekDouble> > Nektar::SolverUtils::FilterHistoryPoints::m_pointDatMap
protected

Definition at line 110 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), and GetPoint().

◆ m_pointNumMap

std::map<LibUtilities::PtsType, Array<OneD, int> > Nektar::SolverUtils::FilterHistoryPoints::m_pointNumMap
protected

Definition at line 111 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), and GetPoint().

◆ m_waveSpace

bool Nektar::SolverUtils::FilterHistoryPoints::m_waveSpace
protected