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::FilterBenchmark Class Reference

Records activation and repolarisation times. More...

#include <FilterBenchmark.h>

Inheritance diagram for Nektar::FilterBenchmark:
Inheritance graph
[legend]
Collaboration diagram for Nektar::FilterBenchmark:
Collaboration graph
[legend]

Public Member Functions

 FilterBenchmark (const LibUtilities::SessionReaderSharedPtr &pSession, const std::map< std::string, std::string > &pParams)
 Construct the benchmark filter.
virtual ~FilterBenchmark ()
 Destructor for the benchmark filter.
- 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 SolverUtils::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
 Name of the class.

Protected Member Functions

virtual void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 Initialises the benchmark filter and allocates storage.
virtual void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 Update recorded times.
virtual void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 Finalises the benchmark filter and write out recorded data.
virtual bool v_IsTimeDependent ()
 Identifies that the benchmark filter is time dependent.

Private Attributes

std::vector< Array< OneD,
NekDouble > > 
m_threshold
 Storage for activation and repolarisation times.
Array< OneD, int > m_idx
 Number of activations and repolarisations detected for each point.
Array< OneD, int > m_polarity
 Indicates if the previous event was an activation or repolarisation.
NekDouble m_startTime
 Time at which to start detecting activations and repolarisations.
NekDouble m_thresholdValue
 Value at which tissue is considered active.
NekDouble m_initialValue
 Initial time to use in storage array.
std::string m_outputFile
 Filename of output files.
LibUtilities::FieldIOSharedPtr m_fld
 FieldIO object used for writing output files.

Friends

class MemoryManager< FilterBenchmark >

Additional Inherited Members

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

Detailed Description

Records activation and repolarisation times.

This class records the sequence of activation and repolarisation times across the entire domain into a two-dimensional storage structure. At each timestep, the voltage at each point in the domain is examined to identify if it has crossed the threshold value. If so, the time of crossing is recorded. Auxiliary arrays hold the current index of each point (i.e. the number of crossings of the threshold) and the type of the last crossing (activation or repolarisation).

Definition at line 45 of file FilterBenchmark.h.

Constructor & Destructor Documentation

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

Construct the benchmark filter.

Parameters
pSessionSession reader for IO
pParamsParameters of filter

Definition at line 61 of file FilterBenchmark.cpp.

References ASSERTL0, m_fld, m_initialValue, m_outputFile, m_startTime, and m_thresholdValue.

: Filter(pSession)
{
ASSERTL0(pParams.find("ThresholdValue") != pParams.end(),
"Missing parameter 'ThresholdValue'.");
m_thresholdValue = atof(pParams.find("ThresholdValue")->second.c_str());
ASSERTL0(pParams.find("InitialValue") != pParams.end(),
"Missing parameter 'InitialValue'.");
m_initialValue = atof(pParams.find("InitialValue")->second.c_str());
ASSERTL0(!(pParams.find("OutputFile")->second.empty()),
"Missing parameter 'OutputFile'.");
m_outputFile = pParams.find("OutputFile")->second;
m_startTime = 0.0;
if (pParams.find("StartTime") != pParams.end())
{
m_startTime = atof(pParams.find("StartTime")->second.c_str());
}
pSession->GetComm());
}
Nektar::FilterBenchmark::~FilterBenchmark ( )
virtual

Destructor for the benchmark filter.

Definition at line 91 of file FilterBenchmark.cpp.

{
}

Member Function Documentation

static SolverUtils::FilterSharedPtr Nektar::FilterBenchmark::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 FilterBenchmark.h.

void Nektar::FilterBenchmark::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
protectedvirtual

Finalises the benchmark filter and write out recorded data.

Writes out the crossings to file.

Parameters
pFieldsField storage expansion list.
timeCurrent time.

Implements Nektar::SolverUtils::Filter.

Definition at line 175 of file FilterBenchmark.cpp.

References m_fld, m_outputFile, and m_threshold.

{
for (int i = 0; i < m_threshold.size() - 1; ++i)
{
std::stringstream vOutputFilename;
vOutputFilename << m_outputFile << "_" << i << ".fld";
std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
= pFields[0]->GetFieldDefinitions();
std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
Array<OneD, NekDouble> vCoeffs(pFields[0]->GetNcoeffs());
pFields[0]->FwdTrans_IterPerExp(m_threshold[i], vCoeffs);
// copy Data into FieldData and set variable
for(int i = 0; i < FieldDef.size(); ++i)
{
// Could do a search here to find correct variable
FieldDef[i]->m_fields.push_back("m");
pFields[0]->AppendFieldData(FieldDef[i], FieldData[i], vCoeffs);
}
m_fld->Write(vOutputFilename.str(),FieldDef,FieldData);
}
}
void Nektar::FilterBenchmark::v_Initialise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
protectedvirtual

Initialises the benchmark filter and allocates storage.

Implements Nektar::SolverUtils::Filter.

Definition at line 102 of file FilterBenchmark.cpp.

References m_idx, m_initialValue, m_polarity, and m_threshold.

{
m_threshold.push_back(Array<OneD, NekDouble>(
pFields[0]->GetNpoints(), m_initialValue));
m_idx = Array<OneD, int> (pFields[0]->GetNpoints(), 0);
m_polarity = Array<OneD, int> (pFields[0]->GetNpoints(), -1);
}
bool Nektar::FilterBenchmark::v_IsTimeDependent ( )
protectedvirtual

Identifies that the benchmark filter is time dependent.

Returns
This filter is time dependent.

Implements Nektar::SolverUtils::Filter.

Definition at line 207 of file FilterBenchmark.cpp.

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

Update recorded times.

Checks each point in the domain to determine if it has crossed the threshold. The direction of crossing is determined. Additional storage is allocated if needed.

Parameters
pFieldsField storage expansion lists
timeCurrent time

Implements Nektar::SolverUtils::Filter.

Definition at line 121 of file FilterBenchmark.cpp.

References m_idx, m_initialValue, m_polarity, m_startTime, m_threshold, m_thresholdValue, Nektar::LibUtilities::ReduceMax, and Vmath::Vmax().

{
// Only proceed if the start time has passed
if (time < m_startTime)
{
return;
}
// Examine each point in turn
for (int i = 0; i < pFields[0]->GetNpoints(); ++i)
{
if ((m_polarity[i] == -1 &&
pFields[0]->GetPhys()[i] > m_thresholdValue) ||
(m_polarity[i] == 1 &&
pFields[0]->GetPhys()[i] < m_thresholdValue))
{
// If APD less than 50ms, remove last activation
if (m_polarity[i] == 1 &&
time - m_threshold[m_idx[i]][i] < 50)
{
m_idx[i]--;
}
else
{
m_threshold[m_idx[i]][i] = time;
m_idx[i]++;
}
// Update polarity of last crossing
m_polarity[i] *= -1;
}
}
// Allocate additional storage if any point has as many crossings as
// current storage permits.
int max_idx = Vmath::Vmax(pFields[0]->GetNpoints(), m_idx, 1);
pFields[0]->GetSession()->GetComm()->AllReduce(max_idx,
if (m_threshold.size() == max_idx)
{
m_threshold.push_back(Array<OneD, NekDouble>(
pFields[0]->GetNpoints(), m_initialValue));
}
}

Friends And Related Function Documentation

friend class MemoryManager< FilterBenchmark >
friend

Definition at line 48 of file FilterBenchmark.h.

Member Data Documentation

std::string Nektar::FilterBenchmark::className
static
Initial value:

Name of the class.

Definition at line 61 of file FilterBenchmark.h.

LibUtilities::FieldIOSharedPtr Nektar::FilterBenchmark::m_fld
private

FieldIO object used for writing output files.

Definition at line 103 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Finalise().

Array<OneD, int> Nektar::FilterBenchmark::m_idx
private

Number of activations and repolarisations detected for each point.

Definition at line 91 of file FilterBenchmark.h.

Referenced by v_Initialise(), and v_Update().

NekDouble Nektar::FilterBenchmark::m_initialValue
private

Initial time to use in storage array.

Definition at line 99 of file FilterBenchmark.h.

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

std::string Nektar::FilterBenchmark::m_outputFile
private

Filename of output files.

Definition at line 101 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Finalise().

Array<OneD, int> Nektar::FilterBenchmark::m_polarity
private

Indicates if the previous event was an activation or repolarisation.

Definition at line 93 of file FilterBenchmark.h.

Referenced by v_Initialise(), and v_Update().

NekDouble Nektar::FilterBenchmark::m_startTime
private

Time at which to start detecting activations and repolarisations.

Definition at line 95 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Update().

std::vector<Array<OneD, NekDouble> > Nektar::FilterBenchmark::m_threshold
private

Storage for activation and repolarisation times.

Definition at line 89 of file FilterBenchmark.h.

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

NekDouble Nektar::FilterBenchmark::m_thresholdValue
private

Value at which tissue is considered active.

Definition at line 97 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Update().