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/algorithm/string/classification.hpp>
38#include <boost/algorithm/string/predicate.hpp>
39#include <boost/algorithm/string/split.hpp>
40#include <boost/program_options.hpp>
41
42namespace Nektar::SolverUtils
43{
47
50 const std::shared_ptr<EquationSystem> &pEquation, const ParamMap &pParams)
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 // Time after which we need to write checkfiles
79 it = pParams.find("OutputStartTime");
80 if (it == pParams.end())
81 {
83 }
84 else
85 {
86 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
88 }
89
90 // Restart file
91 it = pParams.find("RestartFile");
92 if (it == pParams.end())
93 {
94 m_restartFile = "";
95 }
96 else
97 {
98 ASSERTL0(it->second.length() > 0, "Missing parameter 'RestartFile'.");
99 if (it->second.find_last_of('.') != std::string::npos)
100 {
101 m_restartFile = it->second;
102 }
103 else
104 {
105 std::stringstream outname;
106 outname << it->second << ".fld";
107 m_restartFile = outname.str();
108 }
109 }
110
111 // OutputFrequency
112 it = pParams.find("OutputFrequency");
113 if (it == pParams.end())
114 {
115 m_outputFrequency = m_session->GetParameter("NumSteps");
116 }
117 else
118 {
119 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
120 m_outputFrequency = round(equ.Evaluate());
121 }
122
123 // The base class can use SampleFrequency = OutputFrequency
124 // (Derived classes need to override this if needed)
126
127 // Phase sampling option
128 it = pParams.find("PhaseAverage");
129 if (it == pParams.end())
130 {
131 m_phaseSample = false;
132 }
133 else
134 {
135 std::string sOption = it->second.c_str();
136 m_phaseSample = (boost::iequals(sOption, "true")) ||
137 (boost::iequals(sOption, "yes"));
138 }
139
140 if (m_phaseSample)
141 {
142 auto itPeriod = pParams.find("PhaseAveragePeriod");
143 auto itPhase = pParams.find("PhaseAveragePhase");
144
145 // Error if only one of the required params for PhaseAverage is present
146 ASSERTL0(
147 (itPeriod != pParams.end() && itPhase != pParams.end()),
148 "The phase sampling feature requires both 'PhaseAveragePeriod' and "
149 "'PhaseAveragePhase' to be set.");
150
151 LibUtilities::Equation equPeriod(m_session->GetInterpreter(),
152 itPeriod->second);
153 m_phaseSamplePeriod = equPeriod.Evaluate();
154
155 LibUtilities::Equation equPhase(m_session->GetInterpreter(),
156 itPhase->second);
157 m_phaseSamplePhase = equPhase.Evaluate();
158
159 // Check that phase and period are within required limits
161 "PhaseAveragePeriod must be greater than 0.");
163 "PhaseAveragePhase must be between 0 and 1.");
164
165 // Load sampling frequency, overriding the previous value
166 it = pParams.find("SampleFrequency");
167 if (it == pParams.end())
168 {
170 }
171 else
172 {
173 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
174 m_sampleFrequency = round(equ.Evaluate());
175 }
176
177 // Compute tolerance within which sampling occurs.
179
180 // Display worst case scenario sampling tolerance for exact phase, if
181 // verbose option is active
182 if (m_session->GetComm()->GetRank() == 0 &&
183 m_session->DefinesCmdLineArgument("verbose"))
184 {
185 std::cout << "Phase sampling activated with period "
186 << m_phaseSamplePeriod << " and phase "
187 << m_phaseSamplePhase << "." << std::endl
188 << "Sampling within a tolerance of "
189 << std::setprecision(6) << m_phaseTolerance << "."
190 << std::endl;
191 }
192 }
193
194 m_numSamples = 0;
195 m_index = 0;
196 m_outputIndex = 0;
197
198 //
199 // FieldConvert modules
200 //
201 m_f = std::shared_ptr<Field>(new Field());
202 std::vector<std::string> modcmds;
203 // Process modules
204 std::stringstream moduleStream;
205 it = pParams.find("Modules");
206 if (it != pParams.end())
207 {
208 moduleStream.str(it->second);
209 }
210 while (!moduleStream.fail())
211 {
212 std::string sMod;
213 moduleStream >> sMod;
214 if (!moduleStream.fail())
215 {
216 modcmds.push_back(sMod);
217 }
218 }
219 // Output module
220 modcmds.push_back(m_outputFile);
221 // Create modules
222 CreateModules(modcmds);
223 // Strip options from m_outputFile
224 std::vector<std::string> tmp;
225 boost::split(tmp, m_outputFile, boost::is_any_of(":"));
226 m_outputFile = tmp[0];
227
228 // Prevent checking before overwriting
229 it = pParams.find("options");
230 if (it != pParams.end())
231 {
232 int argc = 0;
233 std::string strargv;
234 std::vector<char *> argv;
235
236 strargv = "dummy " + it->second;
237 strargv.push_back((char)0);
238 bool flag = true;
239 for (size_t i = 0; strargv[i]; ++i)
240 {
241 if (strargv[i] != ' ' && flag)
242 {
243 argv.push_back(&strargv[i]);
244 ++argc;
245 flag = false;
246 }
247 if (strargv[i] == ' ')
248 {
249 flag = true;
250 strargv[i] = 0;
251 }
252 }
253
254 po::options_description desc("Available options");
255
256 // clang-format off
257 desc.add_options()
258 ("help,h",
259 "Produce this help message.")
260 ("modules-list,l",
261 "Print the list of available modules.")
262 ("output-points,n", po::value<int>(),
263 "Output at n equipspaced points along the "
264 "collapsed coordinates (for .dat, .vtu).")
265 ("output-points-hom-z", po::value<int>(),
266 "Number of planes in the z-direction for output of "
267 "Homogeneous 1D expansion(for .dat, .vtu).")
268 ("error,e",
269 "Write error of fields for regression checking")
270 ("force-output,f",
271 "Force the output to be written without any checks")
272 ("range,r", po::value<std::string>(),
273 "Define output range i.e. (-r xmin,xmax,ymin,ymax,zmin,zmax) "
274 "in which any vertex is contained.")
275 ("no-equispaced",
276 "Do not use equispaced output.")
277 ("nparts", po::value<int>(),
278 "Define nparts if running serial problem to mimic "
279 "parallel run with many partitions.")
280 ("npz", po::value<int>(),
281 "Used to define number of partitions in z for Homogeneous1D "
282 "expansions for parallel runs.")
283 ("npt", po::value<int>(),
284 "Used to define number of partitions in time for Parareal runs. ")
285 ("onlyshape", po::value<std::string>(),
286 "Only use element with defined shape type i.e. -onlyshape "
287 " Tetrahedron")
288 ("part-only", po::value<int>(),
289 "Partition into specified npart partitions and exit")
290 ("part-only-overlapping", po::value<int>(),
291 "Partition into specified npart overlapping partitions and exit")
292 ("modules-opt,p", po::value<std::string>(),
293 "Print options for a module.")
294 ("module,m", po::value<std::vector<std::string> >(),
295 "Specify modules which are to be used.")
296 ("use-session-variables",
297 "Use variables defined in session for output")
298 ("use-session-expansion",
299 "Use expansion defined in session.")
300 ("verbose,v",
301 "Enable verbose mode.");
302 // clang-format on
303
304 po::options_description hidden("Hidden options");
305
306 // clang-format off
307 hidden.add_options()
308 ("input-file", po::value<std::vector<std::string> >(),
309 "Input filename");
310 // clang-format on
311
312 po::options_description cmdline_options;
313 cmdline_options.add(hidden).add(desc);
314
315 po::options_description visible("Allowed options");
316 visible.add(desc);
317
318 po::positional_options_description p;
319 p.add("input-file", -1);
320
321 try
322 {
323 po::store(po::command_line_parser(argc, &(argv[0]))
324 .options(cmdline_options)
325 .positional(p)
326 .run(),
327 m_vm);
328 po::notify(m_vm);
329 }
330 catch (const std::exception &e)
331 {
332 std::cerr << e.what() << std::endl;
333 std::cerr << desc;
334 }
335 }
336 m_vm.insert(std::make_pair("force-output", po::variable_value()));
337}
338
340{
341}
342
345 const NekDouble &time)
346{
347 v_FillVariablesName(pFields);
348
349 // m_variables need to be filled by a derived class
350 m_outFields.resize(m_variables.size());
351 int nfield;
352
353 for (int n = 0; n < m_variables.size(); ++n)
354 {
355 // if n >= pFields.size() assum we have used n=0 field
356 nfield = (n < pFields.size()) ? n : 0;
357
358 m_outFields[n] =
359 Array<OneD, NekDouble>(pFields[nfield]->GetNcoeffs(), 0.0);
360 }
361
362 m_fieldMetaData["InitialTime"] = boost::lexical_cast<std::string>(time);
363
364 // Load restart file if necessary
365 if (m_restartFile != "")
366 {
367 // Load file
368 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
369 std::vector<std::vector<NekDouble>> fieldData;
370 LibUtilities::FieldMetaDataMap fieldMetaData;
373 fld->Import(m_restartFile, fieldDef, fieldData, fieldMetaData);
374
375 // Extract fields to output
376 int nfield = -1, k;
377 for (int j = 0; j < m_variables.size(); ++j)
378 {
379 // see if m_variables is part of pFields definition and if
380 // so use that field for extract
381 for (k = 0; k < pFields.size(); ++k)
382 {
383 if (pFields[k]->GetSession()->GetVariable(k) == m_variables[j])
384 {
385 nfield = k;
386 break;
387 }
388 }
389 if (nfield == -1)
390 {
391 nfield = 0;
392 }
393
394 for (int i = 0; i < fieldData.size(); ++i)
395 {
396 pFields[nfield]->ExtractDataToCoeffs(
397 fieldDef[i], fieldData[i], m_variables[j], m_outFields[j]);
398 }
399 }
400
401 // Load information for numSamples
402 if (fieldMetaData.count("NumberOfFieldDumps"))
403 {
404 m_numSamples = atoi(fieldMetaData["NumberOfFieldDumps"].c_str());
405 }
406 else
407 {
408 m_numSamples = 1;
409 }
410
411 if (fieldMetaData.count("InitialTime"))
412 {
413 m_fieldMetaData["InitialTime"] = fieldMetaData["InitialTime"];
414 }
415
416 // Load information for outputIndex
417 if (fieldMetaData.count("FilterFileNum"))
418 {
419 m_outputIndex = atoi(fieldMetaData["FilterFileNum"].c_str());
420 }
421
422 // Divide by scale
423 NekDouble scale = v_GetScale();
424 for (int n = 0; n < m_outFields.size(); ++n)
425 {
426 Vmath::Smul(m_outFields[n].size(), 1.0 / scale, m_outFields[n], 1,
427 m_outFields[n], 1);
428 }
429 }
430}
431
434{
435 int nfield = pFields.size();
436 m_variables.resize(pFields.size());
437 for (int n = 0; n < nfield; ++n)
438 {
439 m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
440 }
441
442 // Need to create a dummy coeffs vector to get extra variables names...
443 std::vector<Array<OneD, NekDouble>> coeffs(nfield);
444 for (int n = 0; n < nfield; ++n)
445 {
446 coeffs[n] = pFields[n]->GetCoeffs();
447 }
448 // Get extra variables
449 auto equ = m_equ.lock();
450 ASSERTL0(equ, "Weak pointer expired");
451 equ->ExtraFldOutput(coeffs, m_variables);
452}
453
456 const NekDouble &time)
457{
458 m_index++;
459 if (m_index % m_sampleFrequency > 0 ||
460 (time - m_outputStartTime) < -1.0e-07)
461 {
462 return;
463 }
464
465 // Append extra fields
466 int nfield = pFields.size();
467 std::vector<Array<OneD, NekDouble>> coeffs(nfield);
468 for (int n = 0; n < nfield; ++n)
469 {
470 coeffs[n] = pFields[n]->GetCoeffs();
471 }
472 std::vector<std::string> variables = m_variables;
473 auto equ = m_equ.lock();
474 ASSERTL0(equ, "Weak pointer expired");
475 equ->ExtraFldOutput(coeffs, variables);
476
477 if (m_phaseSample)
478 {
479 // The sample is added to the filter only if the current time
480 // corresponds to the correct phase. Introducing M as number of
481 // cycles and N nondimensional phase (between 0 and 1):
482 // t = M * m_phaseSamplePeriod + N * m_phaseSamplePeriod
483 int currentCycle = floor(time / m_phaseSamplePeriod);
484 NekDouble currentPhase = time / m_phaseSamplePeriod - currentCycle;
485
486 // Evaluate phase relative to the requested value.
487 NekDouble relativePhase = fabs(m_phaseSamplePhase - currentPhase);
488
489 // Check if relative phase is within required tolerance and sample.
490 // Care must be taken to handle the cases at phase 0 as the sample might
491 // have to be taken at the very end of the previous cycle instead.
492 if (relativePhase < m_phaseTolerance ||
493 fabs(relativePhase - 1) < m_phaseTolerance)
494 {
495 m_numSamples++;
496 v_ProcessSample(pFields, coeffs, time);
497 if (m_session->GetComm()->GetRank() == 0 &&
498 m_session->DefinesCmdLineArgument("verbose"))
499 {
500 std::cout << "Sample: " << std::setw(8) << std::left
501 << m_numSamples << "Phase: " << std::setw(8)
502 << std::left << currentPhase << std::endl;
503 }
504 }
505 }
506 else
507 {
508 m_numSamples++;
509 v_ProcessSample(pFields, coeffs, time);
510 }
511
512 if (m_index % m_outputFrequency == 0)
513 {
514 m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
515 v_PrepareOutput(pFields, time);
516 m_fieldMetaData["FilterFileNum"] = std::to_string(++m_outputIndex);
517 OutputField(pFields, m_outputIndex);
518 }
519}
520
523 const NekDouble &time)
524{
525 m_fieldMetaData["FinalTime"] = boost::lexical_cast<std::string>(time);
526 v_PrepareOutput(pFields, time);
527 OutputField(pFields);
528}
529
532 &pFields,
533 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
534 [[maybe_unused]] const NekDouble &time)
535{
536 for (int n = 0; n < m_outFields.size(); ++n)
537 {
538 Vmath::Vcopy(m_outFields[n].size(), fieldcoeffs[n], 1, m_outFields[n],
539 1);
540 }
541}
542
545{
546 NekDouble scale = v_GetScale();
547 for (int n = 0; n < m_outFields.size(); ++n)
548 {
549 Vmath::Smul(m_outFields[n].size(), scale, m_outFields[n], 1,
550 m_outFields[n], 1);
551 }
552
553 CreateFields(pFields);
554
555 // Determine new file name
556 std::stringstream outname;
557 int dot = m_outputFile.find_last_of('.');
558 std::string name = m_outputFile.substr(0, dot);
559 std::string ext = m_outputFile.substr(dot, m_outputFile.length() - dot);
560 std::string suffix = v_GetFileSuffix();
561 if (dump == -1) // final dump
562 {
563 outname << name << suffix << ext;
564 }
565 else
566 {
567 outname << name << "_" << dump << suffix << ext;
568 }
569 m_modules[m_modules.size() - 1]->RegisterConfig("outfile", outname.str());
570
571 // Run field process.
572 for (int n = 0; n < SIZE_ModulePriority; ++n)
573 {
574 ModulePriority priority = static_cast<ModulePriority>(n);
575 for (int i = 0; i < m_modules.size(); ++i)
576 {
577 if (m_modules[i]->GetModulePriority() == priority)
578 {
579 m_modules[i]->Process(m_vm);
580 }
581 }
582 }
583
584 // Empty m_f to save memory
585 m_f->ClearField();
586
587 if (dump != -1) // not final dump so rescale
588 {
589 for (int n = 0; n < m_outFields.size(); ++n)
590 {
591 Vmath::Smul(m_outFields[n].size(), 1.0 / scale, m_outFields[n], 1,
592 m_outFields[n], 1);
593 }
594 }
595}
596
598{
599 return true;
600}
601
602void FilterFieldConvert::CreateModules(std::vector<std::string> &modcmds)
603{
604 for (int i = 0; i < modcmds.size(); ++i)
605 {
606 // First split each command by the colon separator.
607 std::vector<std::string> tmp1;
608 ModuleKey module;
609 int offset = 1;
610
611 boost::split(tmp1, modcmds[i], boost::is_any_of(":"));
612
613 if (i == modcmds.size() - 1)
614 {
615 module.first = eOutputModule;
616
617 // If no colon detected, automatically detect mesh type from
618 // file extension. Otherwise override and use tmp1[1] as the
619 // module to load. This also allows us to pass options to
620 // input/output modules. So, for example, to override
621 // filename.xml to be read as vtk, you use:
622 //
623 // filename.xml:vtk:opt1=arg1:opt2=arg2
624 if (tmp1.size() == 1)
625 {
626 int dot = tmp1[0].find_last_of('.') + 1;
627 std::string ext = tmp1[0].substr(dot, tmp1[0].length() - dot);
628 module.second = ext;
629 tmp1.push_back(std::string("outfile=") + tmp1[0]);
630 }
631 else
632 {
633 module.second = tmp1[1];
634 tmp1.push_back(std::string("outfile=") + tmp1[0]);
635 offset++;
636 }
637 }
638 else
639 {
640 module.first = eProcessModule;
641 module.second = tmp1[0];
642 }
643
644 // Create modules
645 ModuleSharedPtr mod;
646 mod = GetModuleFactory().CreateInstance(module, m_f);
647 m_modules.push_back(mod);
648
649 // Set options for this module.
650 for (int j = offset; j < tmp1.size(); ++j)
651 {
652 std::vector<std::string> tmp2;
653 boost::split(tmp2, tmp1[j], boost::is_any_of("="));
654
655 if (tmp2.size() == 1)
656 {
657 mod->RegisterConfig(tmp2[0]);
658 }
659 else if (tmp2.size() == 2)
660 {
661 mod->RegisterConfig(tmp2[0], tmp2[1]);
662 }
663 else
664 {
665 std::cerr << "ERROR: Invalid module configuration: format is "
666 << "either :arg or :arg=val" << std::endl;
667 abort();
668 }
669 }
670
671 // Ensure configuration options have been set.
672 mod->SetDefaults();
673 }
674
675 // Include equispaced output if needed
676 Array<OneD, int> modulesCount(SIZE_ModulePriority, 0);
677 for (int i = 0; i < m_modules.size(); ++i)
678 {
679 ++modulesCount[m_modules[i]->GetModulePriority()];
680 }
681 if (modulesCount[eModifyPts] != 0 && modulesCount[eCreatePts] == 0 &&
682 modulesCount[eConvertExpToPts] == 0)
683 {
684 ModuleKey module;
685 ModuleSharedPtr mod;
686 module.first = eProcessModule;
687 module.second = std::string("equispacedoutput");
688 mod = GetModuleFactory().CreateInstance(module, m_f);
689 m_modules.insert(m_modules.end() - 1, mod);
690 mod->SetDefaults();
691 }
692
693 // Check if modules provided are compatible
695}
696
699{
700 // Fill some parameters of m_f
701 m_f->m_session = m_session;
702 m_f->m_graph = pFields[0]->GetGraph();
703 m_f->m_comm = m_f->m_session->GetComm();
704 m_f->m_fieldMetaDataMap = m_fieldMetaData;
705 m_f->m_fieldPts = LibUtilities::NullPtsField;
706 // Create m_f->m_exp
707 m_f->m_numHomogeneousDir = 0;
708 if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
709 {
710 m_f->m_numHomogeneousDir = 1;
711 }
712 else if (pFields[0]->GetExpType() == MultiRegions::e3DH2D)
713 {
714 m_f->m_numHomogeneousDir = 2;
715 }
716
717 m_f->m_exp.resize(m_variables.size());
718 m_f->m_exp[0] = pFields[0];
719 int nfield;
720 for (int n = 0; n < m_variables.size(); ++n)
721 {
722 // if n >= pFields.size() assume we have used n=0 field
723 nfield = (n < pFields.size()) ? n : 0;
724
725 m_f->m_exp[n] =
726 m_f->AppendExpList(m_f->m_numHomogeneousDir, m_variables[0]);
727 m_f->m_exp[n]->SetWaveSpace(false);
728
729 ASSERTL1(pFields[nfield]->GetNcoeffs() == m_outFields[n].size(),
730 "pFields[nfield] does not have the "
731 "same number of coefficients as m_outFields[n]");
732
733 m_f->m_exp[n]->ExtractCoeffsToCoeffs(pFields[nfield], m_outFields[n],
734 m_f->m_exp[n]->UpdateCoeffs());
735
736 m_f->m_exp[n]->BwdTrans(m_f->m_exp[n]->GetCoeffs(),
737 m_f->m_exp[n]->UpdatePhys());
738 }
739 m_f->m_variables = m_variables;
740}
741
742// This function checks validity conditions for the list of modules provided
743void FilterFieldConvert::CheckModules(std::vector<ModuleSharedPtr> &modules)
744{
745 // Count number of modules by priority
746 Array<OneD, int> modulesCount(SIZE_ModulePriority, 0);
747 for (int i = 0; i < modules.size(); ++i)
748 {
749 ++modulesCount[modules[i]->GetModulePriority()];
750 }
751
752 // FilterFieldConvert already starts with m_exp, so anything before
753 // eModifyExp is not valid, and also eCreatePts
754 if (modulesCount[eCreateGraph] != 0 ||
755 modulesCount[eCreateFieldData] != 0 ||
756 modulesCount[eModifyFieldData] != 0 || modulesCount[eCreateExp] != 0 ||
757 modulesCount[eFillExp] != 0 || modulesCount[eCreatePts] != 0)
758 {
759 std::stringstream ss;
760 ss << "Module(s): ";
761 for (int i = 0; i < modules.size(); ++i)
762 {
763 if (modules[i]->GetModulePriority() == eCreateGraph ||
764 modules[i]->GetModulePriority() == eCreateFieldData ||
765 modules[i]->GetModulePriority() == eModifyFieldData ||
766 modules[i]->GetModulePriority() == eCreateExp ||
767 modules[i]->GetModulePriority() == eFillExp ||
768 modules[i]->GetModulePriority() == eCreatePts)
769 {
770 ss << modules[i]->GetModuleName() << " ";
771 }
772 }
773 ss << "not compatible with FilterFieldConvert.";
774 ASSERTL0(false, ss.str());
775 }
776
777 // Modules of type eConvertExpToPts are not compatible with eBndExtraction
778 if (modulesCount[eConvertExpToPts] != 0 &&
779 modulesCount[eBndExtraction] != 0)
780 {
781 std::stringstream ss;
782 ss << "Module(s): ";
783 for (int i = 0; i < modules.size(); ++i)
784 {
785 if (modules[i]->GetModulePriority() == eBndExtraction)
786 {
787 ss << modules[i]->GetModuleName() << " ";
788 }
789 }
790 ss << "is not compatible with module(s): ";
791 for (int i = 0; i < modules.size(); ++i)
792 {
793 if (modules[i]->GetModulePriority() == eConvertExpToPts)
794 {
795 ss << modules[i]->GetModuleName() << " ";
796 }
797 }
798 ss << ".";
799 ASSERTL0(false, ss.str());
800 }
801}
802
803} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
SOLVER_UTILS_EXPORT bool v_IsTimeDependent() override
virtual SOLVER_UTILS_EXPORT void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
static std::string className
Name of the class.
void CreateFields(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
void CheckModules(std::vector< ModuleSharedPtr > &modules)
SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
SOLVER_UTILS_EXPORT void v_Update(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)
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
SOLVER_UTILS_EXPORT ~FilterFieldConvert() override
void CreateModules(std::vector< std::string > &modcmds)
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_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::shared_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
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()
std::vector< ModuleSharedPtr > m_modules
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
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
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:185
FilterFactory & GetFilterFactory()
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
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825