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

#include <FilterFieldConvert.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterFieldConvert (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterFieldConvert () 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
 
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
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
 
virtual SOLVER_UTILS_EXPORT void v_ProcessSample (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetScale ()
 
virtual SOLVER_UTILS_EXPORT std::string v_GetFileSuffix ()
 
void OutputField (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int dump=-1)
 
SOLVER_UTILS_EXPORT bool v_IsTimeDependent () override
 
void CreateModules (std::vector< std::string > &modcmds)
 
void CreateFields (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
void CheckModules (std::vector< ModuleSharedPtr > &modules)
 
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

unsigned int m_numSamples
 
unsigned int m_outputFrequency
 
unsigned int m_sampleFrequency
 
std::string m_outputFile
 
std::string m_restartFile
 
unsigned int m_index
 
unsigned int m_outputIndex
 
bool m_phaseSample
 
NekDouble m_phaseSamplePeriod
 
NekDouble m_phaseSamplePhase
 
NekDouble m_phaseTolerance
 
NekDouble m_dt
 
std::vector< ModuleSharedPtrm_modules
 
LibUtilities::FieldMetaDataMap m_fieldMetaData
 
std::vector< Array< OneD, NekDouble > > m_outFields
 
std::vector< std::string > m_variables
 
FieldSharedPtr m_f
 
po::variables_map m_vm
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Friends

class MemoryManager< FilterFieldConvert >
 

Additional Inherited Members

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

Detailed Description

Definition at line 46 of file FilterFieldConvert.h.

Constructor & Destructor Documentation

◆ FilterFieldConvert()

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

Definition at line 48 of file FilterFieldConvert.cpp.

51 : Filter(pSession, pEquation)
52{
53 m_dt = m_session->GetParameter("TimeStep");
54
55 // OutputFile
56 auto it = pParams.find("OutputFile");
57 if (it == pParams.end())
58 {
59 std::stringstream outname;
60 outname << m_session->GetSessionName() << ".fld";
61 m_outputFile = outname.str();
62 }
63 else
64 {
65 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
66 if (it->second.find_last_of('.') != std::string::npos)
67 {
68 m_outputFile = it->second;
69 }
70 else
71 {
72 std::stringstream outname;
73 outname << it->second << ".fld";
74 m_outputFile = outname.str();
75 }
76 }
77
78 // Restart file
79 it = pParams.find("RestartFile");
80 if (it == pParams.end())
81 {
82 m_restartFile = "";
83 }
84 else
85 {
86 ASSERTL0(it->second.length() > 0, "Missing parameter 'RestartFile'.");
87 if (it->second.find_last_of('.') != std::string::npos)
88 {
89 m_restartFile = it->second;
90 }
91 else
92 {
93 std::stringstream outname;
94 outname << it->second << ".fld";
95 m_restartFile = outname.str();
96 }
97 }
98
99 // OutputFrequency
100 it = pParams.find("OutputFrequency");
101 if (it == pParams.end())
102 {
103 m_outputFrequency = m_session->GetParameter("NumSteps");
104 }
105 else
106 {
107 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
108 m_outputFrequency = round(equ.Evaluate());
109 }
110
111 // The base class can use SampleFrequency = OutputFrequency
112 // (Derived classes need to override this if needed)
114
115 // Phase sampling option
116 it = pParams.find("PhaseAverage");
117 if (it == pParams.end())
118 {
119 m_phaseSample = false;
120 }
121 else
122 {
123 std::string sOption = it->second.c_str();
124 m_phaseSample = (boost::iequals(sOption, "true")) ||
125 (boost::iequals(sOption, "yes"));
126 }
127
128 if (m_phaseSample)
129 {
130 auto itPeriod = pParams.find("PhaseAveragePeriod");
131 auto itPhase = pParams.find("PhaseAveragePhase");
132
133 // Error if only one of the required params for PhaseAverage is present
134 ASSERTL0(
135 (itPeriod != pParams.end() && itPhase != pParams.end()),
136 "The phase sampling feature requires both 'PhaseAveragePeriod' and "
137 "'PhaseAveragePhase' to be set.");
138
139 LibUtilities::Equation equPeriod(m_session->GetInterpreter(),
140 itPeriod->second);
141 m_phaseSamplePeriod = equPeriod.Evaluate();
142
143 LibUtilities::Equation equPhase(m_session->GetInterpreter(),
144 itPhase->second);
145 m_phaseSamplePhase = equPhase.Evaluate();
146
147 // Check that phase and period are within required limits
149 "PhaseAveragePeriod must be greater than 0.");
151 "PhaseAveragePhase must be between 0 and 1.");
152
153 // Load sampling frequency, overriding the previous value
154 it = pParams.find("SampleFrequency");
155 if (it == pParams.end())
156 {
158 }
159 else
160 {
161 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
162 m_sampleFrequency = round(equ.Evaluate());
163 }
164
165 // Compute tolerance within which sampling occurs.
167
168 // Display worst case scenario sampling tolerance for exact phase, if
169 // verbose option is active
170 if (m_session->GetComm()->GetRank() == 0 &&
171 m_session->DefinesCmdLineArgument("verbose"))
172 {
173 std::cout << "Phase sampling activated with period "
174 << m_phaseSamplePeriod << " and phase "
175 << m_phaseSamplePhase << "." << std::endl
176 << "Sampling within a tolerance of "
177 << std::setprecision(6) << m_phaseTolerance << "."
178 << std::endl;
179 }
180 }
181
182 m_numSamples = 0;
183 m_index = 0;
184 m_outputIndex = 0;
185
186 //
187 // FieldConvert modules
188 //
189 m_f = std::shared_ptr<Field>(new Field());
190 std::vector<std::string> modcmds;
191 // Process modules
192 std::stringstream moduleStream;
193 it = pParams.find("Modules");
194 if (it != pParams.end())
195 {
196 moduleStream.str(it->second);
197 }
198 while (!moduleStream.fail())
199 {
200 std::string sMod;
201 moduleStream >> sMod;
202 if (!moduleStream.fail())
203 {
204 modcmds.push_back(sMod);
205 }
206 }
207 // Output module
208 modcmds.push_back(m_outputFile);
209 // Create modules
210 CreateModules(modcmds);
211 // Strip options from m_outputFile
212 std::vector<std::string> tmp;
213 boost::split(tmp, m_outputFile, boost::is_any_of(":"));
214 m_outputFile = tmp[0];
215
216 // Prevent checking before overwriting
217 it = pParams.find("options");
218 if (it != pParams.end())
219 {
220 int argc = 0;
221 std::string strargv;
222 std::vector<char *> argv;
223
224 strargv = "dummy " + it->second;
225 strargv.push_back((char)0);
226 bool flag = true;
227 for (size_t i = 0; strargv[i]; ++i)
228 {
229 if (strargv[i] != ' ' && flag)
230 {
231 argv.push_back(&strargv[i]);
232 ++argc;
233 flag = false;
234 }
235 if (strargv[i] == ' ')
236 {
237 flag = true;
238 strargv[i] = 0;
239 }
240 }
241
242 po::options_description desc("Available options");
243
244 // clang-format off
245 desc.add_options()
246 ("help,h",
247 "Produce this help message.")
248 ("modules-list,l",
249 "Print the list of available modules.")
250 ("output-points,n", po::value<int>(),
251 "Output at n equipspaced points along the "
252 "collapsed coordinates (for .dat, .vtu).")
253 ("output-points-hom-z", po::value<int>(),
254 "Number of planes in the z-direction for output of "
255 "Homogeneous 1D expansion(for .dat, .vtu).")
256 ("error,e",
257 "Write error of fields for regression checking")
258 ("force-output,f",
259 "Force the output to be written without any checks")
260 ("range,r", po::value<std::string>(),
261 "Define output range i.e. (-r xmin,xmax,ymin,ymax,zmin,zmax) "
262 "in which any vertex is contained.")
263 ("no-equispaced",
264 "Do not use equispaced output.")
265 ("nparts", po::value<int>(),
266 "Define nparts if running serial problem to mimic "
267 "parallel run with many partitions.")
268 ("npz", po::value<int>(),
269 "Used to define number of partitions in z for Homogeneous1D "
270 "expansions for parallel runs.")
271 ("npt", po::value<int>(),
272 "Used to define number of partitions in time for Parareal runs. ")
273 ("onlyshape", po::value<std::string>(),
274 "Only use element with defined shape type i.e. -onlyshape "
275 " Tetrahedron")
276 ("part-only", po::value<int>(),
277 "Partition into specified npart partitions and exit")
278 ("part-only-overlapping", po::value<int>(),
279 "Partition into specified npart overlapping partitions and exit")
280 ("modules-opt,p", po::value<std::string>(),
281 "Print options for a module.")
282 ("module,m", po::value<std::vector<std::string> >(),
283 "Specify modules which are to be used.")
284 ("use-session-variables",
285 "Use variables defined in session for output")
286 ("use-session-expansion",
287 "Use expansion defined in session.")
288 ("verbose,v",
289 "Enable verbose mode.");
290 // clang-format on
291
292 po::options_description hidden("Hidden options");
293
294 // clang-format off
295 hidden.add_options()
296 ("input-file", po::value<std::vector<std::string> >(),
297 "Input filename");
298 // clang-format on
299
300 po::options_description cmdline_options;
301 cmdline_options.add(hidden).add(desc);
302
303 po::options_description visible("Allowed options");
304 visible.add(desc);
305
306 po::positional_options_description p;
307 p.add("input-file", -1);
308
309 try
310 {
311 po::store(po::command_line_parser(argc, &(argv[0]))
312 .options(cmdline_options)
313 .positional(p)
314 .run(),
315 m_vm);
316 po::notify(m_vm);
317 }
318 catch (const std::exception &e)
319 {
320 std::cerr << e.what() << std::endl;
321 std::cerr << desc;
322 }
323 }
324 m_vm.insert(std::make_pair("force-output", po::variable_value()));
325}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void CreateModules(std::vector< std::string > &modcmds)
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, CreateModules(), Nektar::LibUtilities::Equation::Evaluate(), m_dt, m_f, m_index, m_numSamples, m_outputFile, m_outputFrequency, m_outputIndex, m_phaseSample, m_phaseSamplePeriod, m_phaseSamplePhase, m_phaseTolerance, m_restartFile, m_sampleFrequency, Nektar::SolverUtils::Filter::m_session, m_vm, CellMLToNektar.cellml_metadata::p, and CellMLToNektar.translators::run().

◆ ~FilterFieldConvert()

Nektar::SolverUtils::FilterFieldConvert::~FilterFieldConvert ( )
override

Definition at line 327 of file FilterFieldConvert.cpp.

328{
329}

Member Function Documentation

◆ CheckModules()

void Nektar::SolverUtils::FilterFieldConvert::CheckModules ( std::vector< ModuleSharedPtr > &  modules)
protected

Definition at line 730 of file FilterFieldConvert.cpp.

731{
732 // Count number of modules by priority
733 Array<OneD, int> modulesCount(SIZE_ModulePriority, 0);
734 for (int i = 0; i < modules.size(); ++i)
735 {
736 ++modulesCount[modules[i]->GetModulePriority()];
737 }
738
739 // FilterFieldConvert already starts with m_exp, so anything before
740 // eModifyExp is not valid, and also eCreatePts
741 if (modulesCount[eCreateGraph] != 0 ||
742 modulesCount[eCreateFieldData] != 0 ||
743 modulesCount[eModifyFieldData] != 0 || modulesCount[eCreateExp] != 0 ||
744 modulesCount[eFillExp] != 0 || modulesCount[eCreatePts] != 0)
745 {
746 std::stringstream ss;
747 ss << "Module(s): ";
748 for (int i = 0; i < modules.size(); ++i)
749 {
750 if (modules[i]->GetModulePriority() == eCreateGraph ||
751 modules[i]->GetModulePriority() == eCreateFieldData ||
752 modules[i]->GetModulePriority() == eModifyFieldData ||
753 modules[i]->GetModulePriority() == eCreateExp ||
754 modules[i]->GetModulePriority() == eFillExp ||
755 modules[i]->GetModulePriority() == eCreatePts)
756 {
757 ss << modules[i]->GetModuleName() << " ";
758 }
759 }
760 ss << "not compatible with FilterFieldConvert.";
761 ASSERTL0(false, ss.str());
762 }
763
764 // Modules of type eConvertExpToPts are not compatible with eBndExtraction
765 if (modulesCount[eConvertExpToPts] != 0 &&
766 modulesCount[eBndExtraction] != 0)
767 {
768 std::stringstream ss;
769 ss << "Module(s): ";
770 for (int i = 0; i < modules.size(); ++i)
771 {
772 if (modules[i]->GetModulePriority() == eBndExtraction)
773 {
774 ss << modules[i]->GetModuleName() << " ";
775 }
776 }
777 ss << "is not compatible with module(s): ";
778 for (int i = 0; i < modules.size(); ++i)
779 {
780 if (modules[i]->GetModulePriority() == eConvertExpToPts)
781 {
782 ss << modules[i]->GetModuleName() << " ";
783 }
784 }
785 ss << ".";
786 ASSERTL0(false, ss.str());
787 }
788}

References ASSERTL0, Nektar::FieldUtils::eBndExtraction, Nektar::FieldUtils::eConvertExpToPts, Nektar::FieldUtils::eCreateExp, Nektar::FieldUtils::eCreateFieldData, Nektar::FieldUtils::eCreateGraph, Nektar::FieldUtils::eCreatePts, Nektar::FieldUtils::eFillExp, Nektar::FieldUtils::eModifyFieldData, and Nektar::FieldUtils::SIZE_ModulePriority.

Referenced by CreateModules().

◆ create()

static FilterSharedPtr Nektar::SolverUtils::FilterFieldConvert::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 52 of file FilterFieldConvert.h.

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

◆ CreateFields()

void Nektar::SolverUtils::FilterFieldConvert::CreateFields ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields)
protected

Definition at line 684 of file FilterFieldConvert.cpp.

686{
687 // Fill some parameters of m_f
688 m_f->m_session = m_session;
689 m_f->m_graph = pFields[0]->GetGraph();
690 m_f->m_comm = m_f->m_session->GetComm();
691 m_f->m_fieldMetaDataMap = m_fieldMetaData;
692 m_f->m_fieldPts = LibUtilities::NullPtsField;
693 // Create m_f->m_exp
694 m_f->m_numHomogeneousDir = 0;
695 if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
696 {
697 m_f->m_numHomogeneousDir = 1;
698 }
699 else if (pFields[0]->GetExpType() == MultiRegions::e3DH2D)
700 {
701 m_f->m_numHomogeneousDir = 2;
702 }
703
704 m_f->m_exp.resize(m_variables.size());
705 m_f->m_exp[0] = pFields[0];
706 int nfield;
707 for (int n = 0; n < m_variables.size(); ++n)
708 {
709 // if n >= pFields.size() assume we have used n=0 field
710 nfield = (n < pFields.size()) ? n : 0;
711
712 m_f->m_exp[n] =
713 m_f->AppendExpList(m_f->m_numHomogeneousDir, m_variables[0]);
714 m_f->m_exp[n]->SetWaveSpace(false);
715
716 ASSERTL1(pFields[nfield]->GetNcoeffs() == m_outFields[n].size(),
717 "pFields[nfield] does not have the "
718 "same number of coefficients as m_outFields[n]");
719
720 m_f->m_exp[n]->ExtractCoeffsToCoeffs(pFields[nfield], m_outFields[n],
721 m_f->m_exp[n]->UpdateCoeffs());
722
723 m_f->m_exp[n]->BwdTrans(m_f->m_exp[n]->GetCoeffs(),
724 m_f->m_exp[n]->UpdatePhys());
725 }
726 m_f->m_variables = m_variables;
727}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
LibUtilities::FieldMetaDataMap m_fieldMetaData
std::vector< Array< OneD, NekDouble > > m_outFields
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:185

References ASSERTL1, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_f, m_fieldMetaData, m_outFields, Nektar::SolverUtils::Filter::m_session, m_variables, and Nektar::LibUtilities::NullPtsField.

Referenced by OutputField().

◆ CreateModules()

void Nektar::SolverUtils::FilterFieldConvert::CreateModules ( std::vector< std::string > &  modcmds)
protected

Definition at line 589 of file FilterFieldConvert.cpp.

590{
591 for (int i = 0; i < modcmds.size(); ++i)
592 {
593 // First split each command by the colon separator.
594 std::vector<std::string> tmp1;
595 ModuleKey module;
596 int offset = 1;
597
598 boost::split(tmp1, modcmds[i], boost::is_any_of(":"));
599
600 if (i == modcmds.size() - 1)
601 {
602 module.first = eOutputModule;
603
604 // If no colon detected, automatically detect mesh type from
605 // file extension. Otherwise override and use tmp1[1] as the
606 // module to load. This also allows us to pass options to
607 // input/output modules. So, for example, to override
608 // filename.xml to be read as vtk, you use:
609 //
610 // filename.xml:vtk:opt1=arg1:opt2=arg2
611 if (tmp1.size() == 1)
612 {
613 int dot = tmp1[0].find_last_of('.') + 1;
614 std::string ext = tmp1[0].substr(dot, tmp1[0].length() - dot);
615 module.second = ext;
616 tmp1.push_back(std::string("outfile=") + tmp1[0]);
617 }
618 else
619 {
620 module.second = tmp1[1];
621 tmp1.push_back(std::string("outfile=") + tmp1[0]);
622 offset++;
623 }
624 }
625 else
626 {
627 module.first = eProcessModule;
628 module.second = tmp1[0];
629 }
630
631 // Create modules
632 ModuleSharedPtr mod;
633 mod = GetModuleFactory().CreateInstance(module, m_f);
634 m_modules.push_back(mod);
635
636 // Set options for this module.
637 for (int j = offset; j < tmp1.size(); ++j)
638 {
639 std::vector<std::string> tmp2;
640 boost::split(tmp2, tmp1[j], boost::is_any_of("="));
641
642 if (tmp2.size() == 1)
643 {
644 mod->RegisterConfig(tmp2[0]);
645 }
646 else if (tmp2.size() == 2)
647 {
648 mod->RegisterConfig(tmp2[0], tmp2[1]);
649 }
650 else
651 {
652 std::cerr << "ERROR: Invalid module configuration: format is "
653 << "either :arg or :arg=val" << std::endl;
654 abort();
655 }
656 }
657
658 // Ensure configuration options have been set.
659 mod->SetDefaults();
660 }
661
662 // Include equispaced output if needed
663 Array<OneD, int> modulesCount(SIZE_ModulePriority, 0);
664 for (int i = 0; i < m_modules.size(); ++i)
665 {
666 ++modulesCount[m_modules[i]->GetModulePriority()];
667 }
668 if (modulesCount[eModifyPts] != 0 && modulesCount[eCreatePts] == 0 &&
669 modulesCount[eConvertExpToPts] == 0)
670 {
671 ModuleKey module;
672 ModuleSharedPtr mod;
673 module.first = eProcessModule;
674 module.second = std::string("equispacedoutput");
675 mod = GetModuleFactory().CreateInstance(module, m_f);
676 m_modules.insert(m_modules.end() - 1, mod);
677 mod->SetDefaults();
678 }
679
680 // Check if modules provided are compatible
682}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void CheckModules(std::vector< ModuleSharedPtr > &modules)
std::vector< ModuleSharedPtr > m_modules
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
std::shared_ptr< Module > ModuleSharedPtr
Definition: Module.h:329
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

References CheckModules(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::FieldUtils::eConvertExpToPts, Nektar::FieldUtils::eCreatePts, Nektar::FieldUtils::eModifyPts, Nektar::FieldUtils::eOutputModule, Nektar::FieldUtils::eProcessModule, Nektar::FieldUtils::GetModuleFactory(), m_f, m_modules, and Nektar::FieldUtils::SIZE_ModulePriority.

Referenced by FilterFieldConvert().

◆ OutputField()

void Nektar::SolverUtils::FilterFieldConvert::OutputField ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
int  dump = -1 
)
protected

Definition at line 530 of file FilterFieldConvert.cpp.

532{
533 NekDouble scale = v_GetScale();
534 for (int n = 0; n < m_outFields.size(); ++n)
535 {
536 Vmath::Smul(m_outFields[n].size(), scale, m_outFields[n], 1,
537 m_outFields[n], 1);
538 }
539
540 CreateFields(pFields);
541
542 // Determine new file name
543 std::stringstream outname;
544 int dot = m_outputFile.find_last_of('.');
545 std::string name = m_outputFile.substr(0, dot);
546 std::string ext = m_outputFile.substr(dot, m_outputFile.length() - dot);
547 std::string suffix = v_GetFileSuffix();
548 if (dump == -1) // final dump
549 {
550 outname << name << suffix << ext;
551 }
552 else
553 {
554 outname << name << "_" << dump << suffix << ext;
555 }
556 m_modules[m_modules.size() - 1]->RegisterConfig("outfile", outname.str());
557
558 // Run field process.
559 for (int n = 0; n < SIZE_ModulePriority; ++n)
560 {
561 ModulePriority priority = static_cast<ModulePriority>(n);
562 for (int i = 0; i < m_modules.size(); ++i)
563 {
564 if (m_modules[i]->GetModulePriority() == priority)
565 {
566 m_modules[i]->Process(m_vm);
567 }
568 }
569 }
570
571 // Empty m_f to save memory
572 m_f->ClearField();
573
574 if (dump != -1) // not final dump so rescale
575 {
576 for (int n = 0; n < m_outFields.size(); ++n)
577 {
578 Vmath::Smul(m_outFields[n].size(), 1.0 / scale, m_outFields[n], 1,
579 m_outFields[n], 1);
580 }
581 }
582}
void CreateFields(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
virtual SOLVER_UTILS_EXPORT NekDouble v_GetScale()
virtual SOLVER_UTILS_EXPORT std::string v_GetFileSuffix()
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References CreateFields(), m_f, m_modules, m_outFields, m_outputFile, m_vm, CellMLToNektar.pycml::name, Nektar::FieldUtils::SIZE_ModulePriority, Vmath::Smul(), v_GetFileSuffix(), and v_GetScale().

Referenced by v_Finalise(), and v_Update().

◆ v_FillVariablesName()

void Nektar::SolverUtils::FilterFieldConvert::v_FillVariablesName ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields)
protectedvirtual

Reimplemented in Nektar::SolverUtils::FilterAverageFields, Nektar::SolverUtils::FilterBodyFittedVelocity, Nektar::SolverUtils::FilterMaxMinFields, and Nektar::SolverUtils::FilterReynoldsStresses.

Definition at line 420 of file FilterFieldConvert.cpp.

422{
423 int nfield = pFields.size();
424 m_variables.resize(pFields.size());
425 for (int n = 0; n < nfield; ++n)
426 {
427 m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
428 }
429
430 // Need to create a dummy coeffs vector to get extra variables names...
431 std::vector<Array<OneD, NekDouble>> coeffs(nfield);
432 for (int n = 0; n < nfield; ++n)
433 {
434 coeffs[n] = pFields[n]->GetCoeffs();
435 }
436 // Get extra variables
437 auto equ = m_equ.lock();
438 ASSERTL0(equ, "Weak pointer expired");
439 equ->ExtraFldOutput(coeffs, m_variables);
440}
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84

References ASSERTL0, Nektar::SolverUtils::Filter::m_equ, and m_variables.

Referenced by Nektar::SolverUtils::FilterAverageFields::v_FillVariablesName(), Nektar::SolverUtils::FilterBodyFittedVelocity::v_FillVariablesName(), Nektar::SolverUtils::FilterMaxMinFields::v_FillVariablesName(), and v_Initialise().

◆ v_Finalise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 508 of file FilterFieldConvert.cpp.

511{
512 m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
513 v_PrepareOutput(pFields, time);
514 OutputField(pFields);
515}
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
void OutputField(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int dump=-1)

References m_fieldMetaData, OutputField(), and v_PrepareOutput().

◆ v_GetFileSuffix()

virtual SOLVER_UTILS_EXPORT std::string Nektar::SolverUtils::FilterFieldConvert::v_GetFileSuffix ( )
inlineprotectedvirtual

◆ v_GetScale()

virtual SOLVER_UTILS_EXPORT NekDouble Nektar::SolverUtils::FilterFieldConvert::v_GetScale ( )
inlineprotectedvirtual

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Reimplemented in Nektar::SolverUtils::FilterMaxMinFields, and Nektar::SolverUtils::FilterReynoldsStresses.

Definition at line 331 of file FilterFieldConvert.cpp.

334{
335 v_FillVariablesName(pFields);
336
337 // m_variables need to be filled by a derived class
338 m_outFields.resize(m_variables.size());
339 int nfield;
340
341 for (int n = 0; n < m_variables.size(); ++n)
342 {
343 // if n >= pFields.size() assum we have used n=0 field
344 nfield = (n < pFields.size()) ? n : 0;
345
346 m_outFields[n] =
347 Array<OneD, NekDouble>(pFields[nfield]->GetNcoeffs(), 0.0);
348 }
349
350 m_fieldMetaData["InitialTime"] = boost::lexical_cast<std::string>(time);
351
352 // Load restart file if necessary
353 if (m_restartFile != "")
354 {
355 // Load file
356 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
357 std::vector<std::vector<NekDouble>> fieldData;
358 LibUtilities::FieldMetaDataMap fieldMetaData;
361 fld->Import(m_restartFile, fieldDef, fieldData, fieldMetaData);
362
363 // Extract fields to output
364 int nfield = -1, k;
365 for (int j = 0; j < m_variables.size(); ++j)
366 {
367 // see if m_variables is part of pFields definition and if
368 // so use that field for extract
369 for (k = 0; k < pFields.size(); ++k)
370 {
371 if (pFields[k]->GetSession()->GetVariable(k) == m_variables[j])
372 {
373 nfield = k;
374 break;
375 }
376 }
377 if (nfield == -1)
378 {
379 nfield = 0;
380 }
381
382 for (int i = 0; i < fieldData.size(); ++i)
383 {
384 pFields[nfield]->ExtractDataToCoeffs(
385 fieldDef[i], fieldData[i], m_variables[j], m_outFields[j]);
386 }
387 }
388
389 // Load information for numSamples
390 if (fieldMetaData.count("NumberOfFieldDumps"))
391 {
392 m_numSamples = atoi(fieldMetaData["NumberOfFieldDumps"].c_str());
393 }
394 else
395 {
396 m_numSamples = 1;
397 }
398
399 if (fieldMetaData.count("InitialTime"))
400 {
401 m_fieldMetaData["InitialTime"] = fieldMetaData["InitialTime"];
402 }
403
404 // Load information for outputIndex
405 if (fieldMetaData.count("FilterFileNum"))
406 {
407 m_outputIndex = atoi(fieldMetaData["FilterFileNum"].c_str());
408 }
409
410 // Divide by scale
411 NekDouble scale = v_GetScale();
412 for (int n = 0; n < m_outFields.size(); ++n)
413 {
414 Vmath::Smul(m_outFields[n].size(), 1.0 / scale, m_outFields[n], 1,
415 m_outFields[n], 1);
416 }
417 }
418}
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:224
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50

References Nektar::LibUtilities::FieldIO::CreateForFile(), m_fieldMetaData, m_numSamples, m_outFields, m_outputIndex, m_restartFile, Nektar::SolverUtils::Filter::m_session, m_variables, Vmath::Smul(), v_FillVariablesName(), and v_GetScale().

Referenced by Nektar::SolverUtils::FilterBodyFittedVelocity::v_Initialise(), Nektar::SolverUtils::FilterMaxMinFields::v_Initialise(), and Nektar::SolverUtils::FilterReynoldsStresses::v_Initialise().

◆ v_IsTimeDependent()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 584 of file FilterFieldConvert.cpp.

585{
586 return true;
587}

◆ v_PrepareOutput()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::FilterFieldConvert::v_PrepareOutput ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
inlineprotectedvirtual

◆ v_ProcessSample()

void Nektar::SolverUtils::FilterFieldConvert::v_ProcessSample ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
const NekDouble time 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::FilterAverageFields, Nektar::SolverUtils::FilterBodyFittedVelocity, Nektar::SolverUtils::FilterMaxMinFields, Nektar::SolverUtils::FilterMovingAverage, and Nektar::SolverUtils::FilterReynoldsStresses.

Definition at line 517 of file FilterFieldConvert.cpp.

522{
523 for (int n = 0; n < m_outFields.size(); ++n)
524 {
525 Vmath::Vcopy(m_outFields[n].size(), fieldcoeffs[n], 1, m_outFields[n],
526 1);
527 }
528}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References m_outFields, and Vmath::Vcopy().

Referenced by v_Update().

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 442 of file FilterFieldConvert.cpp.

445{
446 m_index++;
447 if (m_index % m_sampleFrequency > 0)
448 {
449 return;
450 }
451
452 // Append extra fields
453 int nfield = pFields.size();
454 std::vector<Array<OneD, NekDouble>> coeffs(nfield);
455 for (int n = 0; n < nfield; ++n)
456 {
457 coeffs[n] = pFields[n]->GetCoeffs();
458 }
459 std::vector<std::string> variables = m_variables;
460 auto equ = m_equ.lock();
461 ASSERTL0(equ, "Weak pointer expired");
462 equ->ExtraFldOutput(coeffs, variables);
463
464 if (m_phaseSample)
465 {
466 // The sample is added to the filter only if the current time
467 // corresponds to the correct phase. Introducing M as number of
468 // cycles and N nondimensional phase (between 0 and 1):
469 // t = M * m_phaseSamplePeriod + N * m_phaseSamplePeriod
470 int currentCycle = floor(time / m_phaseSamplePeriod);
471 NekDouble currentPhase = time / m_phaseSamplePeriod - currentCycle;
472
473 // Evaluate phase relative to the requested value.
474 NekDouble relativePhase = fabs(m_phaseSamplePhase - currentPhase);
475
476 // Check if relative phase is within required tolerance and sample.
477 // Care must be taken to handle the cases at phase 0 as the sample might
478 // have to be taken at the very end of the previous cycle instead.
479 if (relativePhase < m_phaseTolerance ||
480 fabs(relativePhase - 1) < m_phaseTolerance)
481 {
482 m_numSamples++;
483 v_ProcessSample(pFields, coeffs, time);
484 if (m_session->GetComm()->GetRank() == 0 &&
485 m_session->DefinesCmdLineArgument("verbose"))
486 {
487 std::cout << "Sample: " << std::setw(8) << std::left
488 << m_numSamples << "Phase: " << std::setw(8)
489 << std::left << currentPhase << std::endl;
490 }
491 }
492 }
493 else
494 {
495 m_numSamples++;
496 v_ProcessSample(pFields, coeffs, time);
497 }
498
499 if (m_index % m_outputFrequency == 0)
500 {
501 m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
502 v_PrepareOutput(pFields, time);
503 m_fieldMetaData["FilterFileNum"] = std::to_string(++m_outputIndex);
504 OutputField(pFields, m_outputIndex);
505 }
506}
virtual SOLVER_UTILS_EXPORT void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)

References ASSERTL0, Nektar::SolverUtils::Filter::m_equ, m_fieldMetaData, m_index, m_numSamples, m_outputFrequency, m_outputIndex, m_phaseSample, m_phaseSamplePeriod, m_phaseSamplePhase, m_phaseTolerance, m_sampleFrequency, Nektar::SolverUtils::Filter::m_session, m_variables, OutputField(), v_PrepareOutput(), and v_ProcessSample().

Friends And Related Function Documentation

◆ MemoryManager< FilterFieldConvert >

friend class MemoryManager< FilterFieldConvert >
friend

Definition at line 1 of file FilterFieldConvert.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterFieldConvert::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 64 of file FilterFieldConvert.h.

◆ m_dt

NekDouble Nektar::SolverUtils::FilterFieldConvert::m_dt
protected

Definition at line 129 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert().

◆ m_f

FieldSharedPtr Nektar::SolverUtils::FilterFieldConvert::m_f
protected

Definition at line 135 of file FilterFieldConvert.h.

Referenced by CreateFields(), CreateModules(), FilterFieldConvert(), and OutputField().

◆ m_fieldMetaData

LibUtilities::FieldMetaDataMap Nektar::SolverUtils::FilterFieldConvert::m_fieldMetaData
protected

◆ m_index

unsigned int Nektar::SolverUtils::FilterFieldConvert::m_index
protected

Definition at line 121 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_modules

std::vector<ModuleSharedPtr> Nektar::SolverUtils::FilterFieldConvert::m_modules
protected

Definition at line 131 of file FilterFieldConvert.h.

Referenced by CreateModules(), and OutputField().

◆ m_numSamples

unsigned int Nektar::SolverUtils::FilterFieldConvert::m_numSamples
protected

◆ m_outFields

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterFieldConvert::m_outFields
protected

◆ m_outputFile

std::string Nektar::SolverUtils::FilterFieldConvert::m_outputFile
protected

Definition at line 119 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and OutputField().

◆ m_outputFrequency

unsigned int Nektar::SolverUtils::FilterFieldConvert::m_outputFrequency
protected

Definition at line 117 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_outputIndex

unsigned int Nektar::SolverUtils::FilterFieldConvert::m_outputIndex
protected

Definition at line 122 of file FilterFieldConvert.h.

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

◆ m_phaseSample

bool Nektar::SolverUtils::FilterFieldConvert::m_phaseSample
protected

Definition at line 125 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_phaseSamplePeriod

NekDouble Nektar::SolverUtils::FilterFieldConvert::m_phaseSamplePeriod
protected

Definition at line 126 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_phaseSamplePhase

NekDouble Nektar::SolverUtils::FilterFieldConvert::m_phaseSamplePhase
protected

Definition at line 127 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_phaseTolerance

NekDouble Nektar::SolverUtils::FilterFieldConvert::m_phaseTolerance
protected

Definition at line 128 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and v_Update().

◆ m_restartFile

std::string Nektar::SolverUtils::FilterFieldConvert::m_restartFile
protected

◆ m_sampleFrequency

unsigned int Nektar::SolverUtils::FilterFieldConvert::m_sampleFrequency
protected

◆ m_variables

std::vector<std::string> Nektar::SolverUtils::FilterFieldConvert::m_variables
protected

◆ m_vm

po::variables_map Nektar::SolverUtils::FilterFieldConvert::m_vm
protected

Definition at line 136 of file FilterFieldConvert.h.

Referenced by FilterFieldConvert(), and OutputField().