Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::SolverUtils::FilterHistoryPoints Class Reference

#include <FilterHistoryPoints.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterHistoryPoints (const LibUtilities::SessionReaderSharedPtr &pSession, const std::map< std::string, std::string > &pParams)
SOLVER_UTILS_EXPORT ~FilterHistoryPoints ()
- Public Member Functions inherited from Nektar::SolverUtils::Filter
SOLVER_UTILS_EXPORT Filter (const LibUtilities::SessionReaderSharedPtr &pSession)
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::map< std::string, std::string > &pParams)
 Creates an instance of this class.

Static Public Attributes

static std::string className = GetFilterFactory().RegisterCreatorFunction("HistoryPoints", FilterHistoryPoints::create)
 Name of the class.

Protected Member Functions

virtual void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual bool v_IsTimeDependent ()

Private Attributes

SpatialDomains::PointGeomVector m_historyPoints
unsigned int m_index
unsigned int m_outputFrequency
unsigned int m_outputPlane
bool m_isHomogeneous1D
std::string m_outputFile
std::ofstream m_outputStream
std::stringstream m_historyPointStream
std::list< std::pair
< SpatialDomains::PointGeomSharedPtr,
Array< OneD, NekDouble > > > 
m_historyList
std::map< int, int > m_historyLocalPointMap

Friends

class MemoryManager< FilterHistoryPoints >

Additional Inherited Members

- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session

Detailed Description

Definition at line 45 of file FilterHistoryPoints.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::FilterHistoryPoints::FilterHistoryPoints ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::map< std::string, std::string > &  pParams 
)

Definition at line 49 of file FilterHistoryPoints.cpp.

References ASSERTL0, m_historyPointStream, m_index, m_isHomogeneous1D, m_outputFile, m_outputFrequency, m_outputPlane, and Nektar::SolverUtils::Filter::m_session.

:
Filter(pSession)
{
if (pParams.find("OutputFile") == pParams.end())
{
m_outputFile = m_session->GetSessionName();
}
else
{
ASSERTL0(!(pParams.find("OutputFile")->second.empty()),
"Missing parameter 'OutputFile'.");
m_outputFile = pParams.find("OutputFile")->second;
}
if (!(m_outputFile.length() >= 4
&& m_outputFile.substr(m_outputFile.length() - 4) == ".his"))
{
m_outputFile += ".his";
}
if (pParams.find("OutputFrequency") == pParams.end())
{
}
else
{
m_outputFrequency = atoi(pParams.find("OutputFrequency")->second.c_str());
}
m_session->MatchSolverInfo("Homogeneous","1D",m_isHomogeneous1D,false);
{
if (pParams.find("OutputPlane") == pParams.end())
{
}
else
{
m_outputPlane = atoi(pParams.find("OutputPlane")->second.c_str());
}
}
ASSERTL0(pParams.find("Points") != pParams.end(),
"Missing parameter 'Points'.");
m_historyPointStream.str(pParams.find("Points")->second);
m_index = 0;
}
Nektar::SolverUtils::FilterHistoryPoints::~FilterHistoryPoints ( )

Definition at line 104 of file FilterHistoryPoints.cpp.

{
}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 51 of file FilterHistoryPoints.h.

{
//p->InitObject();
return p;
}
void Nektar::SolverUtils::FilterHistoryPoints::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
protectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 413 of file FilterHistoryPoints.cpp.

References m_outputStream.

{
if (pFields[0]->GetComm()->GetRank() == 0)
{
m_outputStream.close();
}
}
void Nektar::SolverUtils::FilterHistoryPoints::v_Initialise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
protectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 113 of file FilterHistoryPoints.cpp.

References ASSERTL0, Nektar::NekConstants::kGeomFactorsTol, Nektar::NekConstants::kVertexTheSameDouble, m_historyList, m_historyLocalPointMap, m_historyPoints, m_historyPointStream, m_index, m_isHomogeneous1D, m_outputFile, m_outputPlane, m_outputStream, Nektar::SolverUtils::Filter::m_session, Nektar::LibUtilities::ReduceMax, and v_Update().

{
"No history points in stream.");
m_index = 0;
m_historyList.clear();
// Read history points
Array<OneD, NekDouble> gloCoord(3,0.0);
int dim = pFields[0]->GetGraph()->GetSpaceDimension();
int i = 0;
while (!m_historyPointStream.fail())
{
m_historyPointStream >> gloCoord[0] >> gloCoord[1] >> gloCoord[2];
if(m_isHomogeneous1D) // overwrite with plane z
{
NekDouble Z = (pFields[0]->GetHomogeneousBasis()->GetZ())[m_outputPlane];
if(fabs(gloCoord[2] - Z) > NekConstants::kVertexTheSameDouble)
{
cout <<"Reseting History point from " << gloCoord[2] <<
" to " << Z << endl;
}
gloCoord[2] = Z;
}
if (!m_historyPointStream.fail())
{
::AllocateSharedPtr(dim, i, gloCoord[0],
gloCoord[1], gloCoord[2]);
m_historyPoints.push_back(vert);
++i;
}
}
// Determine the unique process responsible for each history point
// For points on a partition boundary, must select a single process
LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
int vRank = vComm->GetRank();
Array<OneD, int> procList(m_historyPoints.size(), -1);
Array<OneD, int> idList(m_historyPoints.size());
std::vector<Array<OneD, NekDouble> > LocCoords;
for (i = 0; i < m_historyPoints.size(); ++i)
{
Array<OneD, NekDouble> locCoords(3);
// Determine the expansion and local coordinates
m_historyPoints[i]->GetCoords( gloCoord[0],
gloCoord[1],
gloCoord[2]);
idList[i] = pFields[0]->GetExpIndex(gloCoord,locCoords,
// Check if the reverse mapping of the local coordinates gives
// the correct coordinates of the history point. This ensures
// that the correct element is chosen in the manifold case.
if (idList[i] != -1)
{
pFields[0]->GetExp(idList[i])->GetGeom();
Array<OneD, NekDouble> coordVals(e->GetTotPoints());
NekDouble distance = 0.0;
for (int j = 0; j < g->GetCoordim(); ++j)
{
e->BwdTrans(g->GetCoeffs(j), coordVals);
NekDouble x = e->PhysEvaluate(locCoords, coordVals)
- gloCoord[j];
distance += x*x;
}
{
idList[i] = -1;
}
}
// Save Local coordinates for later
LocCoords.push_back(locCoords);
// Set element id to Vid of m_history point for later use
m_historyPoints[i]->SetVid(idList[i]);
// If a matching element is found on this process, note the
// process ID
if (idList[i] != -1)
{
{
int j;
Array<OneD, const unsigned int> IDs = pFields[0]->GetZIDs();
for(j = 0; j < IDs.num_elements(); ++j)
{
if(IDs[j] == m_outputPlane)
{
break;
}
}
if(j != IDs.num_elements())
{
procList[i] = vRank;
}
}
else
{
procList[i] = vRank;
}
}
}
// Reduce process IDs for all history points. The process with
// largest rank will handle the history point
vComm->AllReduce(procList, LibUtilities::ReduceMax);
// Determine the element in which each history point resides.
// If point is not in mesh (on this process), id is -1.
for (i = 0; i < m_historyPoints.size(); ++i)
{
// If point lies on partition boundary, only the proc with max
// rank retains possession.
if (procList[i] != vRank)
{
idList[i] = -1;
}
// If the current process owns this history point, add it to its
// local list of history points.
if (idList[i] != -1)
{
m_historyList.push_back(
Array<OneD, NekDouble> >
(m_historyPoints[i], LocCoords[i]));
}
}
// Collate the element ID list across processes and check each
// history point is allocated to a process
vComm->AllReduce(idList, LibUtilities::ReduceMax);
if (vComm->GetRank() == 0)
{
for (i = 0; i < m_historyPoints.size(); ++i)
{
m_historyPoints[i]->GetCoords( gloCoord[0],
gloCoord[1],
gloCoord[2]);
ASSERTL0(idList[i] != -1, "History point " +
boost::lexical_cast<std::string>(gloCoord[0]) + ", " +
boost::lexical_cast<std::string>(gloCoord[1]) + ", " +
boost::lexical_cast<std::string>(gloCoord[2]) +
" cannot be found in the mesh.");
}
// Open output stream
m_outputStream << "# History data for variables (:";
for (i = 0; i < pFields.num_elements(); ++i)
{
m_outputStream << m_session->GetVariable(i) <<",";
}
{
m_outputStream << ") at points:";
}
else
{
m_outputStream << ") at points:" << endl;
}
for (i = 0; i < m_historyPoints.size(); ++i)
{
m_historyPoints[i]->GetCoords( gloCoord[0],
gloCoord[1],
gloCoord[2]);
m_outputStream << "# \t" << i;
m_outputStream.width(8);
m_outputStream << gloCoord[0];
m_outputStream.width(8);
m_outputStream << gloCoord[1];
m_outputStream.width(8);
m_outputStream << gloCoord[2];
m_outputStream << endl;
}
{
m_outputStream << "(in Wavespace)" << endl;
}
}
v_Update(pFields, time);
}
bool Nektar::SolverUtils::FilterHistoryPoints::v_IsTimeDependent ( )
protectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 425 of file FilterHistoryPoints.cpp.

{
return true;
}
void Nektar::SolverUtils::FilterHistoryPoints::v_Update ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
protectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 320 of file FilterHistoryPoints.cpp.

References Nektar::iterator, m_historyList, m_historyLocalPointMap, m_historyPoints, m_index, m_isHomogeneous1D, m_outputFrequency, m_outputPlane, m_outputStream, and Nektar::LibUtilities::ReduceSum.

Referenced by v_Initialise().

{
// Only output every m_outputFrequency.
{
return;
}
int j = 0;
int k = 0;
int numPoints = m_historyPoints.size();
int numFields = pFields.num_elements();
LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
Array<OneD, NekDouble> data(numPoints*numFields, 0.0);
Array<OneD, NekDouble> gloCoord(3, 0.0);
std::list<std::pair<SpatialDomains::PointGeomSharedPtr, Array<OneD, NekDouble> > >::iterator x;
Array<OneD, NekDouble> physvals;
Array<OneD, NekDouble> locCoord;
int expId;
// Pull out data values field by field
for (j = 0; j < numFields; ++j)
{
{
for (k = 0, x = m_historyList.begin(); x != m_historyList.end();
++x, ++k)
{
locCoord = (*x).second;
expId = (*x).first->GetVid();
physvals = pFields[j]->GetPlane(m_outputPlane)->UpdatePhys() + pFields[j]->GetPhys_Offset(expId);
// transform elemental data if required.
if(pFields[j]->GetPhysState() == false)
{
pFields[j]->GetPlane(m_outputPlane)->GetExp(expId)->BwdTrans(pFields[j]->GetPlane(m_outputPlane)->GetCoeffs() + pFields[j]->GetCoeff_Offset(expId),physvals);
}
// interpolate point can do with zero plane methods
data[m_historyLocalPointMap[k]*numFields+j] = pFields[j]->GetExp(expId)->StdPhysEvaluate(locCoord,physvals);
}
}
else
{
for (k = 0, x = m_historyList.begin(); x != m_historyList.end(); ++x, ++k)
{
locCoord = (*x).second;
expId = (*x).first->GetVid();
physvals = pFields[j]->UpdatePhys() + pFields[j]->GetPhys_Offset(expId);
// transform elemental data if required.
if(pFields[j]->GetPhysState() == false)
{
pFields[j]->GetExp(expId)->BwdTrans(pFields[j]->GetCoeffs() + pFields[j]->GetCoeff_Offset(expId),physvals);
}
// interpolate point
data[m_historyLocalPointMap[k]*numFields+j] = pFields[j]->GetExp(expId)->StdPhysEvaluate(locCoord,physvals);
}
}
}
// Exchange history data
// This could be improved to reduce communication but works for now
vComm->AllReduce(data, LibUtilities::ReduceSum);
// Only the root process writes out history data
if (vComm->GetRank() == 0)
{
// Write data values point by point
for (k = 0; k < m_historyPoints.size(); ++k)
{
m_outputStream.width(8);
m_outputStream << setprecision(6) << time;
for (int j = 0; j < numFields; ++j)
{
m_outputStream.width(25);
m_outputStream << setprecision(16) << data[k*numFields+j];
}
m_outputStream << endl;
}
}
}

Friends And Related Function Documentation

friend class MemoryManager< FilterHistoryPoints >
friend

Definition at line 48 of file FilterHistoryPoints.h.

Member Data Documentation

std::string Nektar::SolverUtils::FilterHistoryPoints::className = GetFilterFactory().RegisterCreatorFunction("HistoryPoints", FilterHistoryPoints::create)
static

Name of the class.

Definition at line 60 of file FilterHistoryPoints.h.

std::list<std::pair<SpatialDomains::PointGeomSharedPtr, Array<OneD, NekDouble> > > Nektar::SolverUtils::FilterHistoryPoints::m_historyList
private

Definition at line 83 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), and v_Update().

std::map<int, int > Nektar::SolverUtils::FilterHistoryPoints::m_historyLocalPointMap
private

Definition at line 84 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), and v_Update().

SpatialDomains::PointGeomVector Nektar::SolverUtils::FilterHistoryPoints::m_historyPoints
private

Definition at line 74 of file FilterHistoryPoints.h.

Referenced by v_Initialise(), and v_Update().

std::stringstream Nektar::SolverUtils::FilterHistoryPoints::m_historyPointStream
private

Definition at line 81 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), and v_Initialise().

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_index
private

Definition at line 75 of file FilterHistoryPoints.h.

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

bool Nektar::SolverUtils::FilterHistoryPoints::m_isHomogeneous1D
private

Definition at line 78 of file FilterHistoryPoints.h.

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

std::string Nektar::SolverUtils::FilterHistoryPoints::m_outputFile
private

Definition at line 79 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), and v_Initialise().

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_outputFrequency
private

Definition at line 76 of file FilterHistoryPoints.h.

Referenced by FilterHistoryPoints(), and v_Update().

unsigned int Nektar::SolverUtils::FilterHistoryPoints::m_outputPlane
private

Definition at line 77 of file FilterHistoryPoints.h.

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

std::ofstream Nektar::SolverUtils::FilterHistoryPoints::m_outputStream
private

Definition at line 80 of file FilterHistoryPoints.h.

Referenced by v_Finalise(), v_Initialise(), and v_Update().