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