Nektar++
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:
[legend]

Public Member Functions

 FilterBenchmark (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 Construct the benchmark filter. More...
 
 ~FilterBenchmark () override
 Destructor for the benchmark filter. More...
 
- 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 Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

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

Private Attributes

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

Friends

class MemoryManager< FilterBenchmark >
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Filter
typedef std::map< std::string, std::string > ParamMap
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

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 44 of file FilterBenchmark.h.

Constructor & Destructor Documentation

◆ FilterBenchmark()

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

Construct the benchmark filter.

Parameters
pSessionSession reader for IO
pParamsParameters of filter

Definition at line 59 of file FilterBenchmark.cpp.

63 : Filter(pSession, pEquation)
64{
65 // ThresholdValue
66 auto it = pParams.find("ThresholdValue");
67 ASSERTL0(it != pParams.end(), "Missing parameter 'ThresholdValue'.");
68 LibUtilities::Equation equ1(m_session->GetInterpreter(), it->second);
69 m_thresholdValue = floor(equ1.Evaluate());
70
71 // InitialValue
72 it = pParams.find("InitialValue");
73 ASSERTL0(it != pParams.end(), "Missing parameter 'InitialValue'.");
74 LibUtilities::Equation equ2(m_session->GetInterpreter(), it->second);
75 m_initialValue = floor(equ2.Evaluate());
76
77 // OutputFile
78 it = pParams.find("OutputFile");
79 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
80 m_outputFile = it->second;
81
82 // StartTime
83 m_startTime = 0.0;
84 it = pParams.find("StartTime");
85 if (it != pParams.end())
86 {
87 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
88 m_startTime = floor(equ.Evaluate());
89 }
90
92}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::FieldIOSharedPtr m_fld
FieldIO object used for writing output files.
NekDouble m_thresholdValue
Value at which tissue is considered active.
std::string m_outputFile
Filename of output files.
NekDouble m_initialValue
Initial time to use in storage array.
NekDouble m_startTime
Time at which to start detecting activations and repolarisations.
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:195
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::FieldIO::CreateDefault(), Nektar::LibUtilities::Equation::Evaluate(), m_fld, m_initialValue, m_outputFile, Nektar::SolverUtils::Filter::m_session, m_startTime, and m_thresholdValue.

◆ ~FilterBenchmark()

Nektar::FilterBenchmark::~FilterBenchmark ( )
override

Destructor for the benchmark filter.

Definition at line 97 of file FilterBenchmark.cpp.

98{
99}

Member Function Documentation

◆ create()

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

◆ v_Finalise()

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

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.

178{
179 for (size_t j = 0; j < m_threshold.size() - 1; ++j)
180 {
181 std::stringstream vOutputFilename;
182 vOutputFilename << m_outputFile << "_" << j << ".fld";
183
184 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
185 pFields[0]->GetFieldDefinitions();
186 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
187
188 Array<OneD, NekDouble> vCoeffs(pFields[0]->GetNcoeffs());
189 pFields[0]->FwdTransLocalElmt(m_threshold[j], vCoeffs);
190
191 // copy Data into FieldData and set variable
192 for (size_t i = 0; i < FieldDef.size(); ++i)
193 {
194 // Could do a search here to find correct variable
195 FieldDef[i]->m_fields.push_back("m");
196 pFields[0]->AppendFieldData(FieldDef[i], FieldData[i], vCoeffs);
197 }
198
199 m_fld->Write(vOutputFilename.str(), FieldDef, FieldData);
200 }
201}
std::vector< Array< OneD, NekDouble > > m_threshold
Storage for activation and repolarisation times.

References m_fld, m_outputFile, and m_threshold.

◆ v_Initialise()

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

Initialises the benchmark filter and allocates storage.

Implements Nektar::SolverUtils::Filter.

Definition at line 106 of file FilterBenchmark.cpp.

109{
110 m_threshold.push_back(
111 Array<OneD, NekDouble>(pFields[0]->GetNpoints(), m_initialValue));
112
113 m_idx = Array<OneD, int>(pFields[0]->GetNpoints(), 0);
114 m_polarity = Array<OneD, int>(pFields[0]->GetNpoints(), -1);
115}
Array< OneD, int > m_polarity
Indicates if the previous event was an activation or repolarisation.
Array< OneD, int > m_idx
Number of activations and repolarisations detected for each point.

References m_idx, m_initialValue, m_polarity, and m_threshold.

◆ v_IsTimeDependent()

bool Nektar::FilterBenchmark::v_IsTimeDependent ( )
overrideprotectedvirtual

Identifies that the benchmark filter is time dependent.

Returns
This filter is time dependent.

Implements Nektar::SolverUtils::Filter.

Definition at line 206 of file FilterBenchmark.cpp.

207{
208 return true;
209}

◆ v_Update()

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

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 124 of file FilterBenchmark.cpp.

127{
128 // Only proceed if the start time has passed
129 if (time < m_startTime)
130 {
131 return;
132 }
133
134 // Examine each point in turn
135 size_t nPts = pFields[0]->GetNpoints();
136 for (size_t i = 0; i < nPts; ++i)
137 {
138 if ((m_polarity[i] == -1 &&
139 pFields[0]->GetPhys()[i] > m_thresholdValue) ||
140 (m_polarity[i] == 1 && pFields[0]->GetPhys()[i] < m_thresholdValue))
141 {
142 // If APD less than 50ms, remove last activation
143 if (m_polarity[i] == 1 && time - m_threshold[m_idx[i]][i] < 50)
144 {
145 m_idx[i]--;
147 }
148 else
149 {
150 m_threshold[m_idx[i]][i] = time;
151 m_idx[i]++;
152 }
153 // Update polarity of last crossing
154 m_polarity[i] *= -1;
155 }
156 }
157
158 // Allocate additional storage if any point has as many crossings as
159 // current storage permits.
160 size_t max_idx = Vmath::Vmax(pFields[0]->GetNpoints(), m_idx, 1);
161 pFields[0]->GetSession()->GetComm()->AllReduce(max_idx,
163 if (m_threshold.size() == max_idx)
164 {
165 m_threshold.push_back(
166 Array<OneD, NekDouble>(pFields[0]->GetNpoints(), m_initialValue));
167 }
168}
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644

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

Friends And Related Function Documentation

◆ MemoryManager< FilterBenchmark >

friend class MemoryManager< FilterBenchmark >
friend

Definition at line 1 of file FilterBenchmark.h.

Member Data Documentation

◆ className

std::string Nektar::FilterBenchmark::className
static
Initial value:
=
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 FilterBenchmark.h.

◆ m_fld

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

FieldIO object used for writing output files.

Definition at line 104 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Finalise().

◆ m_idx

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

Number of activations and repolarisations detected for each point.

Definition at line 92 of file FilterBenchmark.h.

Referenced by v_Initialise(), and v_Update().

◆ m_initialValue

NekDouble Nektar::FilterBenchmark::m_initialValue
private

Initial time to use in storage array.

Definition at line 100 of file FilterBenchmark.h.

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

◆ m_outputFile

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

Filename of output files.

Definition at line 102 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Finalise().

◆ m_polarity

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

Indicates if the previous event was an activation or repolarisation.

Definition at line 94 of file FilterBenchmark.h.

Referenced by v_Initialise(), and v_Update().

◆ m_startTime

NekDouble Nektar::FilterBenchmark::m_startTime
private

Time at which to start detecting activations and repolarisations.

Definition at line 96 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Update().

◆ m_threshold

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

Storage for activation and repolarisation times.

Definition at line 90 of file FilterBenchmark.h.

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

◆ m_thresholdValue

NekDouble Nektar::FilterBenchmark::m_thresholdValue
private

Value at which tissue is considered active.

Definition at line 98 of file FilterBenchmark.h.

Referenced by FilterBenchmark(), and v_Update().