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
153  "PhaseAveragePeriod must be greater than 0.");
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.
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.size() assum we have used n=0 field
240  nfield = (n < pFields.size())? 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.size(); ++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  // Load information for outputIndex
305  if (fieldMetaData.count("FilterFileNum"))
306  {
307  m_outputIndex = atoi(fieldMetaData["FilterFileNum"].c_str());
308  }
309 
310  // Divide by scale
311  NekDouble scale = v_GetScale();
312  for (int n = 0; n < m_outFields.size(); ++n)
313  {
314  Vmath::Smul(m_outFields[n].size(),
315  1.0/scale,
316  m_outFields[n],
317  1,
318  m_outFields[n],
319  1);
320  }
321  }
322 }
323 
326 {
327  int nfield = pFields.size();
328  m_variables.resize(pFields.size());
329  for (int n = 0; n < nfield; ++n)
330  {
331  m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
332  }
333 
334  // Need to create a dummy coeffs vector to get extra variables names...
335  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
336  for (int n = 0; n < nfield; ++n)
337  {
338  coeffs[n] = pFields[n]->GetCoeffs();
339  }
340  // Get extra variables
341  auto equ = m_equ.lock();
342  ASSERTL0(equ, "Weak pointer expired");
343  equ->ExtraFldOutput(coeffs, m_variables);
344 }
345 
348  const NekDouble &time)
349 {
350  m_index++;
351  if (m_index % m_sampleFrequency > 0)
352  {
353  return;
354  }
355 
356  // Append extra fields
357  int nfield = pFields.size();
358  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
359  for (int n = 0; n < nfield; ++n)
360  {
361  coeffs[n] = pFields[n]->GetCoeffs();
362  }
363  std::vector<std::string> variables = m_variables;
364  auto equ = m_equ.lock();
365  ASSERTL0(equ, "Weak pointer expired");
366  equ->ExtraFldOutput(coeffs, variables);
367 
368  if(m_phaseSample)
369  {
370  // The sample is added to the filter only if the current time
371  // corresponds to the correct phase. Introducing M as number of
372  // cycles and N nondimensional phase (between 0 and 1):
373  // t = M * m_phaseSamplePeriod + N * m_phaseSamplePeriod
374  int currentCycle = floor(time / m_phaseSamplePeriod);
375  NekDouble currentPhase = time / m_phaseSamplePeriod - currentCycle;
376 
377  // Evaluate phase relative to the requested value.
378  NekDouble relativePhase = fabs(m_phaseSamplePhase - currentPhase);
379 
380  // Check if relative phase is within required tolerance and sample.
381  // Care must be taken to handle the cases at phase 0 as the sample might
382  // have to be taken at the very end of the previous cycle instead.
383  if (relativePhase < m_phaseTolerance ||
384  fabs(relativePhase-1) < m_phaseTolerance)
385  {
386  m_numSamples++;
387  v_ProcessSample(pFields, coeffs, time);
388  if (m_session->GetComm()->GetRank() == 0 &&
389  m_session->DefinesCmdLineArgument("verbose"))
390  {
391  std::cout << "Sample: " << std::setw(8) << std::left
392  << m_numSamples << "Phase: " << std::setw(8)
393  << std::left << currentPhase << std::endl;
394  }
395  }
396  }
397  else
398  {
399  m_numSamples++;
400  v_ProcessSample(pFields, coeffs, time);
401  }
402 
403  if (m_index % m_outputFrequency == 0)
404  {
405  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
406  v_PrepareOutput(pFields, time);
407  m_fieldMetaData["FilterFileNum"] =
408  boost::lexical_cast<std::string>(++m_outputIndex);
409  OutputField(pFields, m_outputIndex);
410  }
411 }
412 
415  const NekDouble &time)
416 {
417  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
418  v_PrepareOutput(pFields, time);
419  OutputField(pFields);
420 }
421 
424  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
425  const NekDouble &time)
426 {
427  boost::ignore_unused(pFields, time);
428 
429  for(int n = 0; n < m_outFields.size(); ++n)
430  {
431  Vmath::Vcopy(m_outFields[n].size(),
432  fieldcoeffs[n],
433  1,
434  m_outFields[n],
435  1);
436  }
437 }
438 
440  const Array<OneD, const MultiRegions::ExpListSharedPtr> &pFields, int dump)
441 {
442  NekDouble scale = v_GetScale();
443  for (int n = 0; n < m_outFields.size(); ++n)
444  {
445  Vmath::Smul(m_outFields[n].size(),
446  scale,
447  m_outFields[n],
448  1,
449  m_outFields[n],
450  1);
451  }
452 
453  CreateFields(pFields);
454 
455  // Determine new file name
456  std::stringstream outname;
457  int dot = m_outputFile.find_last_of('.');
458  std::string name = m_outputFile.substr(0, dot);
459  std::string ext = m_outputFile.substr(dot, m_outputFile.length() - dot);
460  std::string suffix = v_GetFileSuffix();
461  if (dump == -1) // final dump
462  {
463  outname << name << suffix << ext;
464  }
465  else
466  {
467  outname << name << "_" << dump << suffix << ext;
468  }
469  m_modules[m_modules.size()-1]->RegisterConfig("outfile", outname.str());
470 
471  // Prevent checking before overwriting
472  po::options_description desc("Available options");
473  desc.add_options()
474  ("forceoutput,f",
475  "Force the output to be written without any checks");
476  po::variables_map vm;
477  vm.insert(std::make_pair("forceoutput", po::variable_value()));
478 
479  // Run field process.
480  for (int n = 0; n < SIZE_ModulePriority; ++n)
481  {
482  ModulePriority priority = static_cast<ModulePriority>(n);
483  for (int i = 0; i < m_modules.size(); ++i)
484  {
485  if(m_modules[i]->GetModulePriority() == priority)
486  {
487  m_modules[i]->Process(vm);
488  }
489  }
490  }
491 
492  // Empty m_f to save memory
493  m_f->ClearField();
494 
495  if (dump != -1) // not final dump so rescale
496  {
497  for (int n = 0; n < m_outFields.size(); ++n)
498  {
499  Vmath::Smul(m_outFields[n].size(),
500  1.0 / scale,
501  m_outFields[n],
502  1,
503  m_outFields[n],
504  1);
505  }
506  }
507 }
508 
510 {
511  return true;
512 }
513 
514 void FilterFieldConvert::CreateModules(std::vector<std::string> &modcmds)
515 {
516  for (int i = 0; i < modcmds.size(); ++i)
517  {
518  // First split each command by the colon separator.
519  std::vector<std::string> tmp1;
520  ModuleKey module;
521  int offset = 1;
522 
523  boost::split(tmp1, modcmds[i], boost::is_any_of(":"));
524 
525  if (i == modcmds.size() - 1)
526  {
527  module.first = eOutputModule;
528 
529  // If no colon detected, automatically detect mesh type from
530  // file extension. Otherwise override and use tmp1[1] as the
531  // module to load. This also allows us to pass options to
532  // input/output modules. So, for example, to override
533  // filename.xml to be read as vtk, you use:
534  //
535  // filename.xml:vtk:opt1=arg1:opt2=arg2
536  if (tmp1.size() == 1)
537  {
538  int dot = tmp1[0].find_last_of('.') + 1;
539  std::string ext = tmp1[0].substr(dot, tmp1[0].length() - dot);
540 
541  module.second = ext;
542  tmp1.push_back(std::string("outfile=") + tmp1[0]);
543  }
544  else
545  {
546  module.second = tmp1[1];
547  tmp1.push_back(std::string("outfile=") + tmp1[0]);
548  offset++;
549  }
550  }
551  else
552  {
553  module.first = eProcessModule;
554  module.second = tmp1[0];
555  }
556 
557  // Create modules
558  ModuleSharedPtr mod;
559  mod = GetModuleFactory().CreateInstance(module, m_f);
560  m_modules.push_back(mod);
561 
562  // Set options for this module.
563  for (int j = offset; j < tmp1.size(); ++j)
564  {
565  std::vector<std::string> tmp2;
566  boost::split(tmp2, tmp1[j], boost::is_any_of("="));
567 
568  if (tmp2.size() == 1)
569  {
570  mod->RegisterConfig(tmp2[0]);
571  }
572  else if (tmp2.size() == 2)
573  {
574  mod->RegisterConfig(tmp2[0], tmp2[1]);
575  }
576  else
577  {
578  std::cerr << "ERROR: Invalid module configuration: format is "
579  << "either :arg or :arg=val" << std::endl;
580  abort();
581  }
582  }
583 
584  // Ensure configuration options have been set.
585  mod->SetDefaults();
586  }
587 
588  // Include equispaced output if needed
589  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
590  for (int i = 0; i < m_modules.size(); ++i)
591  {
592  ++modulesCount[m_modules[i]->GetModulePriority()];
593  }
594  if( modulesCount[eModifyPts] != 0 &&
595  modulesCount[eCreatePts] == 0 &&
596  modulesCount[eConvertExpToPts] == 0)
597  {
598  ModuleKey module;
599  ModuleSharedPtr mod;
600  module.first = eProcessModule;
601  module.second = std::string("equispacedoutput");
602  mod = GetModuleFactory().CreateInstance(module, m_f);
603  m_modules.insert(m_modules.end()-1, mod);
604  mod->SetDefaults();
605  }
606 
607  // Check if modules provided are compatible
609 }
610 
613 {
614  // Fill some parameters of m_f
615  m_f->m_session = m_session;
616  m_f->m_graph = pFields[0]->GetGraph();
617  m_f->m_comm = m_f->m_session->GetComm();
618  m_f->m_fieldMetaDataMap = m_fieldMetaData;
619  m_f->m_fieldPts = LibUtilities::NullPtsField;
620  // Create m_f->m_exp
621  m_f->m_numHomogeneousDir = 0;
622  if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
623  {
624  m_f->m_numHomogeneousDir = 1;
625  }
626  else if (pFields[0]->GetExpType() == MultiRegions::e3DH2D)
627  {
628  m_f->m_numHomogeneousDir = 2;
629  }
630 
631  m_f->m_exp.resize(m_variables.size());
632  m_f->m_exp[0] = pFields[0];
633  int nfield;
634  for (int n = 0; n < m_variables.size(); ++n)
635  {
636  // if n >= pFields.size() assume we have used n=0 field
637  nfield = (n < pFields.size())? n: 0;
638 
639  m_f->m_exp[n] = m_f->AppendExpList(
640  m_f->m_numHomogeneousDir, m_variables[0]);
641  m_f->m_exp[n]->SetWaveSpace(false);
642 
643  ASSERTL1(pFields[nfield]->GetNcoeffs() == m_outFields[n].size(),
644  "pFields[nfield] does not have the "
645  "same number of coefficients as m_outFields[n]");
646 
647  m_f->m_exp[n]->ExtractCoeffsToCoeffs(pFields[nfield], m_outFields[n],
648  m_f->m_exp[n]->UpdateCoeffs());
649 
650  m_f->m_exp[n]->BwdTrans( m_f->m_exp[n]->GetCoeffs(),
651  m_f->m_exp[n]->UpdatePhys());
652  }
653  m_f->m_variables= m_variables;
654 }
655 
656 // This function checks validity conditions for the list of modules provided
657 void FilterFieldConvert::CheckModules(std::vector<ModuleSharedPtr> &modules)
658 {
659  // Count number of modules by priority
660  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
661  for (int i = 0; i < modules.size(); ++i)
662  {
663  ++modulesCount[modules[i]->GetModulePriority()];
664  }
665 
666  // FilterFieldConvert already starts with m_exp, so anything before
667  // eModifyExp is not valid, and also eCreatePts
668  if( modulesCount[eCreateGraph] != 0 ||
669  modulesCount[eCreateFieldData] != 0 ||
670  modulesCount[eModifyFieldData] != 0 ||
671  modulesCount[eCreateExp] != 0 ||
672  modulesCount[eFillExp] != 0 ||
673  modulesCount[eCreatePts] != 0)
674  {
675  std::stringstream ss;
676  ss << "Module(s): ";
677  for (int i = 0; i < modules.size(); ++i)
678  {
679  if(modules[i]->GetModulePriority() == eCreateGraph ||
680  modules[i]->GetModulePriority() == eCreateFieldData ||
681  modules[i]->GetModulePriority() == eModifyFieldData ||
682  modules[i]->GetModulePriority() == eCreateExp ||
683  modules[i]->GetModulePriority() == eFillExp ||
684  modules[i]->GetModulePriority() == eCreatePts)
685  {
686  ss << modules[i]->GetModuleName()<<" ";
687  }
688  }
689  ss << "not compatible with FilterFieldConvert.";
690  ASSERTL0(false, ss.str());
691  }
692 
693  // Modules of type eConvertExpToPts are not compatible with eBndExtraction
694  if( modulesCount[eConvertExpToPts] != 0 &&
695  modulesCount[eBndExtraction] != 0)
696  {
697  std::stringstream ss;
698  ss << "Module(s): ";
699  for (int i = 0; i < modules.size(); ++i)
700  {
701  if(modules[i]->GetModulePriority() == eBndExtraction)
702  {
703  ss << modules[i]->GetModuleName()<<" ";
704  }
705  }
706  ss << "is not compatible with module(s): ";
707  for (int i = 0; i < modules.size(); ++i)
708  {
709  if(modules[i]->GetModulePriority() == eConvertExpToPts)
710  {
711  ss << modules[i]->GetModuleName()<<" ";
712  }
713  }
714  ss << ".";
715  ASSERTL0(false, ss.str());
716  }
717 
718 }
719 
720 }
721 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
NekDouble Evaluate() const
Definition: Equation.cpp:95
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
virtual SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
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 ~FilterFieldConvert()
static std::string className
Name of the class.
void CreateFields(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
virtual SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
void CheckModules(std::vector< ModuleSharedPtr > &modules)
virtual SOLVER_UTILS_EXPORT bool v_IsTimeDependent()
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
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.
SOLVER_UTILS_EXPORT FilterFieldConvert(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
void CreateModules(std::vector< std::string > &modcmds)
void OutputField(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int dump=-1)
LibUtilities::FieldMetaDataMap m_fieldMetaData
std::vector< Array< OneD, NekDouble > > m_outFields
virtual SOLVER_UTILS_EXPORT NekDouble v_GetScale()
virtual SOLVER_UTILS_EXPORT std::string v_GetFileSuffix()
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
std::vector< ModuleSharedPtr > m_modules
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:86
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:87
std::map< std::string, std::string > ParamMap
Definition: Filter.h:68
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
std::shared_ptr< Module > ModuleSharedPtr
Definition: Module.h:294
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:183
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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.cpp:225
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199