Nektar++
FilterFieldConvert.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File FilterFieldConvert.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Base clase for filters performing operations on samples
32 // of the field.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <boost/core/ignore_unused.hpp>
38 #include <boost/program_options.hpp>
39 #include <boost/algorithm/string/classification.hpp>
40 #include <boost/algorithm/string/split.hpp>
41 #include <boost/algorithm/string/predicate.hpp>
42 
43 namespace Nektar
44 {
45 namespace SolverUtils
46 {
49  "FieldConvert", FilterFieldConvert::create);
50 
53  const std::weak_ptr<EquationSystem> &pEquation,
54  const ParamMap &pParams)
55  : Filter(pSession, pEquation)
56 {
57  m_dt = m_session->GetParameter("TimeStep");
58 
59  // OutputFile
60  auto it = pParams.find("OutputFile");
61  if (it == pParams.end())
62  {
63  std::stringstream outname;
64  outname << m_session->GetSessionName() << ".fld";
65  m_outputFile = outname.str();
66  }
67  else
68  {
69  ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
70  if ( it->second.find_last_of('.') != std::string::npos)
71  {
72  m_outputFile = it->second;
73  }
74  else
75  {
76  std::stringstream outname;
77  outname << it->second << ".fld";
78  m_outputFile = outname.str();
79  }
80  }
81 
82  // Restart file
83  it = pParams.find("RestartFile");
84  if (it == pParams.end())
85  {
86  m_restartFile = "";
87  }
88  else
89  {
90  ASSERTL0(it->second.length() > 0, "Missing parameter 'RestartFile'.");
91  if ( it->second.find_last_of('.') != std::string::npos)
92  {
93  m_restartFile = it->second;
94  }
95  else
96  {
97  std::stringstream outname;
98  outname << it->second << ".fld";
99  m_restartFile = outname.str();
100  }
101  }
102 
103  // OutputFrequency
104  it = pParams.find("OutputFrequency");
105  if (it == pParams.end())
106  {
107  m_outputFrequency = m_session->GetParameter("NumSteps");
108  }
109  else
110  {
112  m_session->GetInterpreter(), it->second);
113  m_outputFrequency = round(equ.Evaluate());
114  }
115 
116  // The base class can use SampleFrequency = OutputFrequency
117  // (Derived classes need to override this if needed)
119 
120  // Phase sampling option
121  it = pParams.find("PhaseAverage");
122  if (it == pParams.end())
123  {
124  m_phaseSample = false;
125  }
126  else
127  {
128  std::string sOption = it->second.c_str();
129  m_phaseSample = (boost::iequals(sOption, "true")) ||
130  (boost::iequals(sOption, "yes"));
131  }
132 
133  if(m_phaseSample)
134  {
135  auto itPeriod = pParams.find("PhaseAveragePeriod");
136  auto itPhase = pParams.find("PhaseAveragePhase");
137 
138  // Error if only one of the required params for PhaseAverage is present
139  ASSERTL0((itPeriod != pParams.end() && itPhase != pParams.end()),
140  "The phase sampling feature requires both 'PhaseAveragePeriod' and "
141  "'PhaseAveragePhase' to be set.");
142 
143  LibUtilities::Equation equPeriod(
144  m_session->GetInterpreter(), itPeriod->second);
145  m_phaseSamplePeriod = equPeriod.Evaluate();
146 
147  LibUtilities::Equation equPhase(
148  m_session->GetInterpreter(), itPhase->second);
149  m_phaseSamplePhase = equPhase.Evaluate();
150 
151  // Check that phase and period are within required limits
152  ASSERTL0(m_phaseSamplePeriod > 0,
153  "PhaseAveragePeriod must be greater than 0.");
154  ASSERTL0(m_phaseSamplePhase >= 0 && m_phaseSamplePhase <= 1,
155  "PhaseAveragePhase must be between 0 and 1.");
156 
157  // Load sampling frequency, overriding the previous value
158  it = pParams.find("SampleFrequency");
159  if (it == pParams.end())
160  {
161  m_sampleFrequency = 1;
162  }
163  else
164  {
166  m_session->GetInterpreter(), it->second);
167  m_sampleFrequency = round(equ.Evaluate());
168  }
169 
170  // Compute tolerance within which sampling occurs.
171  m_phaseTolerance = m_dt * m_sampleFrequency /
172  (m_phaseSamplePeriod * 2);
173 
174  // Display worst case scenario sampling tolerance for exact phase, if
175  // verbose option is active
176  if (m_session->GetComm()->GetRank() == 0 &&
177  m_session->DefinesCmdLineArgument("verbose"))
178  {
179  std::cout << "Phase sampling activated with period "
180  << m_phaseSamplePeriod << " and phase "
181  << m_phaseSamplePhase << "." << std::endl
182  << "Sampling within a tolerance of "
183  << std::setprecision(6) << m_phaseTolerance << "."
184  << std::endl;
185  }
186  }
187 
188  m_numSamples = 0;
189  m_index = 0;
190  m_outputIndex = 0;
191 
192  //
193  // FieldConvert modules
194  //
195  m_f = std::shared_ptr<Field>(new Field());
196  std::vector<std::string> modcmds;
197  // Process modules
198  std::stringstream moduleStream;
199  it = pParams.find("Modules");
200  if (it != pParams.end())
201  {
202  moduleStream.str(it->second);
203  }
204  while (!moduleStream.fail())
205  {
206  std::string sMod;
207  moduleStream >> sMod;
208  if (!moduleStream.fail())
209  {
210  modcmds.push_back(sMod);
211  }
212  }
213  // Output module
214  modcmds.push_back(m_outputFile);
215  // Create modules
216  CreateModules(modcmds);
217  // Strip options from m_outputFile
218  std::vector<std::string> tmp;
219  boost::split(tmp, m_outputFile, boost::is_any_of(":"));
220  m_outputFile = tmp[0];
221 }
222 
224 {
225 }
226 
229  const NekDouble &time)
230 {
231  v_FillVariablesName(pFields);
232 
233  // m_variables need to be filled by a derived class
234  m_outFields.resize(m_variables.size());
235  int nfield;
236 
237  for (int n = 0; n < m_variables.size(); ++n)
238  {
239  // if n >= pFields.num_elements() assum we have used n=0 field
240  nfield = (n < pFields.num_elements())? n: 0;
241 
242  m_outFields[n] =
243  Array<OneD, NekDouble>(pFields[nfield]->GetNcoeffs(), 0.0);
244  }
245 
246  m_fieldMetaData["InitialTime"] = boost::lexical_cast<std::string>(time);
247 
248  // Load restart file if necessary
249  if (m_restartFile != "")
250  {
251  // Load file
252  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
253  std::vector<std::vector<NekDouble> > fieldData;
254  LibUtilities::FieldMetaDataMap fieldMetaData;
257  fld->Import(m_restartFile, fieldDef, fieldData, fieldMetaData);
258 
259  // Extract fields to output
260  int nfield = -1, k;
261  for (int j = 0; j < m_variables.size(); ++j)
262  {
263  // see if m_variables is part of pFields definition and if
264  // so use that field for extract
265  for(k = 0; k < pFields.num_elements(); ++k)
266  {
267  if(pFields[k]->GetSession()->GetVariable(k)
268  == m_variables[j])
269  {
270  nfield = k;
271  break;
272  }
273  }
274  if(nfield == -1)
275  {
276  nfield = 0;
277  }
278 
279  for (int i = 0; i < fieldData.size(); ++i)
280  {
281  pFields[nfield]->ExtractDataToCoeffs(
282  fieldDef[i],
283  fieldData[i],
284  m_variables[j],
285  m_outFields[j]);
286  }
287  }
288 
289  // Load information for numSamples
290  if (fieldMetaData.count("NumberOfFieldDumps"))
291  {
292  m_numSamples = atoi(fieldMetaData["NumberOfFieldDumps"].c_str());
293  }
294  else
295  {
296  m_numSamples = 1;
297  }
298 
299  if(fieldMetaData.count("InitialTime"))
300  {
301  m_fieldMetaData["InitialTime"] = fieldMetaData["InitialTime"];
302  }
303 
304  // Divide by scale
305  NekDouble scale = v_GetScale();
306  for (int n = 0; n < m_outFields.size(); ++n)
307  {
308  Vmath::Smul(m_outFields[n].num_elements(),
309  1.0/scale,
310  m_outFields[n],
311  1,
312  m_outFields[n],
313  1);
314  }
315  }
316 }
317 
320 {
321  int nfield = pFields.num_elements();
322  m_variables.resize(pFields.num_elements());
323  for (int n = 0; n < nfield; ++n)
324  {
325  m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
326  }
327 
328  // Need to create a dummy coeffs vector to get extra variables names...
329  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
330  for (int n = 0; n < nfield; ++n)
331  {
332  coeffs[n] = pFields[n]->GetCoeffs();
333  }
334  // Get extra variables
335  auto equ = m_equ.lock();
336  ASSERTL0(equ, "Weak pointer expired");
337  equ->ExtraFldOutput(coeffs, m_variables);
338 }
339 
342  const NekDouble &time)
343 {
344  m_index++;
345  if (m_index % m_sampleFrequency > 0)
346  {
347  return;
348  }
349 
350  // Append extra fields
351  int nfield = pFields.num_elements();
352  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
353  for (int n = 0; n < nfield; ++n)
354  {
355  coeffs[n] = pFields[n]->GetCoeffs();
356  }
357  std::vector<std::string> variables = m_variables;
358  auto equ = m_equ.lock();
359  ASSERTL0(equ, "Weak pointer expired");
360  equ->ExtraFldOutput(coeffs, variables);
361 
362  if(m_phaseSample)
363  {
364  // The sample is added to the filter only if the current time
365  // corresponds to the correct phase. Introducing M as number of
366  // cycles and N nondimensional phase (between 0 and 1):
367  // t = M * m_phaseSamplePeriod + N * m_phaseSamplePeriod
368  int currentCycle = floor(time / m_phaseSamplePeriod);
369  NekDouble currentPhase = time / m_phaseSamplePeriod - currentCycle;
370 
371  // Evaluate phase relative to the requested value.
372  NekDouble relativePhase = fabs(m_phaseSamplePhase - currentPhase);
373 
374  // Check if relative phase is within required tolerance and sample.
375  // Care must be taken to handle the cases at phase 0 as the sample might
376  // have to be taken at the very end of the previous cycle instead.
377  if (relativePhase < m_phaseTolerance ||
378  fabs(relativePhase-1) < m_phaseTolerance)
379  {
380  m_numSamples++;
381  v_ProcessSample(pFields, coeffs, time);
382  if (m_session->GetComm()->GetRank() == 0 &&
383  m_session->DefinesCmdLineArgument("verbose"))
384  {
385  std::cout << "Sample: " << std::setw(8) << std::left
386  << m_numSamples << "Phase: " << std::setw(8)
387  << std::left << currentPhase << std::endl;
388  }
389  }
390  }
391  else
392  {
393  m_numSamples++;
394  v_ProcessSample(pFields, coeffs, time);
395  }
396 
397  if (m_index % m_outputFrequency == 0)
398  {
399  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
400  v_PrepareOutput(pFields, time);
401  OutputField(pFields, ++m_outputIndex);
402  }
403 }
404 
407  const NekDouble &time)
408 {
409  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
410  v_PrepareOutput(pFields, time);
411  OutputField(pFields);
412 }
413 
416  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
417  const NekDouble &time)
418 {
419  boost::ignore_unused(pFields, time);
420 
421  for(int n = 0; n < m_outFields.size(); ++n)
422  {
423  Vmath::Vcopy(m_outFields[n].num_elements(),
424  fieldcoeffs[n],
425  1,
426  m_outFields[n],
427  1);
428  }
429 }
430 
432  const Array<OneD, const MultiRegions::ExpListSharedPtr> &pFields, int dump)
433 {
434  NekDouble scale = v_GetScale();
435  for (int n = 0; n < m_outFields.size(); ++n)
436  {
437  Vmath::Smul(m_outFields[n].num_elements(),
438  scale,
439  m_outFields[n],
440  1,
441  m_outFields[n],
442  1);
443  }
444 
445  CreateFields(pFields);
446 
447  // Determine new file name
448  std::stringstream outname;
449  int dot = m_outputFile.find_last_of('.');
450  std::string name = m_outputFile.substr(0, dot);
451  std::string ext = m_outputFile.substr(dot, m_outputFile.length() - dot);
452  std::string suffix = v_GetFileSuffix();
453  if (dump == -1) // final dump
454  {
455  outname << name << suffix << ext;
456  }
457  else
458  {
459  outname << name << "_" << dump << suffix << ext;
460  }
461  m_modules[m_modules.size()-1]->RegisterConfig("outfile", outname.str());
462 
463  // Prevent checking before overwriting
464  po::options_description desc("Available options");
465  desc.add_options()
466  ("forceoutput,f",
467  "Force the output to be written without any checks");
468  po::variables_map vm;
469  vm.insert(std::make_pair("forceoutput", po::variable_value()));
470 
471  // Run field process.
472  for (int n = 0; n < SIZE_ModulePriority; ++n)
473  {
474  ModulePriority priority = static_cast<ModulePriority>(n);
475  for (int i = 0; i < m_modules.size(); ++i)
476  {
477  if(m_modules[i]->GetModulePriority() == priority)
478  {
479  m_modules[i]->Process(vm);
480  }
481  }
482  }
483 
484  // Empty m_f to save memory
485  m_f->ClearField();
486 
487  if (dump != -1) // not final dump so rescale
488  {
489  for (int n = 0; n < m_outFields.size(); ++n)
490  {
491  Vmath::Smul(m_outFields[n].num_elements(),
492  1.0 / scale,
493  m_outFields[n],
494  1,
495  m_outFields[n],
496  1);
497  }
498  }
499 }
500 
502 {
503  return true;
504 }
505 
506 void FilterFieldConvert::CreateModules(std::vector<std::string> &modcmds)
507 {
508  for (int i = 0; i < modcmds.size(); ++i)
509  {
510  // First split each command by the colon separator.
511  std::vector<std::string> tmp1;
512  ModuleKey module;
513  int offset = 1;
514 
515  boost::split(tmp1, modcmds[i], boost::is_any_of(":"));
516 
517  if (i == modcmds.size() - 1)
518  {
519  module.first = eOutputModule;
520 
521  // If no colon detected, automatically detect mesh type from
522  // file extension. Otherwise override and use tmp1[1] as the
523  // module to load. This also allows us to pass options to
524  // input/output modules. So, for example, to override
525  // filename.xml to be read as vtk, you use:
526  //
527  // filename.xml:vtk:opt1=arg1:opt2=arg2
528  if (tmp1.size() == 1)
529  {
530  int dot = tmp1[0].find_last_of('.') + 1;
531  std::string ext = tmp1[0].substr(dot, tmp1[0].length() - dot);
532 
533  module.second = ext;
534  tmp1.push_back(std::string("outfile=") + tmp1[0]);
535  }
536  else
537  {
538  module.second = tmp1[1];
539  tmp1.push_back(std::string("outfile=") + tmp1[0]);
540  offset++;
541  }
542  }
543  else
544  {
545  module.first = eProcessModule;
546  module.second = tmp1[0];
547  }
548 
549  // Create modules
550  ModuleSharedPtr mod;
551  mod = GetModuleFactory().CreateInstance(module, m_f);
552  m_modules.push_back(mod);
553 
554  // Set options for this module.
555  for (int j = offset; j < tmp1.size(); ++j)
556  {
557  std::vector<std::string> tmp2;
558  boost::split(tmp2, tmp1[j], boost::is_any_of("="));
559 
560  if (tmp2.size() == 1)
561  {
562  mod->RegisterConfig(tmp2[0]);
563  }
564  else if (tmp2.size() == 2)
565  {
566  mod->RegisterConfig(tmp2[0], tmp2[1]);
567  }
568  else
569  {
570  std::cerr << "ERROR: Invalid module configuration: format is "
571  << "either :arg or :arg=val" << std::endl;
572  abort();
573  }
574  }
575 
576  // Ensure configuration options have been set.
577  mod->SetDefaults();
578  }
579 
580  // Include equispaced output if needed
581  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
582  for (int i = 0; i < m_modules.size(); ++i)
583  {
584  ++modulesCount[m_modules[i]->GetModulePriority()];
585  }
586  if( modulesCount[eModifyPts] != 0 &&
587  modulesCount[eCreatePts] == 0 &&
588  modulesCount[eConvertExpToPts] == 0)
589  {
590  ModuleKey module;
591  ModuleSharedPtr mod;
592  module.first = eProcessModule;
593  module.second = std::string("equispacedoutput");
594  mod = GetModuleFactory().CreateInstance(module, m_f);
595  m_modules.insert(m_modules.end()-1, mod);
596  mod->SetDefaults();
597  }
598 
599  // Check if modules provided are compatible
601 }
602 
605 {
606  // Fill some parameters of m_f
607  m_f->m_session = m_session;
608  m_f->m_graph = pFields[0]->GetGraph();
609  m_f->m_comm = m_f->m_session->GetComm();
610  m_f->m_fieldMetaDataMap = m_fieldMetaData;
611  m_f->m_fieldPts = LibUtilities::NullPtsField;
612  // Create m_f->m_exp
613  m_f->m_numHomogeneousDir = 0;
614  if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
615  {
616  m_f->m_numHomogeneousDir = 1;
617  }
618  else if (pFields[0]->GetExpType() == MultiRegions::e3DH2D)
619  {
620  m_f->m_numHomogeneousDir = 2;
621  }
622 
623  m_f->m_exp.resize(m_variables.size());
624  m_f->m_exp[0] = pFields[0];
625  int nfield;
626  for (int n = 0; n < m_variables.size(); ++n)
627  {
628  // if n >= pFields.num_elements() assume we have used n=0 field
629  nfield = (n < pFields.num_elements())? n: 0;
630 
631  m_f->m_exp[n] = m_f->AppendExpList(
632  m_f->m_numHomogeneousDir, m_variables[0]);
633  m_f->m_exp[n]->SetWaveSpace(false);
634 
635  ASSERTL1(pFields[nfield]->GetNcoeffs() == m_outFields[n].num_elements(),
636  "pFields[nfield] does not have the "
637  "same number of coefficients as m_outFields[n]");
638 
639  m_f->m_exp[n]->ExtractCoeffsToCoeffs(pFields[nfield], m_outFields[n],
640  m_f->m_exp[n]->UpdateCoeffs());
641 
642  m_f->m_exp[n]->BwdTrans( m_f->m_exp[n]->GetCoeffs(),
643  m_f->m_exp[n]->UpdatePhys());
644  }
645  m_f->m_variables= m_variables;
646 }
647 
648 // This function checks validity conditions for the list of modules provided
649 void FilterFieldConvert::CheckModules(std::vector<ModuleSharedPtr> &modules)
650 {
651  // Count number of modules by priority
652  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
653  for (int i = 0; i < modules.size(); ++i)
654  {
655  ++modulesCount[modules[i]->GetModulePriority()];
656  }
657 
658  // FilterFieldConvert already starts with m_exp, so anything before
659  // eModifyExp is not valid, and also eCreatePts
660  if( modulesCount[eCreateGraph] != 0 ||
661  modulesCount[eCreateFieldData] != 0 ||
662  modulesCount[eModifyFieldData] != 0 ||
663  modulesCount[eCreateExp] != 0 ||
664  modulesCount[eFillExp] != 0 ||
665  modulesCount[eCreatePts] != 0)
666  {
667  std::stringstream ss;
668  ss << "Module(s): ";
669  for (int i = 0; i < modules.size(); ++i)
670  {
671  if(modules[i]->GetModulePriority() == eCreateGraph ||
672  modules[i]->GetModulePriority() == eCreateFieldData ||
673  modules[i]->GetModulePriority() == eModifyFieldData ||
674  modules[i]->GetModulePriority() == eCreateExp ||
675  modules[i]->GetModulePriority() == eFillExp ||
676  modules[i]->GetModulePriority() == eCreatePts)
677  {
678  ss << modules[i]->GetModuleName()<<" ";
679  }
680  }
681  ss << "not compatible with FilterFieldConvert.";
682  ASSERTL0(false, ss.str());
683  }
684 
685  // Modules of type eConvertExpToPts are not compatible with eBndExtraction
686  if( modulesCount[eConvertExpToPts] != 0 &&
687  modulesCount[eBndExtraction] != 0)
688  {
689  std::stringstream ss;
690  ss << "Module(s): ";
691  for (int i = 0; i < modules.size(); ++i)
692  {
693  if(modules[i]->GetModulePriority() == eBndExtraction)
694  {
695  ss << modules[i]->GetModuleName()<<" ";
696  }
697  }
698  ss << "is not compatible with module(s): ";
699  for (int i = 0; i < modules.size(); ++i)
700  {
701  if(modules[i]->GetModulePriority() == eConvertExpToPts)
702  {
703  ss << modules[i]->GetModuleName()<<" ";
704  }
705  }
706  ss << ".";
707  ASSERTL0(false, ss.str());
708  }
709 
710 }
711 
712 }
713 }
void OutputField(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int dump=-1)
LibUtilities::FieldMetaDataMap m_fieldMetaData
virtual SOLVER_UTILS_EXPORT ~FilterFieldConvert()
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
static std::string className
Name of the class.
void CreateFields(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
std::vector< Array< OneD, NekDouble > > m_outFields
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
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.
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
virtual SOLVER_UTILS_EXPORT std::string v_GetFileSuffix()
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::pair< ModuleType, std::string > ModuleKey
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
void CreateModules(std::vector< std::string > &modcmds)
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:87
virtual SOLVER_UTILS_EXPORT NekDouble v_GetScale()
double NekDouble
std::map< std::string, std::string > ParamMap
Definition: Filter.h:68
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:86
std::shared_ptr< Module > ModuleSharedPtr
virtual SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:180
void CheckModules(std::vector< ModuleSharedPtr > &modules)
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:226
std::vector< ModuleSharedPtr > m_modules
virtual SOLVER_UTILS_EXPORT void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
virtual SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual SOLVER_UTILS_EXPORT bool v_IsTimeDependent()
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
SOLVER_UTILS_EXPORT FilterFieldConvert(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
ModuleFactory & GetModuleFactory()