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 
191  //
192  // FieldConvert modules
193  //
194  m_f = std::shared_ptr<Field>(new Field());
195  std::vector<std::string> modcmds;
196  // Process modules
197  std::stringstream moduleStream;
198  it = pParams.find("Modules");
199  if (it != pParams.end())
200  {
201  moduleStream.str(it->second);
202  }
203  while (!moduleStream.fail())
204  {
205  std::string sMod;
206  moduleStream >> sMod;
207  if (!moduleStream.fail())
208  {
209  modcmds.push_back(sMod);
210  }
211  }
212  // Output module
213  modcmds.push_back(m_outputFile);
214  // Create modules
215  CreateModules(modcmds);
216  // Strip options from m_outputFile
217  std::vector<std::string> tmp;
218  boost::split(tmp, m_outputFile, boost::is_any_of(":"));
219  m_outputFile = tmp[0];
220 }
221 
223 {
224 }
225 
228  const NekDouble &time)
229 {
230  v_FillVariablesName(pFields);
231 
232  // m_variables need to be filled by a derived class
233  m_outFields.resize(m_variables.size());
234  int nfield;
235 
236  for (int n = 0; n < m_variables.size(); ++n)
237  {
238  // if n >= pFields.num_elements() assum we have used n=0 field
239  nfield = (n < pFields.num_elements())? n: 0;
240 
241  m_outFields[n] =
242  Array<OneD, NekDouble>(pFields[nfield]->GetNcoeffs(), 0.0);
243  }
244 
245  m_fieldMetaData["InitialTime"] = boost::lexical_cast<std::string>(time);
246 
247  // Load restart file if necessary
248  if (m_restartFile != "")
249  {
250  // Load file
251  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
252  std::vector<std::vector<NekDouble> > fieldData;
253  LibUtilities::FieldMetaDataMap fieldMetaData;
256  fld->Import(m_restartFile, fieldDef, fieldData, fieldMetaData);
257 
258  // Extract fields to output
259  int nfield = -1, k;
260  for (int j = 0; j < m_variables.size(); ++j)
261  {
262  // see if m_variables is part of pFields definition and if
263  // so use that field for extract
264  for(k = 0; k < pFields.num_elements(); ++k)
265  {
266  if(pFields[k]->GetSession()->GetVariable(k)
267  == m_variables[j])
268  {
269  nfield = k;
270  break;
271  }
272  }
273  if(nfield == -1)
274  {
275  nfield = 0;
276  }
277 
278  for (int i = 0; i < fieldData.size(); ++i)
279  {
280  pFields[nfield]->ExtractDataToCoeffs(
281  fieldDef[i],
282  fieldData[i],
283  m_variables[j],
284  m_outFields[j]);
285  }
286  }
287 
288  // Load information for numSamples
289  if (fieldMetaData.count("NumberOfFieldDumps"))
290  {
291  m_numSamples = atoi(fieldMetaData["NumberOfFieldDumps"].c_str());
292  }
293  else
294  {
295  m_numSamples = 1;
296  }
297 
298  if(fieldMetaData.count("InitialTime"))
299  {
300  m_fieldMetaData["InitialTime"] = fieldMetaData["InitialTime"];
301  }
302 
303  // Load information for outputIndex
304  if (fieldMetaData.count("FilterFileNum"))
305  {
306  m_outputIndex = atoi(fieldMetaData["FilterFileNum"].c_str());
307  }
308  else
309  {
310  m_outputIndex = 0;
311  }
312 
313  // Divide by scale
314  NekDouble scale = v_GetScale();
315  for (int n = 0; n < m_outFields.size(); ++n)
316  {
317  Vmath::Smul(m_outFields[n].num_elements(),
318  1.0/scale,
319  m_outFields[n],
320  1,
321  m_outFields[n],
322  1);
323  }
324  }
325 }
326 
329 {
330  int nfield = pFields.num_elements();
331  m_variables.resize(pFields.num_elements());
332  for (int n = 0; n < nfield; ++n)
333  {
334  m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
335  }
336 
337  // Need to create a dummy coeffs vector to get extra variables names...
338  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
339  for (int n = 0; n < nfield; ++n)
340  {
341  coeffs[n] = pFields[n]->GetCoeffs();
342  }
343  // Get extra variables
344  auto equ = m_equ.lock();
345  ASSERTL0(equ, "Weak pointer expired");
346  equ->ExtraFldOutput(coeffs, m_variables);
347 }
348 
351  const NekDouble &time)
352 {
353  m_index++;
354  if (m_index % m_sampleFrequency > 0)
355  {
356  return;
357  }
358 
359  // Append extra fields
360  int nfield = pFields.num_elements();
361  std::vector<Array<OneD, NekDouble> > coeffs(nfield);
362  for (int n = 0; n < nfield; ++n)
363  {
364  coeffs[n] = pFields[n]->GetCoeffs();
365  }
366  std::vector<std::string> variables = m_variables;
367  auto equ = m_equ.lock();
368  ASSERTL0(equ, "Weak pointer expired");
369  equ->ExtraFldOutput(coeffs, variables);
370 
371  if(m_phaseSample)
372  {
373  // The sample is added to the filter only if the current time
374  // corresponds to the correct phase. Introducing M as number of
375  // cycles and N nondimensional phase (between 0 and 1):
376  // t = M * m_phaseSamplePeriod + N * m_phaseSamplePeriod
377  int currentCycle = floor(time / m_phaseSamplePeriod);
378  NekDouble currentPhase = time / m_phaseSamplePeriod - currentCycle;
379 
380  // Evaluate phase relative to the requested value.
381  NekDouble relativePhase = fabs(m_phaseSamplePhase - currentPhase);
382 
383  // Check if relative phase is within required tolerance and sample.
384  // Care must be taken to handle the cases at phase 0 as the sample might
385  // have to be taken at the very end of the previous cycle instead.
386  if (relativePhase < m_phaseTolerance ||
387  fabs(relativePhase-1) < m_phaseTolerance)
388  {
389  m_numSamples++;
390  v_ProcessSample(pFields, coeffs, time);
391  if (m_session->GetComm()->GetRank() == 0 &&
392  m_session->DefinesCmdLineArgument("verbose"))
393  {
394  std::cout << "Sample: " << std::setw(8) << std::left
395  << m_numSamples << "Phase: " << std::setw(8)
396  << std::left << currentPhase << std::endl;
397  }
398  }
399  }
400  else
401  {
402  m_numSamples++;
403  v_ProcessSample(pFields, coeffs, time);
404  }
405 
406  if (m_index % m_outputFrequency == 0)
407  {
408  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
409  v_PrepareOutput(pFields, time);
410  m_fieldMetaData["FilterFileNum"] =
411  boost::lexical_cast<std::string>(++m_outputIndex);
412  OutputField(pFields, m_outputIndex);
413  }
414 }
415 
418  const NekDouble &time)
419 {
420  m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
421  v_PrepareOutput(pFields, time);
422  OutputField(pFields);
423 }
424 
427  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
428  const NekDouble &time)
429 {
430  boost::ignore_unused(pFields, time);
431 
432  for(int n = 0; n < m_outFields.size(); ++n)
433  {
434  Vmath::Vcopy(m_outFields[n].num_elements(),
435  fieldcoeffs[n],
436  1,
437  m_outFields[n],
438  1);
439  }
440 }
441 
443  const Array<OneD, const MultiRegions::ExpListSharedPtr> &pFields, int dump)
444 {
445  NekDouble scale = v_GetScale();
446  for (int n = 0; n < m_outFields.size(); ++n)
447  {
448  Vmath::Smul(m_outFields[n].num_elements(),
449  scale,
450  m_outFields[n],
451  1,
452  m_outFields[n],
453  1);
454  }
455 
456  CreateFields(pFields);
457 
458  // Determine new file name
459  std::stringstream outname;
460  int dot = m_outputFile.find_last_of('.');
461  std::string name = m_outputFile.substr(0, dot);
462  std::string ext = m_outputFile.substr(dot, m_outputFile.length() - dot);
463  std::string suffix = v_GetFileSuffix();
464  if (dump == -1) // final dump
465  {
466  outname << name << suffix << ext;
467  }
468  else
469  {
470  outname << name << "_" << dump << suffix << ext;
471  }
472  m_modules[m_modules.size()-1]->RegisterConfig("outfile", outname.str());
473 
474  // Prevent checking before overwriting
475  po::options_description desc("Available options");
476  desc.add_options()
477  ("forceoutput,f",
478  "Force the output to be written without any checks");
479  po::variables_map vm;
480  vm.insert(std::make_pair("forceoutput", po::variable_value()));
481 
482  // Run field process.
483  for (int n = 0; n < SIZE_ModulePriority; ++n)
484  {
485  ModulePriority priority = static_cast<ModulePriority>(n);
486  for (int i = 0; i < m_modules.size(); ++i)
487  {
488  if(m_modules[i]->GetModulePriority() == priority)
489  {
490  m_modules[i]->Process(vm);
491  }
492  }
493  }
494 
495  // Empty m_f to save memory
496  m_f->ClearField();
497 
498  if (dump != -1) // not final dump so rescale
499  {
500  for (int n = 0; n < m_outFields.size(); ++n)
501  {
502  Vmath::Smul(m_outFields[n].num_elements(),
503  1.0 / scale,
504  m_outFields[n],
505  1,
506  m_outFields[n],
507  1);
508  }
509  }
510 }
511 
513 {
514  return true;
515 }
516 
517 void FilterFieldConvert::CreateModules(std::vector<std::string> &modcmds)
518 {
519  for (int i = 0; i < modcmds.size(); ++i)
520  {
521  // First split each command by the colon separator.
522  std::vector<std::string> tmp1;
523  ModuleKey module;
524  int offset = 1;
525 
526  boost::split(tmp1, modcmds[i], boost::is_any_of(":"));
527 
528  if (i == modcmds.size() - 1)
529  {
530  module.first = eOutputModule;
531 
532  // If no colon detected, automatically detect mesh type from
533  // file extension. Otherwise override and use tmp1[1] as the
534  // module to load. This also allows us to pass options to
535  // input/output modules. So, for example, to override
536  // filename.xml to be read as vtk, you use:
537  //
538  // filename.xml:vtk:opt1=arg1:opt2=arg2
539  if (tmp1.size() == 1)
540  {
541  int dot = tmp1[0].find_last_of('.') + 1;
542  std::string ext = tmp1[0].substr(dot, tmp1[0].length() - dot);
543 
544  module.second = ext;
545  tmp1.push_back(std::string("outfile=") + tmp1[0]);
546  }
547  else
548  {
549  module.second = tmp1[1];
550  tmp1.push_back(std::string("outfile=") + tmp1[0]);
551  offset++;
552  }
553  }
554  else
555  {
556  module.first = eProcessModule;
557  module.second = tmp1[0];
558  }
559 
560  // Create modules
561  ModuleSharedPtr mod;
562  mod = GetModuleFactory().CreateInstance(module, m_f);
563  m_modules.push_back(mod);
564 
565  // Set options for this module.
566  for (int j = offset; j < tmp1.size(); ++j)
567  {
568  std::vector<std::string> tmp2;
569  boost::split(tmp2, tmp1[j], boost::is_any_of("="));
570 
571  if (tmp2.size() == 1)
572  {
573  mod->RegisterConfig(tmp2[0]);
574  }
575  else if (tmp2.size() == 2)
576  {
577  mod->RegisterConfig(tmp2[0], tmp2[1]);
578  }
579  else
580  {
581  std::cerr << "ERROR: Invalid module configuration: format is "
582  << "either :arg or :arg=val" << std::endl;
583  abort();
584  }
585  }
586 
587  // Ensure configuration options have been set.
588  mod->SetDefaults();
589  }
590 
591  // Include equispaced output if needed
592  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
593  for (int i = 0; i < m_modules.size(); ++i)
594  {
595  ++modulesCount[m_modules[i]->GetModulePriority()];
596  }
597  if( modulesCount[eModifyPts] != 0 &&
598  modulesCount[eCreatePts] == 0 &&
599  modulesCount[eConvertExpToPts] == 0)
600  {
601  ModuleKey module;
602  ModuleSharedPtr mod;
603  module.first = eProcessModule;
604  module.second = std::string("equispacedoutput");
605  mod = GetModuleFactory().CreateInstance(module, m_f);
606  m_modules.insert(m_modules.end()-1, mod);
607  mod->SetDefaults();
608  }
609 
610  // Check if modules provided are compatible
612 }
613 
616 {
617  // Fill some parameters of m_f
618  m_f->m_session = m_session;
619  m_f->m_graph = pFields[0]->GetGraph();
620  m_f->m_comm = m_f->m_session->GetComm();
621  m_f->m_fieldMetaDataMap = m_fieldMetaData;
622  m_f->m_fieldPts = LibUtilities::NullPtsField;
623  // Create m_f->m_exp
624  m_f->m_numHomogeneousDir = 0;
625  if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
626  {
627  m_f->m_numHomogeneousDir = 1;
628  }
629  else if (pFields[0]->GetExpType() == MultiRegions::e3DH2D)
630  {
631  m_f->m_numHomogeneousDir = 2;
632  }
633 
634  m_f->m_exp.resize(m_variables.size());
635  m_f->m_exp[0] = pFields[0];
636  int nfield;
637  for (int n = 0; n < m_variables.size(); ++n)
638  {
639  // if n >= pFields.num_elements() assume we have used n=0 field
640  nfield = (n < pFields.num_elements())? n: 0;
641 
642  m_f->m_exp[n] = m_f->AppendExpList(
643  m_f->m_numHomogeneousDir, m_variables[0]);
644  m_f->m_exp[n]->SetWaveSpace(false);
645 
646  ASSERTL1(pFields[nfield]->GetNcoeffs() == m_outFields[n].num_elements(),
647  "pFields[nfield] does not have the "
648  "same number of coefficients as m_outFields[n]");
649 
650  m_f->m_exp[n]->ExtractCoeffsToCoeffs(pFields[nfield], m_outFields[n],
651  m_f->m_exp[n]->UpdateCoeffs());
652 
653  m_f->m_exp[n]->BwdTrans( m_f->m_exp[n]->GetCoeffs(),
654  m_f->m_exp[n]->UpdatePhys());
655  }
656  m_f->m_variables= m_variables;
657 }
658 
659 // This function checks validity conditions for the list of modules provided
660 void FilterFieldConvert::CheckModules(std::vector<ModuleSharedPtr> &modules)
661 {
662  // Count number of modules by priority
663  Array< OneD, int> modulesCount(SIZE_ModulePriority,0);
664  for (int i = 0; i < modules.size(); ++i)
665  {
666  ++modulesCount[modules[i]->GetModulePriority()];
667  }
668 
669  // FilterFieldConvert already starts with m_exp, so anything before
670  // eModifyExp is not valid, and also eCreatePts
671  if( modulesCount[eCreateGraph] != 0 ||
672  modulesCount[eCreateFieldData] != 0 ||
673  modulesCount[eModifyFieldData] != 0 ||
674  modulesCount[eCreateExp] != 0 ||
675  modulesCount[eFillExp] != 0 ||
676  modulesCount[eCreatePts] != 0)
677  {
678  std::stringstream ss;
679  ss << "Module(s): ";
680  for (int i = 0; i < modules.size(); ++i)
681  {
682  if(modules[i]->GetModulePriority() == eCreateGraph ||
683  modules[i]->GetModulePriority() == eCreateFieldData ||
684  modules[i]->GetModulePriority() == eModifyFieldData ||
685  modules[i]->GetModulePriority() == eCreateExp ||
686  modules[i]->GetModulePriority() == eFillExp ||
687  modules[i]->GetModulePriority() == eCreatePts)
688  {
689  ss << modules[i]->GetModuleName()<<" ";
690  }
691  }
692  ss << "not compatible with FilterFieldConvert.";
693  ASSERTL0(false, ss.str());
694  }
695 
696  // Modules of type eConvertExpToPts are not compatible with eBndExtraction
697  if( modulesCount[eConvertExpToPts] != 0 &&
698  modulesCount[eBndExtraction] != 0)
699  {
700  std::stringstream ss;
701  ss << "Module(s): ";
702  for (int i = 0; i < modules.size(); ++i)
703  {
704  if(modules[i]->GetModulePriority() == eBndExtraction)
705  {
706  ss << modules[i]->GetModuleName()<<" ";
707  }
708  }
709  ss << "is not compatible with module(s): ";
710  for (int i = 0; i < modules.size(); ++i)
711  {
712  if(modules[i]->GetModulePriority() == eConvertExpToPts)
713  {
714  ss << modules[i]->GetModuleName()<<" ";
715  }
716  }
717  ss << ".";
718  ASSERTL0(false, ss.str());
719  }
720 
721 }
722 
723 }
724 }
#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:199
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
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
std::shared_ptr< Module > ModuleSharedPtr
ModuleFactory & GetModuleFactory()
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:180
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
double NekDouble
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064