Nektar++
EquationSystem.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: EquationSystem.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:
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
39
45
49
52
54
55#include <boost/format.hpp>
56
57#include <iostream>
58#include <string>
59
60using namespace std;
61
62namespace Nektar
63{
64namespace SolverUtils
65{
66
69 LibUtilities::SessionReader::RegisterEnumValue("DEALIASING", "False", 1)};
70
72 LibUtilities::SessionReader::RegisterEnumValue("Projection", "Continuous",
74 LibUtilities::SessionReader::RegisterEnumValue("Projection", "CONTINUOUS",
81 "Projection", "DisContinuous", MultiRegions::eDiscontinuous),
83 "Projection", "Mixed_CG_Discontinuous",
86 "Projection", "MixedCGDG", MultiRegions::eMixed_CG_Discontinuous),
87};
88
89/**
90 * @class EquationSystem
91 *
92 * This class is a base class for all solver implementations. It
93 * provides the underlying generic functionality and interface for
94 * solving equations.
95 *
96 * To solve a steady-state equation, create a derived class from this
97 * class and reimplement the virtual functions to provide custom
98 * implementation for the problem.
99 *
100 * To solve unsteady problems, derive from the UnsteadySystem class
101 * instead which provides general time integration.
102 */
104{
105 static EquationSystemFactory instance;
106 return instance;
107}
108
109/**
110 * This constructor is protected as the objects of this class are never
111 * instantiated directly.
112 * @param pSession The session reader holding problem parameters.
113 */
117 : m_comm(pSession->GetComm()), m_session(pSession), m_graph(pGraph),
118 m_lambda(0), m_infosteps(10),
119 m_fieldMetaDataMap(LibUtilities::NullFieldMetaDataMap)
120{
121 // set up session names in fieldMetaDataMap
122 const vector<std::string> filenames = m_session->GetFilenames();
123
124 for (int i = 0; i < filenames.size(); ++i)
125 {
126 string sessionname = "SessionName";
127 sessionname += boost::lexical_cast<std::string>(i);
128 m_fieldMetaDataMap[sessionname] = filenames[i];
129 m_fieldMetaDataMap["ChkFileNum"] = boost::lexical_cast<std::string>(0);
130 }
131}
132
133/**
134 * @brief Destructor for class EquationSystem.
135 */
137{
139 LocalRegions::MatrixKey::opLess>::ClearManager();
141 LocalRegions::MatrixKey::opLess>::ClearManager();
142}
143
144/**
145 * @brief Initialisation object for EquationSystem.
146 */
147void EquationSystem::v_InitObject(bool DeclareFields)
148{
149 // Save the basename of input file name for output details
150 m_sessionName = m_session->GetSessionName();
151
152 // Instantiate a field reader/writer
154
155 // Also read and store the boundary conditions
159
160 // Set space dimension for use in class
161 m_spacedim = m_graph->GetSpaceDimension();
162
163 // Setting parameteres for homogenous problems
164 m_HomoDirec = 0;
165 m_useFFT = false;
166 m_homogen_dealiasing = false;
167 m_singleMode = false;
168 m_halfMode = false;
169 m_multipleModes = false;
171
172 m_verbose = m_session->DefinesCmdLineArgument("verbose");
173
174 if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
175 {
176 std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
177 m_spacedim = 3;
178
179 if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D") ||
180 (HomoStr == "1D") || (HomoStr == "Homo1D"))
181 {
183 m_session->LoadParameter("LZ", m_LhomZ);
184 m_HomoDirec = 1;
185
186 if (m_session->DefinesSolverInfo("ModeType"))
187 {
188 m_session->MatchSolverInfo("ModeType", "SingleMode",
189 m_singleMode, false);
190 m_session->MatchSolverInfo("ModeType", "HalfMode", m_halfMode,
191 false);
192 m_session->MatchSolverInfo("ModeType", "MultipleModes",
193 m_multipleModes, false);
194 }
195
196 // Stability Analysis flags
197 if (m_session->DefinesSolverInfo("ModeType"))
198 {
199 if (m_singleMode)
200 {
201 m_npointsZ = 2;
202 }
203 else if (m_halfMode)
204 {
205 m_npointsZ = 1;
206 }
207 else if (m_multipleModes)
208 {
209 m_npointsZ = m_session->GetParameter("HomModesZ");
210 }
211 else
212 {
213 ASSERTL0(false, "SolverInfo ModeType not valid");
214 }
215 }
216 else
217 {
218 m_npointsZ = m_session->GetParameter("HomModesZ");
219 }
220 }
221
222 if ((HomoStr == "HOMOGENEOUS2D") || (HomoStr == "Homogeneous2D") ||
223 (HomoStr == "2D") || (HomoStr == "Homo2D"))
224 {
226 m_session->LoadParameter("HomModesY", m_npointsY);
227 m_session->LoadParameter("LY", m_LhomY);
228 m_session->LoadParameter("HomModesZ", m_npointsZ);
229 m_session->LoadParameter("LZ", m_LhomZ);
230 m_HomoDirec = 2;
231 }
232
233 if ((HomoStr == "HOMOGENEOUS3D") || (HomoStr == "Homogeneous3D") ||
234 (HomoStr == "3D") || (HomoStr == "Homo3D"))
235 {
237 m_session->LoadParameter("HomModesY", m_npointsY);
238 m_session->LoadParameter("LY", m_LhomY);
239 m_session->LoadParameter("HomModesZ", m_npointsZ);
240 m_session->LoadParameter("LZ", m_LhomZ);
241 m_HomoDirec = 2;
242 }
243
244 m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
245
246 m_session->MatchSolverInfo("DEALIASING", "True", m_homogen_dealiasing,
247 false);
248 }
249 else
250 {
251 // set to default value so can use to identify 2d or 3D
252 // (homogeneous) expansions
253 m_npointsZ = 1;
254 }
255
256 m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "True",
257 m_specHP_dealiasing, false);
258 if (m_specHP_dealiasing == false)
259 {
260 m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "On",
261 m_specHP_dealiasing, false);
262 }
263
264 // Options to determine type of projection from file or directly
265 // from constructor
266 if (m_session->DefinesSolverInfo("PROJECTION"))
267 {
268 std::string ProjectStr = m_session->GetSolverInfo("PROJECTION");
269
270 if ((ProjectStr == "Continuous") || (ProjectStr == "Galerkin") ||
271 (ProjectStr == "CONTINUOUS") || (ProjectStr == "GALERKIN"))
272 {
274 }
275 else if ((ProjectStr == "MixedCGDG") ||
276 (ProjectStr == "Mixed_CG_Discontinuous"))
277 {
279 }
280 else if (ProjectStr == "DisContinuous")
281 {
283 }
284 else
285 {
286 ASSERTL0(false, "PROJECTION value not recognised");
287 }
288 }
289 else
290 {
291 cerr << "Projection type not specified in SOLVERINFO,"
292 "defaulting to continuous Galerkin"
293 << endl;
295 }
296
297 // Enforce singularity check for some problems
299
300 int i;
301 int nvariables = m_session->GetVariables().size();
302 bool DeclareCoeffPhysArrays = true;
303
305 m_spacedim = m_graph->GetSpaceDimension() + m_HomoDirec;
306 m_expdim = m_graph->GetMeshDimension();
307
308 if (DeclareFields) // declare field if required
309 {
310 /// Continuous field
313 {
314 switch (m_expdim)
315 {
316 case 1:
317 {
320 {
321 const LibUtilities::PointsKey PkeyY(
323 const LibUtilities::BasisKey BkeyY(
325 const LibUtilities::PointsKey PkeyZ(
327 const LibUtilities::BasisKey BkeyZ(
329
330 for (i = 0; i < m_fields.size(); i++)
331 {
334 AllocateSharedPtr(m_session, BkeyY, BkeyZ,
337 m_session->GetVariable(i));
338 }
339 }
340 else
341 {
342 for (i = 0; i < m_fields.size(); i++)
343 {
344 m_fields[i] =
348 m_session->GetVariable(i));
349 }
350 }
351 break;
352 }
353 case 2:
354 {
356 {
357 // Fourier single mode stability analysis
358 if (m_singleMode)
359 {
360 const LibUtilities::PointsKey PkeyZ(
363
364 const LibUtilities::BasisKey BkeyZ(
366 PkeyZ);
367
368 for (i = 0; i < m_fields.size(); i++)
369 {
372 AllocateSharedPtr(
373 m_session, BkeyZ, m_LhomZ, m_useFFT,
375 m_session->GetVariable(i),
377 }
378 }
379 // Half mode stability analysis
380 else if (m_halfMode)
381 {
382 const LibUtilities::PointsKey PkeyZ(
385
386 const LibUtilities::BasisKey BkeyZR(
388 PkeyZ);
389
390 const LibUtilities::BasisKey BkeyZI(
392 PkeyZ);
393
394 for (i = 0; i < m_fields.size(); i++)
395 {
396 if (m_session->GetVariable(i).compare("w") == 0)
397 {
399 MultiRegions::
400 ContField3DHomogeneous1D>::
401 AllocateSharedPtr(
402 m_session, BkeyZI, m_LhomZ,
404 m_graph, m_session->GetVariable(i),
406 }
407 else
408 {
410 MultiRegions::
411 ContField3DHomogeneous1D>::
412 AllocateSharedPtr(
413 m_session, BkeyZR, m_LhomZ,
415 m_graph, m_session->GetVariable(i),
417 }
418 }
419 }
420 // Normal homogeneous 1D
421 else
422 {
423 const LibUtilities::PointsKey PkeyZ(
425 const LibUtilities::BasisKey BkeyZ(
427
428 for (i = 0; i < m_fields.size(); i++)
429 {
432 AllocateSharedPtr(
433 m_session, BkeyZ, m_LhomZ, m_useFFT,
435 m_session->GetVariable(i),
437 }
438 }
439 }
440 else
441 {
442 i = 0;
446 m_session->GetVariable(i),
447 DeclareCoeffPhysArrays,
449 m_fields[0] = firstfield;
450 for (i = 1; i < m_fields.size(); i++)
451 {
452 if (m_graph->SameExpansionInfo(
453 m_session->GetVariable(0),
454 m_session->GetVariable(i)))
455 {
456 m_fields[i] =
459 *firstfield, m_graph,
460 m_session->GetVariable(i),
461 DeclareCoeffPhysArrays,
463 }
464 else
465 {
466 m_fields[i] =
470 m_session->GetVariable(i),
471 DeclareCoeffPhysArrays,
473 }
474 }
475
476 if (m_projectionType ==
478 {
479 /// Setting up the normals
482
483 for (i = 0; i < m_spacedim; ++i)
484 {
485 m_traceNormals[i] =
487 }
488
489 m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
490 }
491 }
492
493 break;
494 }
495 case 3:
496 {
497 i = 0;
501 m_session->GetVariable(i),
502 DeclareCoeffPhysArrays,
504
505 m_fields[0] = firstfield;
506 for (i = 1; i < m_fields.size(); i++)
507 {
508 if (m_graph->SameExpansionInfo(
509 m_session->GetVariable(0),
510 m_session->GetVariable(i)))
511 {
512 m_fields[i] =
515 *firstfield, m_graph,
516 m_session->GetVariable(i),
517 DeclareCoeffPhysArrays,
519 }
520 else
521 {
522 m_fields[i] =
526 m_session->GetVariable(i),
527 DeclareCoeffPhysArrays,
529 }
530 }
531
532 if (m_projectionType ==
534 {
535 /// Setting up the normals
538 for (i = 0; i < m_spacedim; ++i)
539 {
540 m_traceNormals[i] =
542 }
543
544 m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
545 // Call the trace on all fields to ensure DG setup.
546 for (i = 1; i < m_fields.size(); ++i)
547 {
548 m_fields[i]->GetTrace();
549 }
550 }
551 break;
552 }
553 default:
554 ASSERTL0(false, "Expansion dimension not recognised");
555 break;
556 }
557 }
558 // Discontinuous field
559 else
560 {
561 switch (m_expdim)
562 {
563 case 1:
564 {
567 {
568 const LibUtilities::PointsKey PkeyY(
570 const LibUtilities::BasisKey BkeyY(
572 const LibUtilities::PointsKey PkeyZ(
574 const LibUtilities::BasisKey BkeyZ(
576
577 for (i = 0; i < m_fields.size(); i++)
578 {
581 AllocateSharedPtr(m_session, BkeyY, BkeyZ,
584 m_session->GetVariable(i));
585 }
586 }
587 else
588 {
589 for (i = 0; i < m_fields.size(); i++)
590 {
591 m_fields[i] =
595 m_session->GetVariable(i));
596 }
597 }
598
599 break;
600 }
601 case 2:
602 {
604 {
605 const LibUtilities::PointsKey PkeyZ(
607 const LibUtilities::BasisKey BkeyZ(
609
610 for (i = 0; i < m_fields.size(); i++)
611 {
614 AllocateSharedPtr(m_session, BkeyZ, m_LhomZ,
615 m_useFFT,
617 m_session->GetVariable(i));
618 }
619 }
620 else
621 {
622 i = 0;
626 m_session->GetVariable(i));
627 m_fields[0] = firstfield;
628 for (i = 1; i < m_fields.size(); i++)
629 {
630 if (m_graph->SameExpansionInfo(
631 m_session->GetVariable(0),
632 m_session->GetVariable(i)))
633 {
634 m_fields[i] =
637 *firstfield, m_graph,
638 m_session->GetVariable(i));
639 }
640 else
641 {
642 m_fields[i] =
646 m_session->GetVariable(i));
647 }
648 }
649 }
650
651 break;
652 }
653 case 3:
654 {
656 {
657 ASSERTL0(
658 false,
659 "3D fully periodic problems not implemented yet");
660 }
661 else
662 {
663 i = 0;
667 m_session->GetVariable(i));
668 m_fields[0] = firstfield;
669 for (i = 1; i < m_fields.size(); i++)
670 {
671 if (m_graph->SameExpansionInfo(
672 m_session->GetVariable(0),
673 m_session->GetVariable(i)))
674 {
675 m_fields[i] =
678 *firstfield, m_graph,
679 m_session->GetVariable(i));
680 }
681 else
682 {
683 m_fields[i] =
687 m_session->GetVariable(i));
688 }
689 }
690 }
691 break;
692 }
693 default:
694 ASSERTL0(false, "Expansion dimension not recognised");
695 break;
696 }
697
698 // Setting up the normals
700
701 for (i = 0; i < m_spacedim; ++i)
702 {
703 m_traceNormals[i] =
705 }
706
707 m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
708 }
709 // Zero all physical fields initially
711 }
712
713 // Set Default Parameter
714 m_session->LoadParameter("Time", m_time, 0.0);
715 m_session->LoadParameter("TimeStep", m_timestep, 0.0);
716 m_session->LoadParameter("NumSteps", m_steps, 0);
717 m_session->LoadParameter("IO_CheckSteps", m_checksteps, 0);
718 m_session->LoadParameter("IO_CheckTime", m_checktime, 0.0);
719 m_session->LoadParameter("FinTime", m_fintime, 0);
720 m_session->LoadParameter("NumQuadPointsError", m_NumQuadPointsError, 0);
721
722 // Check uniqueness of checkpoint output
723 ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
724 (m_checktime > 0.0 && m_checksteps == 0) ||
725 (m_checktime == 0.0 && m_checksteps > 0),
726 "Only one of IO_CheckTime and IO_CheckSteps "
727 "should be set!");
728 m_session->LoadParameter("TimeIncrementFactor", m_TimeIncrementFactor, 1.0);
729
730 m_nchk = 0;
731 m_iterPIT = 0;
732}
733
734/**
735 *
736 */
738 std::string name, const MultiRegions::ExpListSharedPtr &field, bool cache)
739{
740 MultiRegions::ExpListSharedPtr vField = field;
741 if (!field)
742 {
743 vField = m_fields[0];
744 }
745
746 if (cache)
747 {
748 if ((m_sessionFunctions.find(name) == m_sessionFunctions.end()) ||
749 (m_sessionFunctions[name]->GetSession() != m_session) ||
750 (m_sessionFunctions[name]->GetExpansion() != vField))
751 {
754 m_session, vField, name, cache);
755 }
756
757 return m_sessionFunctions[name];
758 }
759 else
760 {
762 new SessionFunction(m_session, vField, name, cache));
763 }
764}
765
766/**
767 * If boundary conditions are time-dependent, they will be evaluated at
768 * the time specified.
769 * @param time The time at which to evaluate the BCs
770 */
772{
773 std::string varName;
774 int nvariables = m_fields.size();
775 for (int i = 0; i < nvariables; ++i)
776 {
777 varName = m_session->GetVariable(i);
778 m_fields[i]->EvaluateBoundaryConditions(time, varName);
779 }
780}
781
782/**
783 * Compute the error in the L2-norm.
784 * @param field The field to compare.
785 * @param exactsoln The exact solution to compare with.
786 * @param Normalised Normalise L2-error.
787 * @returns Error in the L2-norm.
788 */
790 const Array<OneD, NekDouble> &exactsoln,
791 bool Normalised)
792{
793 NekDouble L2error = -1.0;
794
795 if (m_NumQuadPointsError == 0)
796 {
797 if (m_fields[field]->GetPhysState() == false)
798 {
799 m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
800 m_fields[field]->UpdatePhys());
801 }
802
803 if (exactsoln.size())
804 {
805 L2error =
806 m_fields[field]->L2(m_fields[field]->GetPhys(), exactsoln);
807 }
808 else if (m_session->DefinesFunction("ExactSolution"))
809 {
810 Array<OneD, NekDouble> exactsoln(m_fields[field]->GetNpoints());
811
812 GetFunction("ExactSolution")
813 ->Evaluate(m_session->GetVariable(field), exactsoln, m_time);
814
815 L2error =
816 m_fields[field]->L2(m_fields[field]->GetPhys(), exactsoln);
817 }
818 else
819 {
820 L2error = m_fields[field]->L2(m_fields[field]->GetPhys());
821 }
822
823 if (Normalised == true)
824 {
825 Array<OneD, NekDouble> one(m_fields[field]->GetNpoints(), 1.0);
826
827 NekDouble Vol = m_fields[field]->Integral(one);
828 L2error = sqrt(L2error * L2error / Vol);
829 }
830 }
831 else
832 {
833 Array<OneD, NekDouble> L2INF(2);
834 L2INF = ErrorExtraPoints(field);
835 L2error = L2INF[0];
836 }
837 return L2error;
838}
839
840/**
841 * Compute the error in the L_inf-norm
842 * @param field The field to compare.
843 * @param exactsoln The exact solution to compare with.
844 * @returns Error in the L_inft-norm.
845 */
847 const Array<OneD, NekDouble> &exactsoln)
848{
849 NekDouble Linferror = -1.0;
850
851 if (m_NumQuadPointsError == 0)
852 {
853 if (m_fields[field]->GetPhysState() == false)
854 {
855 m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
856 m_fields[field]->UpdatePhys());
857 }
858
859 if (exactsoln.size())
860 {
861 Linferror =
862 m_fields[field]->Linf(m_fields[field]->GetPhys(), exactsoln);
863 }
864 else if (m_session->DefinesFunction("ExactSolution"))
865 {
866 Array<OneD, NekDouble> exactsoln(m_fields[field]->GetNpoints());
867
868 GetFunction("ExactSolution")
869 ->Evaluate(m_session->GetVariable(field), exactsoln, m_time);
870
871 Linferror =
872 m_fields[field]->Linf(m_fields[field]->GetPhys(), exactsoln);
873 }
874 else
875 {
876 Linferror = m_fields[field]->Linf(m_fields[field]->GetPhys());
877 }
878 }
879 else
880 {
881 Array<OneD, NekDouble> L2INF(2);
882 L2INF = ErrorExtraPoints(field);
883 Linferror = L2INF[1];
884 }
885
886 return Linferror;
887}
888
889/**
890 * Compute the error in the L2-norm, L-inf for a larger number of
891 * quadrature points.
892 * @param field The field to compare.
893 * @returns Error in the L2-norm and L-inf norm.
894 */
896{
897 int NumModes = GetNumExpModes();
898 Array<OneD, NekDouble> L2INF(2);
899
903 LibUtilities::eGaussRadauMAlpha1Beta0);
909 PkeyT1);
911 PkeyT2);
913 PkeyQ1);
915 PkeyQ2);
916
917 LibUtilities::BasisKeyVector Tkeys, Qkeys;
918
919 // make a copy of the ExpansionInfoMap
920 SpatialDomains::ExpansionInfoMap NewExpInfo = m_graph->GetExpansionInfo();
923 NewExpInfo);
924
925 // reset new graph with new keys
926 Tkeys.push_back(BkeyT1);
927 Tkeys.push_back(BkeyT2);
928 m_graph->ResetExpansionInfoToBasisKey(ExpInfo, LibUtilities::eTriangle,
929 Tkeys);
930 Qkeys.push_back(BkeyQ1);
931 Qkeys.push_back(BkeyQ2);
932 m_graph->ResetExpansionInfoToBasisKey(ExpInfo, LibUtilities::eQuadrilateral,
933 Qkeys);
934
937 NewExpInfo);
938
939 int ErrorCoordim = ErrorExp->GetCoordim(0);
940 int ErrorNq = ErrorExp->GetTotPoints();
941
942 Array<OneD, NekDouble> ErrorXc0(ErrorNq, 0.0);
943 Array<OneD, NekDouble> ErrorXc1(ErrorNq, 0.0);
944 Array<OneD, NekDouble> ErrorXc2(ErrorNq, 0.0);
945
946 switch (ErrorCoordim)
947 {
948 case 1:
949 ErrorExp->GetCoords(ErrorXc0);
950 break;
951 case 2:
952 ErrorExp->GetCoords(ErrorXc0, ErrorXc1);
953 break;
954 case 3:
955 ErrorExp->GetCoords(ErrorXc0, ErrorXc1, ErrorXc2);
956 break;
957 }
959 m_session->GetFunction("ExactSolution", field);
960
961 // Evaluate the exact solution
962 Array<OneD, NekDouble> ErrorSol(ErrorNq);
963
964 exSol->Evaluate(ErrorXc0, ErrorXc1, ErrorXc2, m_time, ErrorSol);
965
966 // Calcualte spectral/hp approximation on the quadrature points
967 // of this new expansion basis
968 ErrorExp->BwdTrans(m_fields[field]->GetCoeffs(), ErrorExp->UpdatePhys());
969
970 L2INF[0] = ErrorExp->L2(ErrorExp->GetPhys(), ErrorSol);
971 L2INF[1] = ErrorExp->Linf(ErrorExp->GetPhys(), ErrorSol);
972
973 return L2INF;
974}
975
976/**
977 * Set the physical fields based on a restart file, or a function
978 * describing the initial condition given in the session.
979 * @param initialtime Time at which to evaluate the function.
980 * @param dumpInitialConditions Write the initial condition to file?
981 */
983 bool dumpInitialConditions,
984 const int domain)
985{
986 boost::ignore_unused(initialtime);
987
988 if (m_session->GetComm()->GetRank() == 0)
989 {
990 cout << "Initial Conditions:" << endl;
991 }
992
993 if (m_session->DefinesFunction("InitialConditions"))
994 {
995 GetFunction("InitialConditions")
996 ->Evaluate(m_session->GetVariables(), m_fields, m_time, domain);
997 // Enforce C0 Continutiy of initial condiiton
1000 {
1001 for (int i = 0; i < m_fields.size(); ++i)
1002 {
1003 m_fields[i]->LocalToGlobal();
1004 m_fields[i]->GlobalToLocal();
1005 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
1006 m_fields[i]->UpdatePhys());
1007 }
1008 }
1009
1010 if (m_session->GetComm()->GetRank() == 0)
1011 {
1012 for (int i = 0; i < m_fields.size(); ++i)
1013 {
1014 std::string varName = m_session->GetVariable(i);
1015 cout << " - Field " << varName << ": "
1016 << GetFunction("InitialConditions")
1017 ->Describe(varName, domain)
1018 << endl;
1019 }
1020 }
1021 }
1022 else
1023 {
1024 int nq = m_fields[0]->GetNpoints();
1025 for (int i = 0; i < m_fields.size(); i++)
1026 {
1027 Vmath::Zero(nq, m_fields[i]->UpdatePhys(), 1);
1028 m_fields[i]->SetPhysState(true);
1029 Vmath::Zero(m_fields[i]->GetNcoeffs(), m_fields[i]->UpdateCoeffs(),
1030 1);
1031 if (m_session->GetComm()->GetRank() == 0)
1032 {
1033 cout << " - Field " << m_session->GetVariable(i)
1034 << ": 0 (default)" << endl;
1035 }
1036 }
1037 }
1038
1039 if (dumpInitialConditions && m_checksteps && m_nchk == 0 &&
1040 !ParallelInTime())
1041 {
1043 }
1044 else if (dumpInitialConditions && m_nchk == 0 && ParallelInTime())
1045 {
1046 std::string newdir = m_sessionName + ".pit";
1047 if (!fs::is_directory(newdir))
1048 {
1049 fs::create_directory(newdir);
1050 }
1051 if (m_comm->GetTimeComm()->GetRank() == 0)
1052 {
1053 WriteFld(newdir + "/" + m_sessionName + "_0.fld");
1054 }
1055 }
1056 ++m_nchk;
1057}
1058
1059/**
1060 *
1061 */
1063 Array<OneD, NekDouble> &outfield,
1064 const NekDouble time)
1065{
1066 ASSERTL0(outfield.size() == m_fields[field]->GetNpoints(),
1067 "ExactSolution array size mismatch.");
1068 Vmath::Zero(outfield.size(), outfield, 1);
1069 if (m_session->DefinesFunction("ExactSolution"))
1070 {
1071 GetFunction("ExactSolution")
1072 ->Evaluate(m_session->GetVariable(field), outfield, time);
1073 }
1074}
1075
1076/**
1077 * By default, nothing needs initialising at the EquationSystem level.
1078 */
1079void EquationSystem::v_DoInitialise(bool dumpInitialConditions)
1080{
1081 boost::ignore_unused(dumpInitialConditions);
1082}
1083
1084/**
1085 *
1086 */
1088{
1089}
1090
1091/**
1092 * Virtual function to define if operator in DoSolve is
1093 * negated with regard to the strong form. This is currently
1094 * only used in Arnoldi solves. Default is false.
1095 */
1097{
1098 return false;
1099}
1100
1101/**
1102 *
1103 */
1105{
1106}
1107
1108/**
1109 *
1110 */
1112{
1113}
1114
1115/// Virtual function for generating summary information.
1117{
1118 SessionSummary(l);
1119}
1120
1121/**
1122 * Write the field data to file. The file is named according to the
1123 * session name with the extension .fld appended.
1124 */
1126{
1127 if (!ParallelInTime())
1128 {
1129 // Serial-in-time
1130 WriteFld(m_sessionName + ".fld");
1131 }
1132 else
1133 {
1134 // Parallel-in-time
1135 std::string newdir = m_sessionName + ".pit";
1136 if (!fs::is_directory(newdir))
1137 {
1138 fs::create_directory(newdir);
1139 }
1140 WriteFld(newdir + "/" + m_sessionName + "_" +
1141 boost::lexical_cast<std::string>(
1142 m_windowPIT * m_comm->GetTimeComm()->GetSize() +
1143 m_comm->GetTimeComm()->GetRank() + 1) +
1144 ".fld");
1145 }
1146}
1147
1148/**
1149 * Zero the physical fields.
1150 */
1152{
1153 for (int i = 0; i < m_fields.size(); i++)
1154 {
1155 Vmath::Zero(m_fields[i]->GetNpoints(), m_fields[i]->UpdatePhys(), 1);
1156 }
1157}
1158
1159/**
1160 * FwdTrans the m_fields members
1161 */
1163{
1164 for (int i = 0; i < m_fields.size(); i++)
1165 {
1166 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1167 m_fields[i]->UpdateCoeffs());
1168 m_fields[i]->SetPhysState(false);
1169 }
1170}
1171
1172/**
1173 * Write the n-th checkpoint file.
1174 * @param n The index of the checkpoint file.
1175 */
1177{
1178 if (!ParallelInTime())
1179 {
1180 // Serial-in-time
1181 std::string outname =
1182 m_sessionName + "_" + boost::lexical_cast<std::string>(n);
1183 WriteFld(outname + ".chk");
1184 }
1185 else
1186 {
1187 // Parallel-in-time
1188 std::string paradir = m_sessionName + "_" +
1189 boost::lexical_cast<std::string>(m_iterPIT) +
1190 ".pit";
1191 if (!fs::is_directory(paradir))
1192 {
1193 fs::create_directory(paradir);
1194 }
1195 std::string outname = paradir + "/" + m_sessionName + "_" +
1196 boost::lexical_cast<std::string>(n);
1197 WriteFld(outname + ".chk");
1198 }
1199}
1200
1201/**
1202 * Write the n-th checkpoint file.
1203 * @param n The index of the checkpoint file.
1204 */
1206 const int n, MultiRegions::ExpListSharedPtr &field,
1207 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1208 std::vector<std::string> &variables)
1209{
1210 if (!ParallelInTime())
1211 {
1212 // Serial-in-time
1213 std::string outname =
1214 m_sessionName + "_" + boost::lexical_cast<std::string>(n);
1215 WriteFld(outname, field, fieldcoeffs, variables);
1216 }
1217 else
1218 {
1219 // Parallel-in-time
1220 std::string paradir = m_sessionName + "_" +
1221 boost::lexical_cast<std::string>(m_iterPIT) +
1222 ".pit";
1223 if (!fs::is_directory(paradir))
1224 {
1225 fs::create_directory(paradir);
1226 }
1227 std::string outname = paradir + "/" + m_sessionName + "_" +
1228 boost::lexical_cast<std::string>(n);
1229 WriteFld(outname, field, fieldcoeffs, variables);
1230 }
1231}
1232
1233/**
1234 * Write the n-th base flow into a .chk file
1235 * @param n The index of the base flow file.
1236 */
1238{
1239 std::string outname =
1240 m_sessionName + "_BaseFlow_" + boost::lexical_cast<std::string>(n);
1241
1242 WriteFld(outname + ".chk");
1243}
1244
1245/**
1246 * Writes the field data to a file with the given filename.
1247 * @param outname Filename to write to.
1248 */
1249void EquationSystem::WriteFld(const std::string &outname)
1250{
1251 std::vector<Array<OneD, NekDouble>> fieldcoeffs(m_fields.size());
1252 std::vector<std::string> variables(m_fields.size());
1253
1254 for (int i = 0; i < m_fields.size(); ++i)
1255 {
1256 if (m_fields[i]->GetNcoeffs() == m_fields[0]->GetNcoeffs())
1257 {
1258 fieldcoeffs[i] = m_fields[i]->UpdateCoeffs();
1259 }
1260 else
1261 {
1262 fieldcoeffs[i] = Array<OneD, NekDouble>(m_fields[0]->GetNcoeffs());
1263 m_fields[0]->ExtractCoeffsToCoeffs(
1264 m_fields[i], m_fields[i]->GetCoeffs(), fieldcoeffs[i]);
1265 }
1266 variables[i] = m_boundaryConditions->GetVariable(i);
1267 }
1268
1269 ExtraFldOutput(fieldcoeffs, variables);
1270
1271 WriteFld(outname, m_fields[0], fieldcoeffs, variables);
1272}
1273
1274/**
1275 * Writes the field data to a file with the given filename.
1276 * @param outname Filename to write to.
1277 * @param field ExpList on which data is based.
1278 * @param fieldcoeffs An array of array of expansion coefficients.
1279 * @param variables An array of variable names.
1280 */
1281void EquationSystem::WriteFld(const std::string &outname,
1283 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1284 std::vector<std::string> &variables)
1285{
1286 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
1287 field->GetFieldDefinitions();
1288 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
1289
1290 // Copy Data into FieldData and set variable
1291 for (int j = 0; j < fieldcoeffs.size(); ++j)
1292 {
1293 for (int i = 0; i < FieldDef.size(); ++i)
1294 {
1295 // Could do a search here to find correct variable
1296 FieldDef[i]->m_fields.push_back(variables[j]);
1297 field->AppendFieldData(FieldDef[i], FieldData[i], fieldcoeffs[j]);
1298 }
1299 }
1300
1301 // Update time in field info if required
1302 if (m_fieldMetaDataMap.find("Time") != m_fieldMetaDataMap.end())
1303 {
1304 m_fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(m_time);
1305 }
1306
1307 // Update step in field info if required
1308 if (m_fieldMetaDataMap.find("ChkFileNum") != m_fieldMetaDataMap.end())
1309 {
1310 m_fieldMetaDataMap["ChkFileNum"] =
1311 boost::lexical_cast<std::string>(m_nchk);
1312 }
1313
1314 // If necessary, add mapping information to metadata
1315 // and output mapping coordinates
1317 fields[0] = field;
1321 mapping->Output(fieldMetaDataMap, outname);
1322
1323 // If necessary, add informaton for moving frame reference to metadata
1324 if (m_fieldMetaDataMap.find("Theta_x") != m_fieldMetaDataMap.end())
1325 {
1326 // if one theta exists, add all three thetas
1327 std::vector<std::string> vSuffix = {"_x", "_y", "_z"};
1328 for (int i = 0; i < 3; ++i)
1329 {
1330 std::string sTheta = "Theta" + vSuffix[i];
1331 m_fieldMetaDataMap[sTheta] =
1332 boost::lexical_cast<std::string>(m_movingFrameTheta[i]);
1333 }
1334 }
1335
1336 m_fld->Write(outname, FieldDef, FieldData, fieldMetaDataMap,
1337 m_session->GetBackups());
1338}
1339
1340/**
1341 * Import field from infile and load into \a m_fields. This routine will
1342 * also perform a \a BwdTrans to ensure data is in both the physical and
1343 * coefficient storage.
1344 * @param infile Filename to read.
1345 */
1347 const std::string &infile,
1349{
1350 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1351 std::vector<std::vector<NekDouble>> FieldData;
1354 field_fld->Import(infile, FieldDef, FieldData);
1355
1356 // Copy FieldData into m_fields
1357 for (int j = 0; j < pFields.size(); ++j)
1358 {
1359 Vmath::Zero(pFields[j]->GetNcoeffs(), pFields[j]->UpdateCoeffs(), 1);
1360
1361 for (int i = 0; i < FieldDef.size(); ++i)
1362 {
1363 ASSERTL1(FieldDef[i]->m_fields[j] == m_session->GetVariable(j),
1364 std::string("Order of ") + infile +
1365 std::string(" data and that defined in "
1366 "m_boundaryconditions differs"));
1367
1368 pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1369 FieldDef[i]->m_fields[j],
1370 pFields[j]->UpdateCoeffs());
1371 }
1372 pFields[j]->BwdTrans(pFields[j]->GetCoeffs(), pFields[j]->UpdatePhys());
1373 }
1374}
1375
1376/**
1377 * Import field from infile and load into \a m_fields. This routine will
1378 * also perform a \a BwdTrans to ensure data is in both the physical and
1379 * coefficient storage.
1380 * @param infile Filename to read.
1381 * If optionan \a ndomains is specified it assumes we loop over nodmains
1382 * for each nvariables.
1383 */
1385 const std::string &infile,
1386 Array<OneD, MultiRegions::ExpListSharedPtr> &pFields, const int ndomains)
1387{
1388 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1389 std::vector<std::vector<NekDouble>> FieldData;
1390
1391 LibUtilities::Import(infile, FieldDef, FieldData);
1392
1393 int nvariables = GetNvariables();
1394
1395 ASSERTL0(
1396 ndomains * nvariables == pFields.size(),
1397 "Number of fields does not match the number of variables and domains");
1398
1399 // Copy FieldData into m_fields
1400 for (int j = 0; j < ndomains; ++j)
1401 {
1402 for (int i = 0; i < nvariables; ++i)
1403 {
1404 Vmath::Zero(pFields[j * nvariables + i]->GetNcoeffs(),
1405 pFields[j * nvariables + i]->UpdateCoeffs(), 1);
1406
1407 for (int n = 0; n < FieldDef.size(); ++n)
1408 {
1409 ASSERTL1(FieldDef[n]->m_fields[i] == m_session->GetVariable(i),
1410 std::string("Order of ") + infile +
1411 std::string(" data and that defined in "
1412 "m_boundaryconditions differs"));
1413
1414 pFields[j * nvariables + i]->ExtractDataToCoeffs(
1415 FieldDef[n], FieldData[n], FieldDef[n]->m_fields[i],
1416 pFields[j * nvariables + i]->UpdateCoeffs());
1417 }
1418 pFields[j * nvariables + i]->BwdTrans(
1419 pFields[j * nvariables + i]->GetCoeffs(),
1420 pFields[j * nvariables + i]->UpdatePhys());
1421 }
1422 }
1423}
1424
1425/**
1426 * Import field from infile and load into \a pField. This routine will
1427 * also perform a \a BwdTrans to ensure data is in both the physical and
1428 * coefficient storage.
1429 */
1430void EquationSystem::ImportFld(const std::string &infile,
1432 std::string &pFieldName)
1433{
1434 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1435 std::vector<std::vector<NekDouble>> FieldData;
1436
1439 field_fld->Import(infile, FieldDef, FieldData);
1440 int idx = -1;
1441
1442 Vmath::Zero(pField->GetNcoeffs(), pField->UpdateCoeffs(), 1);
1443
1444 for (int i = 0; i < FieldDef.size(); ++i)
1445 {
1446 // find the index of the required field in the file.
1447 for (int j = 0; j < FieldData.size(); ++j)
1448 {
1449 if (FieldDef[i]->m_fields[j] == pFieldName)
1450 {
1451 idx = j;
1452 }
1453 }
1454 ASSERTL1(idx >= 0, "Field " + pFieldName + " not found.");
1455
1456 pField->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1457 FieldDef[i]->m_fields[idx],
1458 pField->UpdateCoeffs());
1459 }
1460 pField->BwdTrans(pField->GetCoeffs(), pField->UpdatePhys());
1461}
1462
1463/**
1464 * Import field from infile and load into the array \a coeffs.
1465 *
1466 * @param infile Filename to read.
1467 * @param fieldStr an array of string identifying fields to be imported
1468 * @param coeffs array of array of coefficients to store imported data
1469 */
1470void EquationSystem::ImportFld(const std::string &infile,
1471 std::vector<std::string> &fieldStr,
1473{
1474
1475 ASSERTL0(fieldStr.size() <= coeffs.size(),
1476 "length of fieldstr should be the same as coeffs");
1477
1478 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1479 std::vector<std::vector<NekDouble>> FieldData;
1480
1483 field_fld->Import(infile, FieldDef, FieldData);
1484
1485 // Copy FieldData into m_fields
1486 for (int j = 0; j < fieldStr.size(); ++j)
1487 {
1488 Vmath::Zero(coeffs[j].size(), coeffs[j], 1);
1489 for (int i = 0; i < FieldDef.size(); ++i)
1490 {
1491 m_fields[0]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1492 fieldStr[j], coeffs[j]);
1493 }
1494 }
1495}
1496
1497/**
1498 * Write out a summary of the session data.
1499 * @param out Output stream to write data to.
1500 */
1502{
1503 AddSummaryItem(s, "EquationType", m_session->GetSolverInfo("EQTYPE"));
1504 AddSummaryItem(s, "Session Name", m_sessionName);
1505 AddSummaryItem(s, "Spatial Dim.", m_spacedim);
1506 AddSummaryItem(s, "Max SEM Exp. Order",
1507 m_fields[0]->EvalBasisNumModesMax());
1508
1509 if (m_session->GetComm()->GetSize() > 1)
1510 {
1511 AddSummaryItem(s, "Num. Processes", m_session->GetComm()->GetSize());
1512 }
1513
1515 {
1516 AddSummaryItem(s, "Quasi-3D", "Homogeneous in z-direction");
1517 AddSummaryItem(s, "Expansion Dim.", m_expdim + 1);
1518 AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1519 AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1520 AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1521 if (m_halfMode)
1522 {
1523 AddSummaryItem(s, "ModeType", "Half Mode");
1524 }
1525 else if (m_singleMode)
1526 {
1527 AddSummaryItem(s, "ModeType", "Single Mode");
1528 }
1529 else if (m_multipleModes)
1530 {
1531 AddSummaryItem(s, "ModeType", "Multiple Modes");
1532 }
1533 }
1535 {
1536 AddSummaryItem(s, "Quasi-3D", "Homogeneous in yz-plane");
1537 AddSummaryItem(s, "Expansion Dim.", m_expdim + 2);
1538 AddSummaryItem(s, "Num. Hom. Modes (y)", m_npointsY);
1539 AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1540 AddSummaryItem(s, "Hom. length (LY)", m_LhomY);
1541 AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1542 AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1543 }
1544 else
1545 {
1546 AddSummaryItem(s, "Expansion Dim.", m_expdim);
1547 }
1548
1549 if (m_session->DefinesSolverInfo("UpwindType"))
1550 {
1551 AddSummaryItem(s, "Riemann Solver",
1552 m_session->GetSolverInfo("UpwindType"));
1553 }
1554
1555 if (m_session->DefinesSolverInfo("AdvectionType"))
1556 {
1557 std::string AdvectionType;
1558 AdvectionType = m_session->GetSolverInfo("AdvectionType");
1560 s, "Advection Type",
1561 GetAdvectionFactory().GetClassDescription(AdvectionType));
1562 }
1563
1565 {
1566 AddSummaryItem(s, "Projection Type", "Continuous Galerkin");
1567 }
1569 {
1570 AddSummaryItem(s, "Projection Type", "Discontinuous Galerkin");
1571 }
1573 {
1574 AddSummaryItem(s, "Projection Type",
1575 "Mixed Continuous Galerkin and Discontinuous");
1576 }
1577
1578 if (m_session->DefinesSolverInfo("DiffusionType"))
1579 {
1580 std::string DiffusionType;
1581 DiffusionType = m_session->GetSolverInfo("DiffusionType");
1583 s, "Diffusion Type",
1584 GetDiffusionFactory().GetClassDescription(DiffusionType));
1585 }
1586}
1587
1588/**
1589 *
1590 */
1592{
1593 return Array<OneD, bool>(m_session->GetVariables().size(), false);
1594}
1595
1596/**
1597 *
1598 */
1600{
1601 ASSERTL0(false, "This function is not valid for the Base class");
1603 return null;
1604}
1605
1606/**
1607 *
1608 */
1610 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1611 std::vector<std::string> &variables)
1612{
1613 boost::ignore_unused(fieldcoeffs, variables);
1614}
1615
1616} // namespace SolverUtils
1617} // 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 GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:272
Describes the specification for a Basis.
Definition: Basis.h:47
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
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:197
Provides a generic Factory class.
Definition: NekFactory.hpp:105
Defines a specification for a set of points.
Definition: Points.h:55
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SOLVER_UTILS_EXPORT void ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
bool m_useFFT
Flag to determine if FFT is used for homogeneous transform.
SOLVER_UTILS_EXPORT int GetNvariables()
int m_expdim
Expansion dimension.
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
SOLVER_UTILS_EXPORT bool ParallelInTime()
Check if solver use Parallel-in-Time.
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
virtual SOLVER_UTILS_EXPORT ~EquationSystem()
Destructor.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT void ImportFld(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Input field data from the given file.
NekDouble m_time
Current time of simulation.
int m_iterPIT
Number of parallel-in-time time iteration.
bool m_multipleModes
Flag to determine if use multiple homogenenous modes are used.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
int m_windowPIT
Index of windows for parallel-in-time time iteration.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff()
Virtual function for transformation to coefficient space.
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp(void)
Virtual function to identify if operator is negated in DoSolve.
NekDouble m_fintime
Finish time of the simulation.
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
Virtual function for the L_2 error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure(void)
static std::string projectionTypeLookupIds[]
int m_npointsZ
number of points in Z direction (if homogeneous)
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
std::map< std::string, SolverUtils::SessionFunctionSharedPtr > m_sessionFunctions
Map of known SessionFunctions.
NekDouble m_checktime
Time between checkpoints.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
Virtual function for the L_inf error computation between fields and a given exact solution.
static std::string equationSystemTypeLookupIds[]
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
virtual SOLVER_UTILS_EXPORT void v_DoInitialise(bool dumpInitialConditions=true)
Virtual function for initialisation implementation.
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
Input field data from the given file to multiple domains.
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
int m_npointsY
number of points in Y direction (if homogeneous)
std::string m_sessionName
Name of the session.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n)
Write base flow file of m_fields.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
bool m_singleMode
Flag to determine if single homogeneous mode is used.
SOLVER_UTILS_EXPORT void ZeroPhysFields()
int m_HomoDirec
number of homogenous directions
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
enum HomogeneousType m_HomogeneousType
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > ErrorExtraPoints(unsigned int field)
Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, bool > m_checkIfSystemSingular
Flag to indicate if the fields should be checked for singularity.
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
SOLVER_UTILS_EXPORT int GetNcoeffs()
int m_nchk
Number of checkpoints written so far.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Virtual function for solve implementation.
virtual SOLVER_UTILS_EXPORT Array< OneD, bool > v_GetSystemSingularChecks()
int m_steps
Number of steps to take.
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
Array< OneD, NekDouble > m_movingFrameTheta
Moving frame of reference angles with respect to the.
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises EquationSystem class members.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SOLVER_UTILS_EXPORT int GetNumExpModes()
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
bool m_halfMode
Flag to determine if half homogeneous mode is used.
SOLVER_UTILS_EXPORT void FwdTransFields()
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
virtual SOLVER_UTILS_EXPORT void v_Output(void)
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys()
Virtual function for transformation to physical space.
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:53
std::vector< BasisKey > BasisKeyVector
Name for a vector of BasisKeys.
void Import(const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, FieldMetaDataMap &fieldinfomap, const Array< OneD, int > &ElementIDs)
This function allows for data to be imported from an FLD file when a session and/or communicator is n...
Definition: FieldIO.cpp:290
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
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:76
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eFourierSingleModeSpaced
1D Non Evenly-spaced points for Single Mode analysis
Definition: PointsType.h:77
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:66
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:70
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:341
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:270
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr
std::shared_ptr< ExpansionInfoMap > ExpansionInfoMapShPtr
Definition: MeshGraph.h:145
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:68
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:71