Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
Nektar::LibUtilities::SessionReader Class Reference

Reads and parses information from a Nektar++ XML session file. More...

#include <SessionReader.h>

Public Member Functions

 SessionReader (int argc, char *argv[], const std::vector< std::string > &pFilenames, const CommSharedPtr &pComm, const int &timelevel)
 
 ~SessionReader ()
 Destructor.
 
void InitSession (const std::vector< std::string > &filenames=std::vector< std::string >())
 
TiXmlDocument & GetDocument ()
 Provides direct access to the TiXmlDocument object.
 
TiXmlElement * GetElement (const std::string &pPath)
 Provides direct access to the TiXmlElement specified.
 
bool DefinesElement (const std::string &pPath) const
 Tests if a specified element is defined in the XML document.
 
const std::vector< std::string > & GetFilenames () const
 Returns the filename of the loaded XML document.
 
const std::string & GetSessionName () const
 Returns the session name of the loaded XML document.
 
CommSharedPtr GetComm ()
 Returns the communication object.
 
bool GetSharedFilesystem ()
 Returns if file system shared.
 
void Finalise ()
 Finalises the session.
 
bool DefinesParameter (const std::string &name) const
 Checks if a parameter is specified in the XML document.
 
const NekDoubleGetParameter (const std::string &pName) const
 Returns the value of the specified parameter.
 
const ParameterMapGetParameters ()
 
void LoadParameter (const std::string &name, int &var) const
 Load an integer parameter.
 
void LoadParameter (const std::string &name, unsigned int &var) const
 Load an unsigned integer parameter.
 
void LoadParameter (const std::string &name, size_t &var) const
 Load an size_t parameter.
 
void LoadParameter (const std::string &name, int &var, const int &def) const
 Check for and load an integer parameter.
 
void LoadParameter (const std::string &name, unsigned int &var, const unsigned int &def) const
 Check for and load an unsigned integer parameter.
 
void LoadParameter (const std::string &name, size_t &var, const size_t &def) const
 Check for and load an size_t parameter.
 
void LoadParameter (const std::string &name, NekDouble &var) const
 Load a double precision parameter.
 
void LoadParameter (const std::string &name, NekDouble &var, const NekDouble &def) const
 Check for and load a double-precision parameter.
 
void SetParameter (const std::string &name, int &var)
 Set an integer parameter.
 
void SetParameter (const std::string &name, unsigned int &var)
 Set an unsigned integer parameter.
 
void SetParameter (const std::string &name, size_t &var)
 Set an size_t parameter.
 
void SetParameter (const std::string &name, NekDouble &var)
 Set a double precision parameter.
 
bool DefinesSolverInfo (const std::string &name) const
 Checks if a solver info property is specified.
 
const std::string & GetSolverInfo (const std::string &pProperty) const
 Returns the value of the specified solver info property.
 
void SetSolverInfo (const std::string &pProperty, const std::string &pValue)
 Sets the value of the specified solver info property.
 
template<typename T >
const T GetSolverInfoAsEnum (const std::string &pName) const
 Returns the value of the specified solver info property as enum.
 
template<typename T >
const T GetValueAsEnum (const std::string &pName, const std::string &vValue) const
 Returns the value of the specified property and value as enum.
 
void LoadSolverInfo (const std::string &name, std::string &var, const std::string &def="") const
 Check for and load a solver info property.
 
void MatchSolverInfo (const std::string &name, const std::string &trueval, bool &var, const bool &def=false) const
 Check if the value of a solver info property matches.
 
bool MatchSolverInfo (const std::string &name, const std::string &trueval) const
 Check if the value of a solver info property matches.
 
template<typename T >
bool MatchSolverInfoAsEnum (const std::string &name, const T &trueval) const
 Check if the value of a solver info property matches.
 
bool GetBackups () const
 Returns the backups.
 
bool DefinesGlobalSysSolnInfo (const std::string &variable, const std::string &property) const
 
const std::string & GetGlobalSysSolnInfo (const std::string &variable, const std::string &property) const
 
bool DefinesTimeIntScheme () const
 Returns true if the TIMEINTEGRATIONSCHEME section is defined in the session file.
 
const TimeIntSchemeGetTimeIntScheme () const
 Returns the time integration scheme structure m_timeIntScheme from the session file.
 
std::string GetGeometryType () const
 
const std::string & GetVariable (const unsigned int &idx) const
 Returns the name of the variable specified by the given index.
 
void SetVariable (const unsigned int &idx, std::string newname)
 
std::vector< std::string > GetVariables () const
 Returns the names of all variables.
 
bool DefinesFunction (const std::string &name) const
 Checks if a specified function is defined in the XML document.
 
bool DefinesFunction (const std::string &name, const std::string &variable, const int pDomain=0) const
 Checks if a specified function has a given variable defined.
 
EquationSharedPtr GetFunction (const std::string &name, const std::string &variable, const int pDomain=0) const
 Returns an EquationSharedPtr to a given function variable.
 
EquationSharedPtr GetFunction (const std::string &name, const unsigned int &var, const int pDomain=0) const
 Returns an EquationSharedPtr to a given function variable index.
 
enum FunctionType GetFunctionType (const std::string &name, const std::string &variable, const int pDomain=0) const
 Returns the type of a given function variable.
 
enum FunctionType GetFunctionType (const std::string &pName, const unsigned int &pVar, const int pDomain=0) const
 Returns the type of a given function variable index.
 
std::string GetFunctionFilename (const std::string &name, const std::string &variable, const int pDomain=0) const
 Returns the filename to be loaded for a given variable.
 
std::string GetFunctionFilename (const std::string &name, const unsigned int &var, const int pDomain=0) const
 Returns the filename to be loaded for a given variable index.
 
std::string GetFunctionFilenameVariable (const std::string &name, const std::string &variable, const int pDomain=0) const
 Returns the filename variable to be loaded for a given variable index.
 
InterpreterSharedPtr GetInterpreter ()
 Returns the instance of the Interpreter specific to this session.
 
bool DefinesTag (const std::string &pName) const
 Checks if a specified tag is defined.
 
void SetTag (const std::string &pName, const std::string &pValue)
 Sets a specified tag.
 
const std::string & GetTag (const std::string &pName) const
 Returns the value of a specified tag.
 
const FilterMapGetFilters () const
 
bool DefinesCmdLineArgument (const std::string &pName) const
 Checks if a specified cmdline argument has been given.
 
template<typename T >
GetCmdLineArgument (const std::string &pName) const
 Retrieves a command-line argument value.
 
bool GetUpdateOptFile () const
 Get bool to update optimisation file.
 
void SetUpdateOptFile (bool flag)
 Set bool to update optimisation file.
 
size_t GetTimeLevel (void) const
 Get time level (Parallel-in-Time)
 

Static Public Member Functions

static SessionReaderSharedPtr CreateInstance (int argc, char *argv[])
 Creates an instance of the SessionReader class.
 
static SessionReaderSharedPtr CreateInstance (int argc, char *argv[], std::vector< std::string > &pFilenames, const CommSharedPtr &pComm=CommSharedPtr(), const int &timelevel=0)
 Creates an instance of the SessionReader class initialised using a separate list of XML documents.
 
static std::string RegisterEnumValue (std::string pEnum, std::string pString, int pEnumValue)
 Registers an enumeration value.
 
static std::string RegisterDefaultSolverInfo (const std::string &pName, const std::string &pValue)
 Registers the default string value of a solver info property.
 
static std::string RegisterCmdLineArgument (const std::string &pName, const std::string &pShortName, const std::string &pDescription)
 Registers a command-line argument with the session reader.
 
static std::string RegisterCmdLineFlag (const std::string &pName, const std::string &pShortName, const std::string &pDescription)
 Registers a command-line flag with the session reader.
 
static void GetXMLElementTimeLevel (TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
 Get XML elment time level (Parallel-in-Time)
 

Private Member Functions

 SessionReader (int argc, char *argv[])
 Main constructor.
 
void TestSharedFilesystem ()
 
std::vector< std::string > ParseCommandLineArguments (int argc, char *argv[])
 Parse the program arguments and fill m_cmdLineOptions.
 
std::string ParseSessionName (std::vector< std::string > &filenames)
 Parse the session name.
 
void LoadDoc (const std::string &pFilename, TiXmlDocument *pDoc) const
 Loads an xml file into a tinyxml doc and decompresses if needed.
 
TiXmlDocument * MergeDoc (const std::vector< std::string > &pFilenames) const
 Creates an XML document from a list of input files.
 
void ParseDocument ()
 Loads and parses the specified file.
 
void CreateComm (int &argc, char *argv[])
 Loads the given XML document and instantiates an appropriate communication object.
 
void PartitionComm ()
 Partitions the comm object based on session parameters.
 
void ReadParameters (TiXmlElement *conditions)
 Reads the PARAMETERS section of the XML document.
 
void ReadSolverInfo (TiXmlElement *conditions)
 Reads the SOLVERINFO section of the XML document.
 
void ReadGlobalSysSolnInfo (TiXmlElement *conditions)
 Reads the GLOBALSYSSOLNINFO section of the XML document.
 
void ReadTimeIntScheme (TiXmlElement *conditions)
 Reads the TIMEINTEGRATIONSCHEME section of the XML document.
 
void ReadVariables (TiXmlElement *conditions)
 Reads the VARIABLES section of the XML document.
 
void ReadFunctions (TiXmlElement *conditions)
 Reads the FUNCTIONS section of the XML document.
 
void ReadFilters (TiXmlElement *filters)
 Reads the FILTERS section of the XML document.
 
void CmdLineOverride ()
 Enforce parameters from command line arguments.
 
void VerifySolverInfo ()
 Check values of solver info options are valid.
 
void ParseEquals (const std::string &line, std::string &lhs, std::string &rhs)
 Parse a string in the form lhs = rhs.
 

Static Private Member Functions

static EnumMapListGetSolverInfoEnums ()
 String to enumeration map for Solver Info parameters.
 
static SolverInfoMapGetSolverInfoDefaults ()
 Default solver info options.
 
static GloSysSolnInfoListGetGloSysSolnList ()
 GlobalSysSoln Info map.
 
static CmdLineArgMapGetCmdLineArgMap ()
 CmdLine argument map.
 

Private Attributes

boost::program_options::variables_map m_cmdLineOptions
 
CommSharedPtr m_comm
 Communication object.
 
std::vector< std::string > m_filenames
 Filenames.
 
std::string m_sessionName
 Session name of the loaded XML document (filename minus ext).
 
TiXmlDocument * m_xmlDoc
 Pointer to the loaded XML document.
 
ParameterMap m_parameters
 Parameters.
 
SolverInfoMap m_solverInfo
 Solver information properties.
 
GeometricInfoMap m_geometricInfo
 Geometric information properties.
 
InterpreterSharedPtr m_interpreter
 Interpreter instance.
 
FunctionMap m_functions
 Functions.
 
VariableList m_variables
 Variables.
 
TagMap m_tags
 Custom tags.
 
FilterMap m_filters
 Filters map.
 
TimeIntScheme m_timeIntScheme
 Time integration scheme information.
 
size_t m_timeLevel = 0
 Time level.
 
bool m_verbose
 Be verbose.
 
bool m_sharedFilesystem
 Running on a shared filesystem.
 
bool m_backups = true
 Backups.
 
bool m_updateOptFile
 Update optimisation file.
 

Friends

class MemoryManager< SessionReader >
 Support creation through MemoryManager.
 

Detailed Description

Reads and parses information from a Nektar++ XML session file.

This class provides an interface to Nektar++-specific content in a supplied XML document. It also initialises a Nektar++ session including setting up communication for parallel execution and where necessary partitioning the supplied mesh for running across multiple processes.

A session should be initialised at the beginning of a user's application by passing the command-line arguments. This not only allows the SessionReader to extract the name of the XML document to load containing Nektar++ session information, but also supplies the MPI arguments necessary for setting up parallel communication. The SessionReader should be initialised using the CreateInstance function:

static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
std::shared_ptr< SessionReader > SessionReaderSharedPtr

The instance vSession can now be passed to other key Nektar++ components during their construction.

Note
At the end of the user application, it is important to call the Finalise routine in order to finalise any MPI communication and correctly free resources.

The SessionReader class provides streamlined, validated access to session parameters, solver information and functions defined within a Nektar++ XML document. The available routines and their usage is documented below.

In the case of solver information properties, the classes to which these parameters are pertinent may register with the SessionReader class the set of valid values for a given property. Such values may also be associated with an enumeration value for more transparent use of the property values in code.

Definition at line 122 of file SessionReader.h.

Constructor & Destructor Documentation

◆ SessionReader() [1/2]

Nektar::LibUtilities::SessionReader::SessionReader ( int  argc,
char *  argv[],
const std::vector< std::string > &  pFilenames,
const CommSharedPtr pComm,
const int &  timelevel 
)

Definition at line 215 of file BasicUtils/SessionReader.cpp.

218{
219 ASSERTL0(pFilenames.size() > 0, "No filenames specified.");
220
221 ParseCommandLineArguments(argc, argv);
222 m_xmlDoc = nullptr;
223 m_filenames = pFilenames;
224
226
227 // Create communicator
228 if (!pComm.get())
229 {
230 CreateComm(argc, argv);
231 }
232 else
233 {
234 m_comm = pComm;
235 }
236
238
239 // Split up the communicator
241
242 // If running in parallel change the default global sys solution
243 // type.
244 if (m_comm->GetSpaceComm()->GetSize() > 1)
245 {
246 GetSolverInfoDefaults()["GLOBALSYSSOLN"] = "IterativeStaticCond";
247 }
248
250 m_interpreter->SetRandomSeed((m_comm->GetSpaceComm()->GetRank() + 1) *
251 (unsigned int)time(nullptr));
252
253 // Set time level (Parallel-in-Time)
254 m_timeLevel = timelevel;
255}
#define ASSERTL0(condition, msg)
std::vector< std::string > ParseCommandLineArguments(int argc, char *argv[])
Parse the program arguments and fill m_cmdLineOptions.
static SolverInfoMap & GetSolverInfoDefaults()
Default solver info options.
std::string ParseSessionName(std::vector< std::string > &filenames)
Parse the session name.
InterpreterSharedPtr m_interpreter
Interpreter instance.
std::vector< std::string > m_filenames
Filenames.
std::string m_sessionName
Session name of the loaded XML document (filename minus ext).
CommSharedPtr m_comm
Communication object.
void CreateComm(int &argc, char *argv[])
Loads the given XML document and instantiates an appropriate communication object.
void PartitionComm()
Partitions the comm object based on session parameters.
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, CreateComm(), GetSolverInfoDefaults(), m_comm, m_filenames, m_interpreter, m_sessionName, m_timeLevel, m_xmlDoc, ParseCommandLineArguments(), ParseSessionName(), PartitionComm(), and TestSharedFilesystem().

◆ ~SessionReader()

Nektar::LibUtilities::SessionReader::~SessionReader ( )

Destructor.

Definition at line 260 of file BasicUtils/SessionReader.cpp.

261{
262 if (m_xmlDoc)
263 {
264 delete m_xmlDoc;
265 }
266}

References m_xmlDoc.

◆ SessionReader() [2/2]

Nektar::LibUtilities::SessionReader::SessionReader ( int  argc,
char *  argv[] 
)
private

Main constructor.

This constructor parses the command-line arguments given to the user application to set up any MPI communication, read supplied XML session files, and partition meshes where necessary.

Parameters
argcNumber of command-line arguments
argvArray of command-line arguments

Definition at line 183 of file BasicUtils/SessionReader.cpp.

184{
185 m_xmlDoc = nullptr;
187
188 ASSERTL0(m_filenames.size() > 0, "No session file(s) given.");
189
191
192 // Create communicator
193 CreateComm(argc, argv);
194
196
197 // Split up the communicator
199
200 // If running in parallel change the default global sys solution
201 // type.
202 if (m_comm->GetSpaceComm()->GetSize() > 1)
203 {
204 GetSolverInfoDefaults()["GLOBALSYSSOLN"] = "IterativeStaticCond";
205 }
206
208 m_interpreter->SetRandomSeed((m_comm->GetSpaceComm()->GetRank() + 1) *
209 (unsigned int)time(nullptr));
210}

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, CreateComm(), GetSolverInfoDefaults(), m_comm, m_filenames, m_interpreter, m_sessionName, m_xmlDoc, ParseCommandLineArguments(), ParseSessionName(), PartitionComm(), and TestSharedFilesystem().

Member Function Documentation

◆ CmdLineOverride()

void Nektar::LibUtilities::SessionReader::CmdLineOverride ( )
private

Enforce parameters from command line arguments.

Definition at line 2543 of file BasicUtils/SessionReader.cpp.

2544{
2545 // Parse solver info overrides
2546 if (m_cmdLineOptions.count("solverinfo"))
2547 {
2548 std::vector<std::string> solverInfoList =
2549 m_cmdLineOptions["solverinfo"].as<std::vector<std::string>>();
2550
2551 for (size_t i = 0; i < solverInfoList.size(); ++i)
2552 {
2553 std::string lhs, rhs;
2554
2555 try
2556 {
2557 ParseEquals(solverInfoList[i], lhs, rhs);
2558 }
2559 catch (...)
2560 {
2561 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2562 "option: " +
2563 solverInfoList[i]);
2564 }
2565
2566 std::string lhsUpper = boost::to_upper_copy(lhs);
2567 m_solverInfo[lhsUpper] = rhs;
2568 }
2569 }
2570
2571 if (m_cmdLineOptions.count("parameter"))
2572 {
2573 std::vector<std::string> parametersList =
2574 m_cmdLineOptions["parameter"].as<std::vector<std::string>>();
2575
2576 for (size_t i = 0; i < parametersList.size(); ++i)
2577 {
2578 std::string lhs, rhs;
2579
2580 try
2581 {
2582 ParseEquals(parametersList[i], lhs, rhs);
2583 }
2584 catch (...)
2585 {
2586 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2587 "option: " +
2588 parametersList[i]);
2589 }
2590
2591 std::string lhsUpper = boost::to_upper_copy(lhs);
2592
2593 try
2594 {
2595 m_parameters[lhsUpper] = std::stod(rhs);
2596 }
2597 catch (...)
2598 {
2599 NEKERROR(ErrorUtil::efatal, "Unable to convert string: " + rhs +
2600 "to double value.");
2601 }
2602 }
2603 }
2604}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
SolverInfoMap m_solverInfo
Solver information properties.
ParameterMap m_parameters
Parameters.
boost::program_options::variables_map m_cmdLineOptions
void ParseEquals(const std::string &line, std::string &lhs, std::string &rhs)
Parse a string in the form lhs = rhs.

References Nektar::ErrorUtil::efatal, m_cmdLineOptions, m_parameters, m_solverInfo, NEKERROR, and ParseEquals().

Referenced by InitSession().

◆ CreateComm()

void Nektar::LibUtilities::SessionReader::CreateComm ( int &  argc,
char *  argv[] 
)
private

Loads the given XML document and instantiates an appropriate communication object.

Definition at line 1661 of file BasicUtils/SessionReader.cpp.

1662{
1663 if (argc == 0)
1664 {
1665 m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1666 }
1667 else
1668 {
1669 string vCommModule("Serial");
1670 if (GetCommFactory().ModuleExists("ParallelMPI"))
1671 {
1672 vCommModule = "ParallelMPI";
1673 }
1674 if (m_cmdLineOptions.count("cwipi") &&
1675 GetCommFactory().ModuleExists("CWIPI"))
1676 {
1677 vCommModule = "CWIPI";
1678 }
1679
1680 m_comm = GetCommFactory().CreateInstance(vCommModule, argc, argv);
1681 }
1682}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
CommFactory & GetCommFactory()

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::GetCommFactory(), m_cmdLineOptions, and m_comm.

Referenced by SessionReader(), and SessionReader().

◆ CreateInstance() [1/2]

static SessionReaderSharedPtr Nektar::LibUtilities::SessionReader::CreateInstance ( int  argc,
char *  argv[] 
)
inlinestatic

Creates an instance of the SessionReader class.

This function should be used by an application to instantiate the session reader. It should be called at the very beginning of the application before any other processing of command-line arguments. After instantiating the class and setting up any parallel communication, it also calls the main initialisation of the object.

Definition at line 138 of file SessionReader.h.

140 {
143 argv);
144 return p;
145 }
std::vector< double > p(NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

Referenced by Nektar::MovementTests::CreateSession(), Diffusion::Diffusion(), main(), SessionReader_CreateInstance(), Nektar::SolverUtils::DriverParallelInTime::SetParallelInTimeEquationSystem(), Nektar::VarcoeffHashingTest::setupContFieldSolve(), Nektar::SolverUtils::Driver::v_InitObject(), Nektar::FieldUtils::InputXml::v_Process(), Nektar::FieldUtils::ProcessDisplacement::v_Process(), Nektar::FieldUtils::ProcessHomogeneousPlane::v_Process(), Nektar::FieldUtils::ProcessInterpField::v_Process(), Nektar::FieldUtils::ProcessInterpPoints::v_Process(), and Nektar::VortexWaveInteraction::VortexWaveInteraction().

◆ CreateInstance() [2/2]

static SessionReaderSharedPtr Nektar::LibUtilities::SessionReader::CreateInstance ( int  argc,
char *  argv[],
std::vector< std::string > &  pFilenames,
const CommSharedPtr pComm = CommSharedPtr(),
const int &  timelevel = 0 
)
inlinestatic

Creates an instance of the SessionReader class initialised using a separate list of XML documents.

This function should be used by an application to instantiate the session reader. It may be called after processing of command-line arguments. After instantiating the class and setting up any parallel communication, it also calls the main initialisation of the object.

Definition at line 157 of file SessionReader.h.

160 {
163 argc, argv, pFilenames, pComm, timelevel);
164 return p;
165 }

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ DefinesCmdLineArgument()

bool Nektar::LibUtilities::SessionReader::DefinesCmdLineArgument ( const std::string &  pName) const

Checks if a specified cmdline argument has been given.

Definition at line 1466 of file BasicUtils/SessionReader.cpp.

1467{
1468 return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1469}

References m_cmdLineOptions.

Referenced by InitSession(), and PartitionComm().

◆ DefinesElement()

bool Nektar::LibUtilities::SessionReader::DefinesElement ( const std::string &  pPath) const

Tests if a specified element is defined in the XML document.

Definition at line 722 of file BasicUtils/SessionReader.cpp.

723{
724 std::string vPath = boost::to_upper_copy(pPath);
725 std::vector<std::string> st;
726 boost::split(st, vPath, boost::is_any_of("\\/ "));
727 ASSERTL0(st.size() > 0, "No path given in XML element request.");
728
729 TiXmlElement *vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
730 ASSERTL0(vReturn,
731 std::string("Cannot find element '") + st[0] + std::string("'."));
732 for (int i = 1; i < st.size(); ++i)
733 {
734 vReturn = vReturn->FirstChildElement(st[i].c_str());
735 if (!vReturn)
736 {
737 return false;
738 }
739 }
740 return true;
741}

References ASSERTL0, and m_xmlDoc.

◆ DefinesFunction() [1/2]

bool Nektar::LibUtilities::SessionReader::DefinesFunction ( const std::string &  name) const

Checks if a specified function is defined in the XML document.

Definition at line 1214 of file BasicUtils/SessionReader.cpp.

1215{
1216 std::string vName = boost::to_upper_copy(pName);
1217 return m_functions.find(vName) != m_functions.end();
1218}

References m_functions.

◆ DefinesFunction() [2/2]

bool Nektar::LibUtilities::SessionReader::DefinesFunction ( const std::string &  name,
const std::string &  variable,
const int  pDomain = 0 
) const

Checks if a specified function has a given variable defined.

Definition at line 1223 of file BasicUtils/SessionReader.cpp.

1226{
1227 std::string vName = boost::to_upper_copy(pName);
1228
1229 // Check function exists
1230 auto it1 = m_functions.find(vName);
1231 if (it1 != m_functions.end())
1232 {
1233 pair<std::string, int> key(pVariable, pDomain);
1234 pair<std::string, int> defkey("*", pDomain);
1235 bool varExists = it1->second.find(key) != it1->second.end() ||
1236 it1->second.find(defkey) != it1->second.end();
1237 return varExists;
1238 }
1239 return false;
1240}

References m_functions.

◆ DefinesGlobalSysSolnInfo()

bool Nektar::LibUtilities::SessionReader::DefinesGlobalSysSolnInfo ( const std::string &  variable,
const std::string &  property 
) const

Definition at line 1082 of file BasicUtils/SessionReader.cpp.

1084{
1085 auto iter = GetGloSysSolnList().find(pVariable);
1086 if (iter == GetGloSysSolnList().end())
1087 {
1088 return false;
1089 }
1090
1091 std::string vProperty = boost::to_upper_copy(pProperty);
1092
1093 auto iter1 = iter->second.find(vProperty);
1094 if (iter1 == iter->second.end())
1095 {
1096 return false;
1097 }
1098
1099 return true;
1100}
static GloSysSolnInfoList & GetGloSysSolnList()
GlobalSysSoln Info map.

References GetGloSysSolnList().

◆ DefinesParameter()

bool Nektar::LibUtilities::SessionReader::DefinesParameter ( const std::string &  name) const

Checks if a parameter is specified in the XML document.

Definition at line 786 of file BasicUtils/SessionReader.cpp.

787{
788 std::string vName = boost::to_upper_copy(pName);
789 return m_parameters.find(vName) != m_parameters.end();
790}

References m_parameters.

Referenced by export_SessionReader().

◆ DefinesSolverInfo()

bool Nektar::LibUtilities::SessionReader::DefinesSolverInfo ( const std::string &  name) const

Checks if a solver info property is specified.

Definition at line 987 of file BasicUtils/SessionReader.cpp.

988{
989 std::string vName = boost::to_upper_copy(pName);
990 auto infoIter = m_solverInfo.find(vName);
991 return (infoIter != m_solverInfo.end());
992}

References m_solverInfo.

Referenced by export_SessionReader(), GetSolverInfoAsEnum(), and MatchSolverInfo().

◆ DefinesTag()

bool Nektar::LibUtilities::SessionReader::DefinesTag ( const std::string &  pName) const

Checks if a specified tag is defined.

Definition at line 1429 of file BasicUtils/SessionReader.cpp.

1430{
1431 std::string vName = boost::to_upper_copy(pName);
1432 return m_tags.find(vName) != m_tags.end();
1433}

References m_tags.

◆ DefinesTimeIntScheme()

bool Nektar::LibUtilities::SessionReader::DefinesTimeIntScheme ( ) const

Returns true if the TIMEINTEGRATIONSCHEME section is defined in the session file.

Definition at line 1126 of file BasicUtils/SessionReader.cpp.

1127{
1128 return m_timeIntScheme.method != "";
1129}
TimeIntScheme m_timeIntScheme
Time integration scheme information.

References m_timeIntScheme, and Nektar::LibUtilities::TimeIntScheme::method.

◆ Finalise()

void Nektar::LibUtilities::SessionReader::Finalise ( )

Finalises the session.

This routine finalises any parallel communication.

Note
This routine should be called at the very end of a users application.

Definition at line 778 of file BasicUtils/SessionReader.cpp.

779{
780 m_comm->Finalise();
781}

References m_comm.

Referenced by export_SessionReader().

◆ GetBackups()

bool Nektar::LibUtilities::SessionReader::GetBackups ( ) const

Returns the backups.

Definition at line 1206 of file BasicUtils/SessionReader.cpp.

1207{
1208 return m_backups;
1209}

References m_backups.

◆ GetCmdLineArgMap()

CmdLineArgMap & Nektar::LibUtilities::SessionReader::GetCmdLineArgMap ( )
staticprivate

CmdLine argument map.

Lists the possible command-line argument which can be specified for this executable.

This list is populated through the RegisterCmdLineArgument static member function which is called statically from various classes to register command-line arguments they need.

Definition at line 169 of file BasicUtils/SessionReader.cpp.

170{
171 static CmdLineArgMap cmdLineArguments;
172 return cmdLineArguments;
173}
std::map< std::string, CmdLineArg > CmdLineArgMap

Referenced by ParseCommandLineArguments(), RegisterCmdLineArgument(), and RegisterCmdLineFlag().

◆ GetCmdLineArgument()

template<typename T >
T Nektar::LibUtilities::SessionReader::GetCmdLineArgument ( const std::string &  pName) const
inline

Retrieves a command-line argument value.

Definition at line 364 of file SessionReader.h.

365 {
366 return m_cmdLineOptions.find(pName)->second.as<T>();
367 }

References m_cmdLineOptions.

◆ GetComm()

CommSharedPtr Nektar::LibUtilities::SessionReader::GetComm ( )

Returns the communication object.

Definition at line 762 of file BasicUtils/SessionReader.cpp.

763{
764 return m_comm;
765}

References m_comm.

Referenced by export_SessionReader().

◆ GetDocument()

TiXmlDocument & Nektar::LibUtilities::SessionReader::GetDocument ( )

Provides direct access to the TiXmlDocument object.

Definition at line 672 of file BasicUtils/SessionReader.cpp.

673{
674 ASSERTL1(m_xmlDoc, "XML Document not defined.");
675 return *m_xmlDoc;
676}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....

References ASSERTL1, and m_xmlDoc.

◆ GetElement()

TiXmlElement * Nektar::LibUtilities::SessionReader::GetElement ( const std::string &  pPath)

Provides direct access to the TiXmlElement specified.

The single parameter specifies a path to the requested element in a similar format to the filesystem path. Given the following XML:

<NEKTAR>
<CONDITIONS>
<PARAMETERS>
...
</PARAMETERS>
</CONDITIONS>
</NEKTAR>

the PARAMETERS element would be retrieved by requesting the path:

Nektar/Conditions/Parameters
Note
Paths are case-insensitive.
Parameters
pPathPath to requested element.
Returns
Direct pointer to requested XML Element.

Definition at line 700 of file BasicUtils/SessionReader.cpp.

701{
702 std::string vPath = boost::to_upper_copy(pPath);
703 std::vector<std::string> st;
704 boost::split(st, vPath, boost::is_any_of("\\/ "));
705 ASSERTL0(st.size() > 0, "No path given in XML element request.");
706
707 TiXmlElement *vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
708 ASSERTL0(vReturn,
709 std::string("Cannot find element '") + st[0] + std::string("'."));
710 for (int i = 1; i < st.size(); ++i)
711 {
712 vReturn = vReturn->FirstChildElement(st[i].c_str());
713 ASSERTL0(vReturn, std::string("Cannot find element '") + st[i] +
714 std::string("'."));
715 }
716 return vReturn;
717}

References ASSERTL0, and m_xmlDoc.

◆ GetFilenames()

const std::vector< std::string > & Nektar::LibUtilities::SessionReader::GetFilenames ( ) const

Returns the filename of the loaded XML document.

Definition at line 746 of file BasicUtils/SessionReader.cpp.

747{
748 return m_filenames;
749}

References m_filenames.

◆ GetFilters()

const FilterMap & Nektar::LibUtilities::SessionReader::GetFilters ( ) const

Definition at line 1458 of file BasicUtils/SessionReader.cpp.

1459{
1460 return m_filters;
1461}

References m_filters.

◆ GetFunction() [1/2]

EquationSharedPtr Nektar::LibUtilities::SessionReader::GetFunction ( const std::string &  name,
const std::string &  variable,
const int  pDomain = 0 
) const

Returns an EquationSharedPtr to a given function variable.

Definition at line 1245 of file BasicUtils/SessionReader.cpp.

1248{
1249 std::string vName = boost::to_upper_copy(pName);
1250 auto it1 = m_functions.find(vName);
1251
1252 ASSERTL0(it1 != m_functions.end(),
1253 std::string("No such function '") + pName +
1254 std::string("' has been defined in the session file."));
1255
1256 // Check for specific and wildcard definitions
1257 pair<std::string, int> key(pVariable, pDomain);
1258 pair<std::string, int> defkey("*", pDomain);
1259
1260 auto it2 = it1->second.find(key);
1261 auto it3 = it1->second.find(defkey);
1262 bool specific = it2 != it1->second.end();
1263 bool wildcard = it3 != it1->second.end();
1264
1265 // Check function is defined somewhere
1266 ASSERTL0(specific || wildcard, "No such variable " + pVariable +
1267 " in domain " + std::to_string(pDomain) +
1268 " defined for function " + pName +
1269 " in session file.");
1270
1271 // If not specific, must be wildcard
1272 if (!specific)
1273 {
1274 it2 = it3;
1275 }
1276
1277 ASSERTL0((it2->second.m_type == eFunctionTypeExpression),
1278 std::string("Function is defined by a file."));
1279 return it2->second.m_expression;
1280}

References ASSERTL0, Nektar::LibUtilities::eFunctionTypeExpression, and m_functions.

Referenced by export_SessionReader(), and GetFunction().

◆ GetFunction() [2/2]

EquationSharedPtr Nektar::LibUtilities::SessionReader::GetFunction ( const std::string &  name,
const unsigned int &  var,
const int  pDomain = 0 
) const

Returns an EquationSharedPtr to a given function variable index.

Definition at line 1285 of file BasicUtils/SessionReader.cpp.

1288{
1289 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1290 return GetFunction(pName, m_variables[pVar], pDomain);
1291}
VariableList m_variables
Variables.
EquationSharedPtr GetFunction(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns an EquationSharedPtr to a given function variable.

References ASSERTL0, GetFunction(), and m_variables.

◆ GetFunctionFilename() [1/2]

std::string Nektar::LibUtilities::SessionReader::GetFunctionFilename ( const std::string &  name,
const std::string &  variable,
const int  pDomain = 0 
) const

Returns the filename to be loaded for a given variable.

Definition at line 1344 of file BasicUtils/SessionReader.cpp.

1347{
1348 std::string vName = boost::to_upper_copy(pName);
1349 auto it1 = m_functions.find(vName);
1350
1351 ASSERTL0(it1 != m_functions.end(),
1352 std::string("Function '") + pName + std::string("' not found."));
1353
1354 // Check for specific and wildcard definitions
1355 pair<std::string, int> key(pVariable, pDomain);
1356 pair<std::string, int> defkey("*", pDomain);
1357
1358 auto it2 = it1->second.find(key);
1359 auto it3 = it1->second.find(defkey);
1360 bool specific = it2 != it1->second.end();
1361 bool wildcard = it3 != it1->second.end();
1362
1363 // Check function is defined somewhere
1364 ASSERTL0(specific || wildcard, "No such variable " + pVariable +
1365 " in domain " + std::to_string(pDomain) +
1366 " defined for function " + pName +
1367 " in session file.");
1368
1369 // If not specific, must be wildcard
1370 if (!specific)
1371 {
1372 it2 = it3;
1373 }
1374
1375 return it2->second.m_filename;
1376}

References ASSERTL0, and m_functions.

Referenced by GetFunctionFilename().

◆ GetFunctionFilename() [2/2]

std::string Nektar::LibUtilities::SessionReader::GetFunctionFilename ( const std::string &  name,
const unsigned int &  var,
const int  pDomain = 0 
) const

Returns the filename to be loaded for a given variable index.

Definition at line 1381 of file BasicUtils/SessionReader.cpp.

1384{
1385 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1386 return GetFunctionFilename(pName, m_variables[pVar], pDomain);
1387}
std::string GetFunctionFilename(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns the filename to be loaded for a given variable.

References ASSERTL0, GetFunctionFilename(), and m_variables.

◆ GetFunctionFilenameVariable()

std::string Nektar::LibUtilities::SessionReader::GetFunctionFilenameVariable ( const std::string &  name,
const std::string &  variable,
const int  pDomain = 0 
) const

Returns the filename variable to be loaded for a given variable index.

Definition at line 1392 of file BasicUtils/SessionReader.cpp.

1395{
1396 std::string vName = boost::to_upper_copy(pName);
1397 auto it1 = m_functions.find(vName);
1398
1399 ASSERTL0(it1 != m_functions.end(),
1400 std::string("Function '") + pName + std::string("' not found."));
1401
1402 // Check for specific and wildcard definitions
1403 pair<std::string, int> key(pVariable, pDomain);
1404 pair<std::string, int> defkey("*", pDomain);
1405
1406 auto it2 = it1->second.find(key);
1407 auto it3 = it1->second.find(defkey);
1408 bool specific = it2 != it1->second.end();
1409 bool wildcard = it3 != it1->second.end();
1410
1411 // Check function is defined somewhere
1412 ASSERTL0(specific || wildcard, "No such variable " + pVariable +
1413 " in domain " + std::to_string(pDomain) +
1414 " defined for function " + pName +
1415 " in session file.");
1416
1417 // If not specific, must be wildcard
1418 if (!specific)
1419 {
1420 it2 = it3;
1421 }
1422
1423 return it2->second.m_fileVariable;
1424}

References ASSERTL0, and m_functions.

◆ GetFunctionType() [1/2]

enum FunctionType Nektar::LibUtilities::SessionReader::GetFunctionType ( const std::string &  name,
const std::string &  variable,
const int  pDomain = 0 
) const

Returns the type of a given function variable.

Definition at line 1296 of file BasicUtils/SessionReader.cpp.

1299{
1300 std::string vName = boost::to_upper_copy(pName);
1301 auto it1 = m_functions.find(vName);
1302
1303 ASSERTL0(it1 != m_functions.end(),
1304 std::string("Function '") + pName + std::string("' not found."));
1305
1306 // Check for specific and wildcard definitions
1307 pair<std::string, int> key(pVariable, pDomain);
1308 pair<std::string, int> defkey("*", pDomain);
1309
1310 auto it2 = it1->second.find(key);
1311 auto it3 = it1->second.find(defkey);
1312 bool specific = it2 != it1->second.end();
1313 bool wildcard = it3 != it1->second.end();
1314
1315 // Check function is defined somewhere
1316 ASSERTL0(specific || wildcard, "No such variable " + pVariable +
1317 " in domain " + std::to_string(pDomain) +
1318 " defined for function " + pName +
1319 " in session file.");
1320
1321 // If not specific, must be wildcard
1322 if (!specific)
1323 {
1324 it2 = it3;
1325 }
1326
1327 return it2->second.m_type;
1328}

References ASSERTL0, and m_functions.

Referenced by GetFunctionType().

◆ GetFunctionType() [2/2]

enum FunctionType Nektar::LibUtilities::SessionReader::GetFunctionType ( const std::string &  pName,
const unsigned int &  pVar,
const int  pDomain = 0 
) const

Returns the type of a given function variable index.

Definition at line 1333 of file BasicUtils/SessionReader.cpp.

1336{
1337 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1338 return GetFunctionType(pName, m_variables[pVar], pDomain);
1339}
enum FunctionType GetFunctionType(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns the type of a given function variable.

References ASSERTL0, GetFunctionType(), and m_variables.

◆ GetGeometryType()

std::string Nektar::LibUtilities::SessionReader::GetGeometryType ( ) const

Definition at line 1143 of file BasicUtils/SessionReader.cpp.

1144{
1145 TiXmlElement *xmlGeom =
1146 m_xmlDoc->FirstChildElement("NEKTAR")->FirstChildElement("GEOMETRY");
1147 ASSERTL1(xmlGeom, "Failed to find a GEOMETRY section in m_xmlDoc");
1148
1149 TiXmlAttribute *attr = xmlGeom->FirstAttribute();
1150 while (attr)
1151 {
1152 std::string attrName(attr->Name());
1153 if (attrName == "HDF5FILE")
1154 {
1155 // there is a file pointer, therefore is HDF5
1156 return "HDF5";
1157 }
1158 // Get the next attribute.
1159 attr = attr->Next();
1160 }
1161
1162 // Check the VERTEX block. If this is compressed, assume the file is
1163 // compressed, otherwise assume uncompressed.
1164 TiXmlElement *element = xmlGeom->FirstChildElement("VERTEX");
1165 string IsCompressed;
1166 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
1167
1168 if (IsCompressed.size() > 0)
1169 {
1170 return "XmlCompressed";
1171 }
1172
1173 // no file pointer or compressed, just standard xml
1174 return "Xml";
1175}

References ASSERTL1, and m_xmlDoc.

◆ GetGlobalSysSolnInfo()

const std::string & Nektar::LibUtilities::SessionReader::GetGlobalSysSolnInfo ( const std::string &  variable,
const std::string &  property 
) const

Definition at line 1105 of file BasicUtils/SessionReader.cpp.

1107{
1108 auto iter = GetGloSysSolnList().find(pVariable);
1109 ASSERTL0(iter != GetGloSysSolnList().end(),
1110 "Failed to find variable in GlobalSysSolnInfoList");
1111
1112 std::string vProperty = boost::to_upper_copy(pProperty);
1113 auto iter1 = iter->second.find(vProperty);
1114
1115 ASSERTL0(iter1 != iter->second.end(),
1116 "Failed to find property: " + vProperty +
1117 " in GlobalSysSolnInfoList");
1118
1119 return iter1->second;
1120}

References ASSERTL0, and GetGloSysSolnList().

◆ GetGloSysSolnList()

GloSysSolnInfoList & Nektar::LibUtilities::SessionReader::GetGloSysSolnList ( )
staticprivate

GlobalSysSoln Info map.

List of values for GlobalSysSoln parameters to be used to override details given in SolverInfo

This list is populated by ReadGlobalSysSolnInfo if the GLOBALSYSSOLNINFO section is defined in the input file. This List allows for details to define for the Global Sys solver for each variable.

Definition at line 155 of file BasicUtils/SessionReader.cpp.

156{
157 static GloSysSolnInfoList gloSysSolnInfoList;
158 return gloSysSolnInfoList;
159}
std::map< std::string, GloSysInfoMap > GloSysSolnInfoList

Referenced by DefinesGlobalSysSolnInfo(), GetGlobalSysSolnInfo(), and ReadGlobalSysSolnInfo().

◆ GetInterpreter()

InterpreterSharedPtr Nektar::LibUtilities::SessionReader::GetInterpreter ( )

Returns the instance of the Interpreter specific to this session.

Definition at line 2631 of file BasicUtils/SessionReader.cpp.

2632{
2633 return m_interpreter;
2634}

References m_interpreter.

◆ GetParameter()

const NekDouble & Nektar::LibUtilities::SessionReader::GetParameter ( const std::string &  pName) const

Returns the value of the specified parameter.

If the parameter is not defined, termination occurs. Therefore, the parameters existence should be tested for using DefinesParameter before calling this function.

Parameters
pNameThe name of a floating-point parameter.
Returns
The value of the floating-point parameter.

Definition at line 800 of file BasicUtils/SessionReader.cpp.

801{
802 std::string vName = boost::to_upper_copy(pName);
803 auto paramIter = m_parameters.find(vName);
804
805 ASSERTL0(paramIter != m_parameters.end(),
806 "Unable to find requested parameter: " + pName);
807
808 return paramIter->second;
809}

References ASSERTL0, and m_parameters.

Referenced by export_SessionReader().

◆ GetParameters()

const ParameterMap & Nektar::LibUtilities::SessionReader::GetParameters ( )

Getter for the session's parameters

Returns
A reference to the parameter map

Definition at line 815 of file BasicUtils/SessionReader.cpp.

816{
817 return m_parameters;
818}

References m_parameters.

Referenced by export_SessionReader().

◆ GetSessionName()

const std::string & Nektar::LibUtilities::SessionReader::GetSessionName ( ) const

Returns the session name of the loaded XML document.

Definition at line 754 of file BasicUtils/SessionReader.cpp.

755{
756 return m_sessionName;
757}

References m_sessionName.

Referenced by export_SessionReader().

◆ GetSharedFilesystem()

bool Nektar::LibUtilities::SessionReader::GetSharedFilesystem ( )

Returns if file system shared.

Definition at line 767 of file BasicUtils/SessionReader.cpp.

768{
769 return m_sharedFilesystem;
770}
bool m_sharedFilesystem
Running on a shared filesystem.

References m_sharedFilesystem.

Referenced by export_SessionReader().

◆ GetSolverInfo()

const std::string & Nektar::LibUtilities::SessionReader::GetSolverInfo ( const std::string &  pProperty) const

Returns the value of the specified solver info property.

Definition at line 997 of file BasicUtils/SessionReader.cpp.

999{
1000 std::string vProperty = boost::to_upper_copy(pProperty);
1001 auto iter = m_solverInfo.find(vProperty);
1002
1003 ASSERTL1(iter != m_solverInfo.end(),
1004 "Unable to find requested property: " + pProperty);
1005
1006 return iter->second;
1007}

References ASSERTL1, and m_solverInfo.

Referenced by export_SessionReader(), and GetSolverInfoAsEnum().

◆ GetSolverInfoAsEnum()

template<typename T >
const T Nektar::LibUtilities::SessionReader::GetSolverInfoAsEnum ( const std::string &  pName) const
inline

Returns the value of the specified solver info property as enum.

Definition at line 513 of file SessionReader.h.

515{
516 std::string vName = boost::to_upper_copy(pName);
518 "Solver info '" + pName + "' not defined.");
519
520 std::string vValue = GetSolverInfo(vName);
521 auto x = GetSolverInfoEnums().find(vName);
522 ASSERTL0(x != GetSolverInfoEnums().end(),
523 "Enum for SolverInfo property '" + pName + "' not found.");
524
525 auto y = x->second.find(vValue);
526 ASSERTL0(y != x->second.end(),
527 "Value of SolverInfo property '" + pName + "' is invalid.");
528
529 return T(y->second);
530}
bool DefinesSolverInfo(const std::string &name) const
Checks if a solver info property is specified.
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.
const std::string & GetSolverInfo(const std::string &pProperty) const
Returns the value of the specified solver info property.

References ASSERTL0, DefinesSolverInfo(), GetSolverInfo(), and GetSolverInfoEnums().

◆ GetSolverInfoDefaults()

SolverInfoMap & Nektar::LibUtilities::SessionReader::GetSolverInfoDefaults ( )
staticprivate

Default solver info options.

List of default values for solver information parameters to be used in the case of them not being provided.

This list is populated through the RegisterDefaultSolverInfo static member variable which is called statically from various classes to register the default value for a given parameter.

Definition at line 140 of file BasicUtils/SessionReader.cpp.

141{
142 static SolverInfoMap solverInfoMap;
143 return solverInfoMap;
144}
std::map< std::string, std::string > SolverInfoMap

Referenced by ReadSolverInfo(), RegisterDefaultSolverInfo(), SessionReader(), and SessionReader().

◆ GetSolverInfoEnums()

EnumMapList & Nektar::LibUtilities::SessionReader::GetSolverInfoEnums ( )
staticprivate

String to enumeration map for Solver Info parameters.

This map of maps stores the list of valid string values for a number of solver information parameters. The top level map connects different parameter names to their list of possible values. The list of possible values is also a map, mapping a valid string to a corresponding enum value.

This list is populated through the RegisterEnumValue static member function which is called statically from various classes to register the valid values for solver info parameters associated with them. The map is therefore fully populated before the SessionReader class is instantiated and a file is read in and parsed.

Definition at line 126 of file BasicUtils/SessionReader.cpp.

127{
128 static EnumMapList solverInfoEnums;
129 return solverInfoEnums;
130}
std::map< std::string, EnumMap > EnumMapList

Referenced by GetSolverInfoAsEnum(), GetValueAsEnum(), RegisterEnumValue(), and VerifySolverInfo().

◆ GetTag()

const std::string & Nektar::LibUtilities::SessionReader::GetTag ( const std::string &  pName) const

Returns the value of a specified tag.

Definition at line 1447 of file BasicUtils/SessionReader.cpp.

1448{
1449 std::string vName = boost::to_upper_copy(pName);
1450 auto vTagIterator = m_tags.find(vName);
1451 ASSERTL0(vTagIterator != m_tags.end(), "Requested tag does not exist.");
1452 return vTagIterator->second;
1453}

References ASSERTL0, and m_tags.

◆ GetTimeIntScheme()

const TimeIntScheme & Nektar::LibUtilities::SessionReader::GetTimeIntScheme ( ) const

Returns the time integration scheme structure m_timeIntScheme from the session file.

Definition at line 1135 of file BasicUtils/SessionReader.cpp.

1136{
1137 return m_timeIntScheme;
1138}

References m_timeIntScheme.

◆ GetTimeLevel()

size_t Nektar::LibUtilities::SessionReader::GetTimeLevel ( void  ) const
inline

Get time level (Parallel-in-Time)

Definition at line 390 of file SessionReader.h.

391 {
392 return m_timeLevel;
393 }

References m_timeLevel.

◆ GetUpdateOptFile()

bool Nektar::LibUtilities::SessionReader::GetUpdateOptFile ( ) const
inline

Get bool to update optimisation file.

Definition at line 378 of file SessionReader.h.

379 {
380 return m_updateOptFile;
381 }
bool m_updateOptFile
Update optimisation file.

References m_updateOptFile.

◆ GetValueAsEnum()

template<typename T >
const T Nektar::LibUtilities::SessionReader::GetValueAsEnum ( const std::string &  pName,
const std::string &  vValue 
) const
inline

Returns the value of the specified property and value as enum.

Definition at line 536 of file SessionReader.h.

538{
539 std::string vName = boost::to_upper_copy(pName);
540
541 auto x = GetSolverInfoEnums().find(vName);
542 ASSERTL0(x != GetSolverInfoEnums().end(),
543 "Enum for property '" + pName + "' not found.");
544
545 auto y = x->second.find(pValue);
546 ASSERTL0(y != x->second.end(),
547 "Value of property '" + pValue + "' is invalid.");
548 return T(y->second);
549}

References ASSERTL0, and GetSolverInfoEnums().

◆ GetVariable()

const std::string & Nektar::LibUtilities::SessionReader::GetVariable ( const unsigned int &  idx) const

Returns the name of the variable specified by the given index.

Definition at line 1180 of file BasicUtils/SessionReader.cpp.

1181{
1182 ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1183 return m_variables[idx];
1184}

References ASSERTL0, and m_variables.

Referenced by export_SessionReader().

◆ GetVariables()

std::vector< std::string > Nektar::LibUtilities::SessionReader::GetVariables ( ) const

Returns the names of all variables.

Definition at line 1198 of file BasicUtils/SessionReader.cpp.

1199{
1200 return m_variables;
1201}

References m_variables.

Referenced by export_SessionReader().

◆ GetXMLElementTimeLevel()

void Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel ( TiXmlElement *&  element,
const size_t  timeLevel,
const bool  enableCheck = true 
)
static

Get XML elment time level (Parallel-in-Time)

Definition at line 1474 of file BasicUtils/SessionReader.cpp.

1477{
1478 if (Element && Element->FirstChildElement("TIMELEVEL"))
1479 {
1480 Element = Element->FirstChildElement("TIMELEVEL");
1481 std::string timeLevelStr;
1482 while (Element)
1483 {
1484 std::stringstream tagcontent;
1485 tagcontent << *Element;
1486 ASSERTL0(Element->Attribute("VALUE"),
1487 "Missing LEVEL attribute in solver info "
1488 "XML element: \n\t'" +
1489 tagcontent.str() + "'");
1490 timeLevelStr = Element->Attribute("VALUE");
1491 ASSERTL0(!timeLevelStr.empty(),
1492 "LEVEL attribute must be non-empty in XML "
1493 "element: \n\t'" +
1494 tagcontent.str() + "'");
1495 if (stoi(timeLevelStr) == timeLevel)
1496 {
1497 break;
1498 }
1499 Element = Element->NextSiblingElement("TIMELEVEL");
1500 }
1501 if (enableCheck)
1502 {
1503 ASSERTL0(stoi(timeLevelStr) == timeLevel,
1504 "TIMELEVEL value " + std::to_string(timeLevel) +
1505 " not found in solver info "
1506 "XML element: \n\t'");
1507 }
1508 }
1509}

References ASSERTL0.

Referenced by Nektar::Collections::CollectionOptimisation::CollectionOptimisation(), Nektar::GlobalMapping::Mapping::Load(), Nektar::SolverUtils::Forcing::Load(), Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryRegions(), Nektar::SpatialDomains::MeshPartition::ReadConditions(), Nektar::SpatialDomains::MeshGraph::ReadExpansionInfo(), Nektar::SpatialDomains::MeshPartition::ReadExpansions(), ReadFilters(), ReadFunctions(), ReadGlobalSysSolnInfo(), ReadParameters(), ReadSolverInfo(), ReadTimeIntScheme(), ReadVariables(), Nektar::GlobalMapping::Mapping::ReplaceField(), Nektar::Collections::CollectionOptimisation::UpdateOptFile(), Nektar::SpatialDomains::MeshGraphIOHDF5::v_PartitionMesh(), and Nektar::SpatialDomains::MeshGraphIOXml::WriteXMLGeometry().

◆ InitSession()

void Nektar::LibUtilities::SessionReader::InitSession ( const std::vector< std::string > &  filenames = std::vector<std::string>())

Performs the main initialisation of the object. The XML file provided on the command-line is loaded and any mesh partitioning is done. The resulting process-specific XML file (containing the process's geometry partition) is then reloaded and parsed.

Definition at line 274 of file BasicUtils/SessionReader.cpp.

275{
276 // Re-load filenames for session if required.
277 if (filenames.size() > 0)
278 {
279 m_filenames = filenames;
280 }
281
282 // check specified opt file
283 std::string optfile;
284 int exists;
285
286 if (DefinesCmdLineArgument("use-opt-file"))
287 {
288 optfile =
289 m_cmdLineOptions.find("use-opt-file")->second.as<std::string>();
290 exists = fs::exists(optfile.c_str());
291 ASSERTL0(exists, "A valid .opt file was not specified "
292 "with the --use-opt-file command line option");
293
294 m_filenames.push_back(optfile);
295
296 // put opt file at beginning
297 std::rotate(m_filenames.rbegin(), m_filenames.rbegin() + 1,
298 m_filenames.rend());
299 }
300 else // check for write-opt-file
301 {
302 // check for opt file
303 optfile = m_sessionName.substr(0, m_sessionName.find("_xml/")) + ".opt";
304 exists = fs::exists(optfile.c_str());
305
306 // For Paralell-in-Time
307 if (exists && m_comm->IsParallelInTime())
308 {
309 TiXmlDocument doc;
310 doc.LoadFile(optfile);
311 TiXmlElement *xmlTag = doc.FirstChildElement("NEKTAR")
312 ->FirstChildElement("COLLECTIONS")
313 ->FirstChildElement("TIMELEVEL");
314 if (xmlTag)
315 {
316 // if there is a TIMELEVEL tag, then turn existing flag to false
317 // and check if the required time level is specified in the
318 // optfile.
319 exists = false;
320 while (xmlTag)
321 {
322 if (m_timeLevel == stoi(xmlTag->Attribute("VALUE")))
323 {
324 exists = true;
325 break;
326 }
327 xmlTag = xmlTag->NextSiblingElement();
328 }
329 }
330 }
331
332 if (exists)
333 {
334 m_filenames.push_back(optfile);
335 // rotate order so opt file can be overwritten by
336 // direct choice in xml file
337 std::rotate(m_filenames.rbegin(), m_filenames.rbegin() + 1,
338 m_filenames.rend());
339 }
340 else
341 {
342 m_updateOptFile = true;
343 }
344 }
345
346 // Merge document if required.
347 if (m_xmlDoc)
348 {
349 delete m_xmlDoc;
350 }
351
353
354 // Parse the XML data in #m_xmlDoc
356
357 // Override SOLVERINFO and parameters with any specified on the
358 // command line.
360
361 // Verify SOLVERINFO values
363
364 // Disable backups if NEKTAR_DISABLE_BACKUPS is set.
365 if (std::getenv("NEKTAR_DISABLE_BACKUPS") != nullptr)
366 {
367 m_backups = false;
368 }
369
370 // In verbose mode, print out parameters and solver info sections
371 if (m_verbose && m_comm)
372 {
373 if (m_comm->TreatAsRankZero() && m_parameters.size() > 0)
374 {
375 cout << "Parameters:" << endl;
376 for (auto &x : m_parameters)
377 {
378 cout << "\t" << x.first << " = " << x.second << endl;
379 }
380 cout << endl;
381 }
382
383 if (m_comm->TreatAsRankZero() && m_solverInfo.size() > 0)
384 {
385 cout << "Solver Info:" << endl;
386 for (auto &x : m_solverInfo)
387 {
388 cout << "\t" << x.first << " = " << x.second << endl;
389 }
390 cout << endl;
391 }
392 }
393}
TiXmlDocument * MergeDoc(const std::vector< std::string > &pFilenames) const
Creates an XML document from a list of input files.
void CmdLineOverride()
Enforce parameters from command line arguments.
void ParseDocument()
Loads and parses the specified file.
void VerifySolverInfo()
Check values of solver info options are valid.
bool DefinesCmdLineArgument(const std::string &pName) const
Checks if a specified cmdline argument has been given.

References ASSERTL0, CmdLineOverride(), DefinesCmdLineArgument(), m_backups, m_cmdLineOptions, m_comm, m_filenames, m_parameters, m_sessionName, m_solverInfo, m_timeLevel, m_updateOptFile, m_verbose, m_xmlDoc, MergeDoc(), ParseDocument(), and VerifySolverInfo().

Referenced by export_SessionReader().

◆ LoadDoc()

void Nektar::LibUtilities::SessionReader::LoadDoc ( const std::string &  pFilename,
TiXmlDocument *  pDoc 
) const
private

Loads an xml file into a tinyxml doc and decompresses if needed.

Definition at line 1514 of file BasicUtils/SessionReader.cpp.

1516{
1517 if (pFilename.size() > 3 &&
1518 pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1519 {
1520 ifstream file(pFilename.c_str(), ios_base::in | ios_base::binary);
1521 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1522 stringstream ss;
1523 io::filtering_streambuf<io::input> in;
1524 in.push(io::gzip_decompressor());
1525 in.push(file);
1526 try
1527 {
1528 io::copy(in, ss);
1529 ss >> (*pDoc);
1530 }
1531 catch (io::gzip_error &)
1532 {
1534 "Error: File '" + pFilename + "' is corrupt.");
1535 }
1536 }
1537 else if (pFilename.size() > 4 &&
1538 pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1539 {
1540 fs::path pdirname(pFilename);
1541 boost::format pad("P%1$07d.xml");
1542 pad % m_comm->GetSpaceComm()->GetRank();
1543 fs::path pRankFilename(pad.str());
1544 fs::path fullpath = pdirname / pRankFilename;
1545
1546 ifstream file(PortablePath(fullpath).c_str());
1547 ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1548 file >> (*pDoc);
1549 }
1550 else
1551 {
1552 ifstream file(pFilename.c_str());
1553 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1554 file >> (*pDoc);
1555 }
1556}
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.

References ASSERTL0, Nektar::ErrorUtil::efatal, m_comm, NEKERROR, and Nektar::LibUtilities::PortablePath().

Referenced by MergeDoc().

◆ LoadParameter() [1/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
int &  var 
) const

Load an integer parameter.

Definition at line 823 of file BasicUtils/SessionReader.cpp.

824{
825 std::string vName = boost::to_upper_copy(pName);
826 auto paramIter = m_parameters.find(vName);
827 ASSERTL0(paramIter != m_parameters.end(),
828 "Required parameter '" + pName + "' not specified in session.");
829 NekDouble param = round(paramIter->second);
830 pVar = checked_cast<int>(param);
831}

References ASSERTL0, and m_parameters.

◆ LoadParameter() [2/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
int &  var,
const int &  def 
) const

Check for and load an integer parameter.

Definition at line 836 of file BasicUtils/SessionReader.cpp.

838{
839 std::string vName = boost::to_upper_copy(pName);
840 auto paramIter = m_parameters.find(vName);
841 if (paramIter != m_parameters.end())
842 {
843 NekDouble param = round(paramIter->second);
844 pVar = checked_cast<int>(param);
845 }
846 else
847 {
848 pVar = pDefault;
849 }
850}

References m_parameters.

◆ LoadParameter() [3/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
NekDouble var 
) const

Load a double precision parameter.

Definition at line 920 of file BasicUtils/SessionReader.cpp.

922{
923 std::string vName = boost::to_upper_copy(pName);
924 auto paramIter = m_parameters.find(vName);
925 ASSERTL0(paramIter != m_parameters.end(),
926 "Required parameter '" + pName + "' not specified in session.");
927 pVar = paramIter->second;
928}

References ASSERTL0, and m_parameters.

◆ LoadParameter() [4/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
NekDouble var,
const NekDouble def 
) const

Check for and load a double-precision parameter.

Definition at line 933 of file BasicUtils/SessionReader.cpp.

935{
936 std::string vName = boost::to_upper_copy(pName);
937 auto paramIter = m_parameters.find(vName);
938 if (paramIter != m_parameters.end())
939 {
940 pVar = paramIter->second;
941 }
942 else
943 {
944 pVar = pDefault;
945 }
946}

References m_parameters.

◆ LoadParameter() [5/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
size_t &  var 
) const

Load an size_t parameter.

Definition at line 888 of file BasicUtils/SessionReader.cpp.

889{
890 std::string vName = boost::to_upper_copy(pName);
891 auto paramIter = m_parameters.find(vName);
892 ASSERTL0(paramIter != m_parameters.end(),
893 "Required parameter '" + pName + "' not specified in session.");
894 NekDouble param = round(paramIter->second);
895 pVar = checked_cast<int>(param);
896}

References ASSERTL0, and m_parameters.

◆ LoadParameter() [6/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
size_t &  var,
const size_t &  def 
) const

Check for and load an size_t parameter.

Definition at line 901 of file BasicUtils/SessionReader.cpp.

903{
904 std::string vName = boost::to_upper_copy(pName);
905 auto paramIter = m_parameters.find(vName);
906 if (paramIter != m_parameters.end())
907 {
908 NekDouble param = round(paramIter->second);
909 pVar = checked_cast<int>(param);
910 }
911 else
912 {
913 pVar = pDefault;
914 }
915}

References m_parameters.

◆ LoadParameter() [7/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
unsigned int &  var 
) const

Load an unsigned integer parameter.

Definition at line 855 of file BasicUtils/SessionReader.cpp.

857{
858 std::string vName = boost::to_upper_copy(pName);
859 auto paramIter = m_parameters.find(vName);
860 ASSERTL0(paramIter != m_parameters.end(),
861 "Required parameter '" + pName + "' not specified in session.");
862 NekDouble param = round(paramIter->second);
863 pVar = checked_cast<unsigned int>(param);
864}

References ASSERTL0, and m_parameters.

◆ LoadParameter() [8/8]

void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
unsigned int &  var,
const unsigned int &  def 
) const

Check for and load an unsigned integer parameter.

Definition at line 869 of file BasicUtils/SessionReader.cpp.

871{
872 std::string vName = boost::to_upper_copy(pName);
873 auto paramIter = m_parameters.find(vName);
874 if (paramIter != m_parameters.end())
875 {
876 NekDouble param = round(paramIter->second);
877 pVar = checked_cast<unsigned int>(param);
878 }
879 else
880 {
881 pVar = pDefault;
882 }
883}

References m_parameters.

◆ LoadSolverInfo()

void Nektar::LibUtilities::SessionReader::LoadSolverInfo ( const std::string &  name,
std::string &  var,
const std::string &  def = "" 
) const

Check for and load a solver info property.

Definition at line 1027 of file BasicUtils/SessionReader.cpp.

1029{
1030 std::string vName = boost::to_upper_copy(pName);
1031 auto infoIter = m_solverInfo.find(vName);
1032 if (infoIter != m_solverInfo.end())
1033 {
1034 pVar = infoIter->second;
1035 }
1036 else
1037 {
1038 pVar = pDefault;
1039 }
1040}

References m_solverInfo.

◆ MatchSolverInfo() [1/2]

bool Nektar::LibUtilities::SessionReader::MatchSolverInfo ( const std::string &  name,
const std::string &  trueval 
) const

Check if the value of a solver info property matches.

Definition at line 1064 of file BasicUtils/SessionReader.cpp.

1066{
1067 if (DefinesSolverInfo(pName))
1068 {
1069 std::string vName = boost::to_upper_copy(pName);
1070 auto iter = m_solverInfo.find(vName);
1071 if (iter != m_solverInfo.end())
1072 {
1073 return boost::iequals(iter->second, pTrueVal);
1074 }
1075 }
1076 return false;
1077}

References DefinesSolverInfo(), and m_solverInfo.

◆ MatchSolverInfo() [2/2]

void Nektar::LibUtilities::SessionReader::MatchSolverInfo ( const std::string &  name,
const std::string &  trueval,
bool &  var,
const bool &  def = false 
) const

Check if the value of a solver info property matches.

Definition at line 1045 of file BasicUtils/SessionReader.cpp.

1048{
1049 std::string vName = boost::to_upper_copy(pName);
1050 auto infoIter = m_solverInfo.find(vName);
1051 if (infoIter != m_solverInfo.end())
1052 {
1053 pVar = boost::iequals(infoIter->second, pTrueVal);
1054 }
1055 else
1056 {
1057 pVar = pDefault;
1058 }
1059}

References m_solverInfo.

◆ MatchSolverInfoAsEnum()

template<typename T >
bool Nektar::LibUtilities::SessionReader::MatchSolverInfoAsEnum ( const std::string &  name,
const T &  trueval 
) const
inline

Check if the value of a solver info property matches.

Definition at line 503 of file SessionReader.h.

505{
506 return (GetSolverInfoAsEnum<T>(name) == trueval);
507}

◆ MergeDoc()

TiXmlDocument * Nektar::LibUtilities::SessionReader::MergeDoc ( const std::vector< std::string > &  pFilenames) const
private

Creates an XML document from a list of input files.

Definition at line 1561 of file BasicUtils/SessionReader.cpp.

1563{
1564 ASSERTL0(pFilenames.size() > 0, "No filenames for merging.");
1565
1566 // Read the first document
1567 TiXmlDocument *vMainDoc = new TiXmlDocument;
1568 LoadDoc(pFilenames[0], vMainDoc);
1569
1570 TiXmlHandle vMainHandle(vMainDoc);
1571 TiXmlElement *vMainNektar =
1572 GetChildElementOrThrow(pFilenames[0], "NEKTAR", vMainHandle);
1573
1574 // Read all subsequent XML documents.
1575 // For each element within the NEKTAR tag, use it to replace the
1576 // version already present in the loaded XML data.
1577 for (int i = 1; i < pFilenames.size(); ++i)
1578 {
1579 if ((pFilenames[i].compare(pFilenames[i].size() - 3, 3, "xml") == 0) ||
1580 (pFilenames[i].compare(pFilenames[i].size() - 6, 6, "xml.gz") ==
1581 0) ||
1582 (pFilenames[i].compare(pFilenames[i].size() - 3, 3, "opt") == 0))
1583 {
1584 TiXmlDocument *vTempDoc = new TiXmlDocument;
1585 LoadDoc(pFilenames[i], vTempDoc);
1586
1587 TiXmlHandle docHandle(vTempDoc);
1588 TiXmlElement *vTempNektar =
1589 GetChildElementOrThrow(pFilenames[i], "NEKTAR", docHandle);
1590 TiXmlElement *p = vTempNektar->FirstChildElement();
1591
1592 while (p)
1593 {
1594 TiXmlElement *vMainEntry =
1595 vMainNektar->FirstChildElement(p->Value());
1596
1597 // First check if the new item is in fact blank
1598 // replace if it is a COLLECTIONS section however.
1599 if (!p->FirstChild() && vMainEntry &&
1600 !boost::iequals(p->Value(), "COLLECTIONS"))
1601 {
1602 std::string warningmsg =
1603 "File " + pFilenames[i] + " contains " +
1604 "an empty XML element " + std::string(p->Value()) +
1605 " which will be ignored.";
1606 NEKERROR(ErrorUtil::ewarning, warningmsg.c_str());
1607 }
1608 else
1609 {
1610 if (vMainEntry)
1611 {
1612 vMainNektar->RemoveChild(vMainEntry);
1613 }
1614 TiXmlElement *q = new TiXmlElement(*p);
1615 vMainNektar->LinkEndChild(q);
1616 }
1617 p = p->NextSiblingElement();
1618 }
1619 delete vTempDoc;
1620 }
1621 }
1622 return vMainDoc;
1623}
void LoadDoc(const std::string &pFilename, TiXmlDocument *pDoc) const
Loads an xml file into a tinyxml doc and decompresses if needed.
TiXmlElement * GetChildElementOrThrow(const std::string &filename, std::string elementName, const TiXmlHandle &docHandle)
std::vector< double > q(NPUPPER *NPUPPER)

References ASSERTL0, Nektar::ErrorUtil::ewarning, Nektar::LibUtilities::GetChildElementOrThrow(), LoadDoc(), and NEKERROR.

Referenced by InitSession().

◆ ParseCommandLineArguments()

std::vector< std::string > Nektar::LibUtilities::SessionReader::ParseCommandLineArguments ( int  argc,
char *  argv[] 
)
private

Parse the program arguments and fill m_cmdLineOptions.

Parses the command-line arguments for known options and filenames.

Definition at line 439 of file BasicUtils/SessionReader.cpp.

441{
442 // List the publically visible options (listed using --help).
443 po::options_description desc("Allowed options");
444 po::options_description dep("Deprecated options");
445
446 // clang-format off
447 desc.add_options()
448 ("force-output,f","disables backups files and forces output to be "
449 "written without any checks")
450 ("help,h", "print this help message")
451 ("solverinfo,I", po::value<vector<std::string>>(),
452 "override a SOLVERINFO property")
453 ("parameter,P", po::value<vector<std::string>>(),
454 "override a parameter")
455 ("verbose,v", "be verbose")
456 ("version,V", "print version information")
457 ("no-exp-opt", "Do not use expansion optimisation for collections")
458 ("npx", po::value<int>(),
459 "number of procs in X-dir")
460 ("npy", po::value<int>(),
461 "number of procs in Y-dir")
462 ("npz", po::value<int>(),
463 "number of procs in Z-dir")
464 ("nsz", po::value<int>(),
465 "number of slices in Z-dir")
466 ("npt", po::value<int>(),
467 "number of procs in T-dir (parareal)")
468 ("part-only", po::value<int>(),
469 "only partition mesh into N partitions.")
470 ("part-only-overlapping", po::value<int>(),
471 "only partition mesh into N overlapping partitions.")
472 ("part-info", "output partition information")
473 ("write-opt-file","write an optimisation file")
474 ("use-opt-file", po::value<std::string>(),
475 "use an optimisation file");
476 // clang-format on
477
478#ifdef NEKTAR_USE_CWIPI
479 desc.add_options()("cwipi", po::value<std::string>(), "set CWIPI name");
480#endif
481
482 for (auto &cmdIt : GetCmdLineArgMap())
483 {
484 std::string names = cmdIt.first;
485 if (cmdIt.second.shortName != "")
486 {
487 names += "," + cmdIt.second.shortName;
488 }
489 if (cmdIt.second.isFlag)
490 {
491 desc.add_options()(names.c_str(), cmdIt.second.description.c_str());
492 }
493 else
494 {
495 desc.add_options()(names.c_str(), po::value<std::string>(),
496 cmdIt.second.description.c_str());
497 }
498 }
499
500 // Deprecated options: introduced in 5.4.0 to homogenise command-line
501 // options to use '-' instead of camelCase or no spaces.
502 std::map<std::string, std::string> deprecated = {
503 {"forceoutput", "force-output"},
504 {"writeoptfile", "write-opt-file"},
505 {"useoptfile", "use-opt-file"}};
506
507 for (auto &d : deprecated)
508 {
509 std::string description = "Deprecated: use --" + d.second;
510 dep.add_options()(d.first.c_str(), description.c_str());
511 }
512
513 // List hidden options (e.g. session file arguments are not actually
514 // specified using the input-file option by the user).
515 po::options_description hidden("Hidden options");
516
517 hidden.add_options()("input-file", po::value<vector<string>>(),
518 "input filename");
519
520 // Combine all options for the parser
521 po::options_description all("All options");
522 all.add(desc).add(dep).add(hidden);
523
524 // Session file is a positional option
525 po::positional_options_description p;
526 p.add("input-file", -1);
527
528 // Parse the command-line options
529 po::parsed_options parsed = po::command_line_parser(argc, argv)
530 .options(all)
531 .positional(p)
532 .allow_unregistered()
533 .run();
534
535 // Extract known options to map and update
536 po::store(parsed, m_cmdLineOptions);
537 po::notify(m_cmdLineOptions);
538
539 // Help message
540 if (m_cmdLineOptions.count("help"))
541 {
542 cout << desc;
543 exit(0);
544 }
545
546 // Version information
547 if (m_cmdLineOptions.count("version"))
548 {
549 cout << "Nektar++ version " << NEKTAR_VERSION;
550
551 if (NekConstants::kGitSha1 != "GITDIR-NOTFOUND")
552 {
554 string branch(NekConstants::kGitBranch);
555 boost::replace_all(branch, "refs/heads/", "");
556
557 cout << " (git changeset " << sha1.substr(0, 8) << ", ";
558
559 if (branch == "")
560 {
561 cout << "detached head";
562 }
563 else
564 {
565 cout << "head " << branch;
566 }
567
568 cout << ")";
569 }
570
571 cout << endl;
572 exit(0);
573 }
574
575 // Deal with deprecated options.
576 for (auto &d : deprecated)
577 {
578 if (m_cmdLineOptions.count(d.first))
579 {
580 std::cerr << "Warning: --" << d.first << " deprecated: use --"
581 << d.second << std::endl;
582 m_cmdLineOptions.emplace(
583 d.second, po::variable_value(m_cmdLineOptions[d.first]));
584 }
585 }
586
587 // Enable verbose mode
588 if (m_cmdLineOptions.count("verbose"))
589 {
590 m_verbose = true;
591 }
592 else
593 {
594 m_verbose = false;
595 }
596
597 // Disable backups
598 if (m_cmdLineOptions.count("force-output"))
599 {
600 m_backups = false;
601 }
602 else
603 {
604 m_backups = true;
605 }
606
607 // Enable update optimisation file
608 if (m_cmdLineOptions.count("write-opt-file"))
609 {
610 m_updateOptFile = true;
611 }
612 else
613 {
614 m_updateOptFile = false;
615 }
616
617 // Print a warning for unknown options
618 for (auto &x : parsed.options)
619 {
620 if (x.unregistered)
621 {
622 cout << "Warning: Unknown option: " << x.string_key << endl;
623 }
624 }
625
626 // Return the vector of filename(s) given as positional options
627 if (m_cmdLineOptions.count("input-file"))
628 {
629 return m_cmdLineOptions["input-file"].as<std::vector<std::string>>();
630 }
631 else
632 {
633 return std::vector<std::string>();
634 }
635}
#define NEKTAR_VERSION
static CmdLineArgMap & GetCmdLineArgMap()
CmdLine argument map.
const std::string kGitBranch
Definition GitRevision.h:46
const std::string kGitSha1
Definition GitRevision.h:45
std::vector< double > d(NPUPPER *NPUPPER)
Definition sha1.cpp:72

References GetCmdLineArgMap(), Nektar::NekConstants::kGitBranch, Nektar::NekConstants::kGitSha1, m_backups, m_cmdLineOptions, m_updateOptFile, m_verbose, and NEKTAR_VERSION.

Referenced by SessionReader(), and SessionReader().

◆ ParseDocument()

void Nektar::LibUtilities::SessionReader::ParseDocument ( )
private

Loads and parses the specified file.

Definition at line 1628 of file BasicUtils/SessionReader.cpp.

1629{
1630 // Check we actually have a document loaded.
1631 ASSERTL0(m_xmlDoc, "No XML document loaded.");
1632
1633 TiXmlHandle docHandle(m_xmlDoc);
1634 TiXmlElement *e;
1635
1636 // Look for all data in CONDITIONS block.
1637 e = docHandle.FirstChildElement("NEKTAR")
1638 .FirstChildElement("CONDITIONS")
1639 .Element();
1640
1641 // Read the various sections of the CONDITIONS block
1642 ReadParameters(e);
1643 ReadSolverInfo(e);
1646 ReadVariables(e);
1647 ReadFunctions(e);
1648
1649 // Look for all data in FILTERS block.
1650 e = docHandle.FirstChildElement("NEKTAR")
1651 .FirstChildElement("FILTERS")
1652 .Element();
1653
1654 // Read the various sections of the FILTERS block
1655 ReadFilters(e);
1656}
void ReadSolverInfo(TiXmlElement *conditions)
Reads the SOLVERINFO section of the XML document.
void ReadVariables(TiXmlElement *conditions)
Reads the VARIABLES section of the XML document.
void ReadTimeIntScheme(TiXmlElement *conditions)
Reads the TIMEINTEGRATIONSCHEME section of the XML document.
void ReadFilters(TiXmlElement *filters)
Reads the FILTERS section of the XML document.
void ReadGlobalSysSolnInfo(TiXmlElement *conditions)
Reads the GLOBALSYSSOLNINFO section of the XML document.
void ReadFunctions(TiXmlElement *conditions)
Reads the FUNCTIONS section of the XML document.
void ReadParameters(TiXmlElement *conditions)
Reads the PARAMETERS section of the XML document.

References ASSERTL0, m_xmlDoc, ReadFilters(), ReadFunctions(), ReadGlobalSysSolnInfo(), ReadParameters(), ReadSolverInfo(), ReadTimeIntScheme(), and ReadVariables().

Referenced by InitSession().

◆ ParseEquals()

void Nektar::LibUtilities::SessionReader::ParseEquals ( const std::string &  line,
std::string &  lhs,
std::string &  rhs 
)
private

Parse a string in the form lhs = rhs.

Pull out lhs and rhs and eliminate any spaces.

Definition at line 2511 of file BasicUtils/SessionReader.cpp.

2513{
2514 /// Pull out lhs and rhs and eliminate any spaces.
2515 size_t beg = line.find_first_not_of(" ");
2516 size_t end = line.find_first_of("=");
2517 // Check for no parameter name
2518 if (beg == end)
2519 {
2520 throw 1;
2521 }
2522 // Check for no parameter value
2523 if (end != line.find_last_of("="))
2524 {
2525 throw 1;
2526 }
2527 // Check for no equals sign
2528 if (end == std::string::npos)
2529 {
2530 throw 1;
2531 }
2532
2533 lhs = line.substr(line.find_first_not_of(" "), end - beg);
2534 lhs = lhs.substr(0, lhs.find_last_not_of(" ") + 1);
2535 rhs = line.substr(line.find_last_of("=") + 1);
2536 rhs = rhs.substr(rhs.find_first_not_of(" "));
2537 rhs = rhs.substr(0, rhs.find_last_not_of(" ") + 1);
2538}

Referenced by CmdLineOverride(), and ReadParameters().

◆ ParseSessionName()

std::string Nektar::LibUtilities::SessionReader::ParseSessionName ( std::vector< std::string > &  filenames)
private

Parse the session name.

Definition at line 640 of file BasicUtils/SessionReader.cpp.

641{
642 ASSERTL0(!filenames.empty(), "At least one filename expected.");
643
644 std::string retval = "";
645
646 // First input file defines the session name
647 std::string fname = filenames[0];
648
649 // If loading a pre-partitioned mesh, remove _xml extension
650 if (fname.size() > 4 && fname.substr(fname.size() - 4, 4) == "_xml")
651 {
652 retval = fname.substr(0, fname.find_last_of("_"));
653 }
654 // otherwise remove the .xml extension
655 else if (fname.size() > 4 && fname.substr(fname.size() - 4, 4) == ".xml")
656 {
657 retval = fname.substr(0, fname.find_last_of("."));
658 }
659 // If compressed .xml.gz, remove both extensions
660 else if (fname.size() > 7 && fname.substr(fname.size() - 7, 7) == ".xml.gz")
661 {
662 retval = fname.substr(0, fname.find_last_of("."));
663 retval = retval.substr(0, retval.find_last_of("."));
664 }
665
666 return retval;
667}

References ASSERTL0.

Referenced by SessionReader(), and SessionReader().

◆ PartitionComm()

void Nektar::LibUtilities::SessionReader::PartitionComm ( )
private

Partitions the comm object based on session parameters.

Splits the processes into a cartesian grid and creates communicators for each row and column of the grid.

Definition at line 1688 of file BasicUtils/SessionReader.cpp.

1689{
1690 if (m_comm->GetSize() > 1)
1691 {
1692 int nProcZ = 1;
1693 int nProcY = 1;
1694 int nProcX = 1;
1695 int nStripZ = 1;
1696 int nTime = 0;
1697 if (DefinesCmdLineArgument("npx"))
1698 {
1699 nProcX = GetCmdLineArgument<int>("npx");
1700 }
1701 if (DefinesCmdLineArgument("npy"))
1702 {
1703 nProcY = GetCmdLineArgument<int>("npy");
1704 }
1705 if (DefinesCmdLineArgument("npz"))
1706 {
1707 nProcZ = GetCmdLineArgument<int>("npz");
1708 }
1709 if (DefinesCmdLineArgument("nsz"))
1710 {
1711 nStripZ = GetCmdLineArgument<int>("nsz");
1712 }
1713 if (DefinesCmdLineArgument("npt"))
1714 {
1715 nTime = GetCmdLineArgument<int>("npt");
1716 }
1717
1718 if (nTime)
1719 {
1720 ASSERTL0(m_comm->GetSize() % nTime == 0,
1721 "Cannot exactly partition time using npt value.");
1722 ASSERTL0((m_comm->GetSize() / nTime) % (nProcZ * nProcY * nProcX) ==
1723 0,
1724 "Cannot exactly partition using PROC_Z value.");
1725 ASSERTL0(nProcZ % nProcY == 0,
1726 "Cannot exactly partition using PROC_Y value.");
1727 ASSERTL0(nProcY % nProcX == 0,
1728 "Cannot exactly partition using PROC_X value.");
1729 }
1730 else
1731 {
1732 ASSERTL0(m_comm->GetSize() % (nProcZ * nProcY * nProcX) == 0,
1733 "Cannot exactly partition using PROC_Z value.");
1734 ASSERTL0(nProcZ % nProcY == 0,
1735 "Cannot exactly partition using PROC_Y value.");
1736 ASSERTL0(nProcY % nProcX == 0,
1737 "Cannot exactly partition using PROC_X value.");
1738 }
1739
1740 // Number of processes associated with the spectral method
1741 int nProcSm = nProcZ * nProcY * nProcX;
1742
1743 // Number of processes associated with the spectral element
1744 // method.
1745 int nProcSem = nTime ? m_comm->GetSize() / nTime / nProcSm
1746 : m_comm->GetSize() / nProcSm;
1747
1748 m_comm->SplitComm(nProcSm, nProcSem, nTime);
1749 m_comm->GetColumnComm()->SplitComm(nProcZ / nStripZ, nStripZ);
1750 m_comm->GetColumnComm()->GetColumnComm()->SplitComm((nProcY * nProcX),
1751 nProcZ / nStripZ);
1752 m_comm->GetColumnComm()->GetColumnComm()->GetColumnComm()->SplitComm(
1753 nProcX, nProcY);
1754 }
1755}

References ASSERTL0, DefinesCmdLineArgument(), and m_comm.

Referenced by SessionReader(), and SessionReader().

◆ ReadFilters()

void Nektar::LibUtilities::SessionReader::ReadFilters ( TiXmlElement *  filters)
private

Reads the FILTERS section of the XML document.

Definition at line 2452 of file BasicUtils/SessionReader.cpp.

2453{
2454 if (!filters)
2455 {
2456 return;
2457 }
2458
2459 m_filters.clear();
2460
2461 TiXmlElement *filter = filters->FirstChildElement("FILTER");
2462
2463 while (filter)
2464 {
2465 ASSERTL0(filter->Attribute("TYPE"),
2466 "Missing attribute 'TYPE' for filter.");
2467 std::string typeStr = filter->Attribute("TYPE");
2468
2469 int domainID = -1;
2470 if (filter->Attribute("DOMAIN"))
2471 {
2472 // domainID = int(filter->Attribute("DOMAIN"));
2473 int err = filter->QueryIntAttribute("DOMAIN", &domainID);
2474 ASSERTL0(err == TIXML_SUCCESS,
2475 "Unable to read attribute DOMAIN in filter.");
2476 }
2477
2478 std::map<std::string, std::string> vParams;
2479
2480 TiXmlElement *element = filter;
2481 GetXMLElementTimeLevel(element, m_timeLevel, false);
2482 TiXmlElement *param = element->FirstChildElement("PARAM");
2483 while (param)
2484 {
2485 ASSERTL0(param->Attribute("NAME"),
2486 "Missing attribute 'NAME' for parameter in filter " +
2487 typeStr + "'.");
2488 std::string nameStr = param->Attribute("NAME");
2489
2490 ASSERTL0(param->GetText(), "Empty value string for param.");
2491 std::string valueStr = param->GetText();
2492
2493 vParams[nameStr] = valueStr;
2494
2495 param = param->NextSiblingElement("PARAM");
2496 }
2497
2498 FilterDefinition filterDef;
2499 filterDef.domain = domainID;
2500 filterDef.name = typeStr;
2501 filterDef.params = vParams;
2502 m_filters.push_back(filterDef);
2503
2504 filter = filter->NextSiblingElement("FILTER");
2505 }
2506}
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)

References ASSERTL0, Nektar::LibUtilities::FilterDefinition::domain, GetXMLElementTimeLevel(), m_filters, m_timeLevel, Nektar::LibUtilities::FilterDefinition::name, and Nektar::LibUtilities::FilterDefinition::params.

Referenced by ParseDocument().

◆ ReadFunctions()

void Nektar::LibUtilities::SessionReader::ReadFunctions ( TiXmlElement *  conditions)
private

Reads the FUNCTIONS section of the XML document.

Definition at line 2227 of file BasicUtils/SessionReader.cpp.

2228{
2229 m_functions.clear();
2230
2231 if (!conditions)
2232 {
2233 return;
2234 }
2235
2236 // Scan through conditions section looking for functions.
2237 TiXmlElement *function = conditions->FirstChildElement("FUNCTION");
2238
2239 while (function)
2240 {
2241 stringstream tagcontent;
2242 tagcontent << *function;
2243
2244 // Every function must have a NAME attribute
2245 ASSERTL0(function->Attribute("NAME"),
2246 "Functions must have a NAME attribute defined in XML "
2247 "element: \n\t'" +
2248 tagcontent.str() + "'");
2249 std::string functionStr = function->Attribute("NAME");
2250 ASSERTL0(!functionStr.empty(),
2251 "Functions must have a non-empty name in XML "
2252 "element: \n\t'" +
2253 tagcontent.str() + "'");
2254
2255 // Store function names in uppercase to remain case-insensitive.
2256 boost::to_upper(functionStr);
2257
2258 // Retrieve first entry (variable, or file)
2259 TiXmlElement *element = function;
2260 GetXMLElementTimeLevel(element, m_timeLevel, false);
2261 TiXmlElement *variable = element->FirstChildElement();
2262
2263 // Create new function structure with default type of none.
2264 FunctionVariableMap functionVarMap;
2265
2266 // Process all entries in the function block
2267 while (variable)
2268 {
2269 FunctionVariableDefinition funcDef;
2270 std::string conditionType = variable->Value();
2271
2272 // If no var is specified, assume wildcard
2273 std::string variableStr;
2274 if (!variable->Attribute("VAR"))
2275 {
2276 variableStr = "*";
2277 }
2278 else
2279 {
2280 variableStr = variable->Attribute("VAR");
2281 }
2282
2283 // Parse list of variables
2284 std::vector<std::string> variableList;
2285 ParseUtils::GenerateVector(variableStr, variableList);
2286
2287 // If no domain is specified, put to 0
2288 std::string domainStr;
2289 if (!variable->Attribute("DOMAIN"))
2290 {
2291 domainStr = "0";
2292 }
2293 else
2294 {
2295 domainStr = variable->Attribute("DOMAIN");
2296 }
2297
2298 // Parse list of domains
2299 std::vector<std::string> varSplit;
2300 std::vector<unsigned int> domainList;
2301 ParseUtils::GenerateSeqVector(domainStr, domainList);
2302
2303 // if no evars is specified, put "x y z t"
2304 std::string evarsStr = "x y z t";
2305 if (variable->Attribute("EVARS"))
2306 {
2307 evarsStr =
2308 evarsStr + std::string(" ") + variable->Attribute("EVARS");
2309 }
2310
2311 // Expressions are denoted by E
2312 if (conditionType == "E")
2313 {
2314 funcDef.m_type = eFunctionTypeExpression;
2315
2316 // Expression must have a VALUE.
2317 ASSERTL0(variable->Attribute("VALUE"),
2318 "Attribute VALUE expected for function '" +
2319 functionStr + "'.");
2320 std::string fcnStr = variable->Attribute("VALUE");
2321 ASSERTL0(!fcnStr.empty(),
2322 (std::string("Expression for var: ") + variableStr +
2323 std::string(" must be specified."))
2324 .c_str());
2325
2326 // set expression
2327 funcDef.m_expression =
2329 m_interpreter, fcnStr, evarsStr);
2330 }
2331
2332 // Files are denoted by F
2333 else if (conditionType == "F")
2334 {
2335 // Check if transient or not
2336 if (variable->Attribute("TIMEDEPENDENT") &&
2337 boost::lexical_cast<bool>(
2338 variable->Attribute("TIMEDEPENDENT")))
2339 {
2340 funcDef.m_type = eFunctionTypeTransientFile;
2341 }
2342 else
2343 {
2344 funcDef.m_type = eFunctionTypeFile;
2345 }
2346
2347 // File must have a FILE.
2348 ASSERTL0(variable->Attribute("FILE"),
2349 "Attribute FILE expected for function '" +
2350 functionStr + "'.");
2351 std::string filenameStr = variable->Attribute("FILE");
2352 ASSERTL0(!filenameStr.empty(),
2353 "A filename must be specified for the FILE "
2354 "attribute of function '" +
2355 functionStr + "'.");
2356
2357 std::vector<std::string> fSplit;
2358 boost::split(fSplit, filenameStr, boost::is_any_of(":"));
2359 ASSERTL0(fSplit.size() == 1 || fSplit.size() == 2,
2360 "Incorrect filename specification in function " +
2361 functionStr +
2362 "'. "
2363 "Specify variables inside file as: "
2364 "filename:var1,var2");
2365
2366 // set the filename
2367 fs::path fullpath = fSplit[0];
2368 fs::path ftype = fullpath.extension();
2369 if (fullpath.parent_path().extension() == ".pit")
2370 {
2371 string filename = fullpath.stem().string();
2372 fullpath = fullpath.parent_path();
2373 size_t start = filename.find_last_of("_") + 1;
2374 int index =
2375 atoi(filename.substr(start, filename.size()).c_str());
2376 fullpath /= filename.substr(0, start) +
2377 std::to_string(
2378 index + m_comm->GetTimeComm()->GetRank()) +
2379 ftype.string();
2380 }
2381 funcDef.m_filename = fullpath.string();
2382
2383 if (fSplit.size() == 2)
2384 {
2385 ASSERTL0(variableList[0] != "*",
2386 "Filename variable mapping not valid "
2387 "when using * as a variable inside "
2388 "function '" +
2389 functionStr + "'.");
2390
2391 boost::split(varSplit, fSplit[1], boost::is_any_of(","));
2392 ASSERTL0(varSplit.size() == variableList.size(),
2393 "Filename variables should contain the "
2394 "same number of variables defined in "
2395 "VAR in function " +
2396 functionStr + "'.");
2397 }
2398 }
2399
2400 // Nothing else supported so throw an error
2401 else
2402 {
2403 stringstream tagcontent;
2404 tagcontent << *variable;
2405
2407 "Identifier " + conditionType + " in function " +
2408 std::string(function->Attribute("NAME")) +
2409 " is not recognised in XML element: \n\t'" +
2410 tagcontent.str() + "'");
2411 }
2412
2413 // Add variables to function
2414 for (unsigned int i = 0; i < variableList.size(); ++i)
2415 {
2416 for (unsigned int j = 0; j < domainList.size(); ++j)
2417 {
2418 // Check it has not already been defined
2419 pair<std::string, int> key(variableList[i], domainList[j]);
2420 auto fcnsIter = functionVarMap.find(key);
2421 ASSERTL0(fcnsIter == functionVarMap.end(),
2422 "Error setting expression '" + variableList[i] +
2423 " in domain " + std::to_string(domainList[j]) +
2424 "' in function '" + functionStr +
2425 "'. "
2426 "Expression has already been defined.");
2427
2428 if (varSplit.size() > 0)
2429 {
2430 FunctionVariableDefinition funcDef2 = funcDef;
2431 funcDef2.m_fileVariable = varSplit[i];
2432 functionVarMap[key] = funcDef2;
2433 }
2434 else
2435 {
2436 functionVarMap[key] = funcDef;
2437 }
2438 }
2439 }
2440 variable = variable->NextSiblingElement();
2441 }
2442
2443 // Add function definition to map
2444 m_functions[functionStr] = functionVarMap;
2445 function = function->NextSiblingElement("FUNCTION");
2446 }
2447}
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
std::map< std::pair< std::string, int >, FunctionVariableDefinition > FunctionVariableMap

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eFunctionTypeExpression, Nektar::LibUtilities::eFunctionTypeFile, Nektar::LibUtilities::eFunctionTypeTransientFile, Nektar::ParseUtils::GenerateSeqVector(), Nektar::ParseUtils::GenerateVector(), GetXMLElementTimeLevel(), m_comm, Nektar::LibUtilities::FunctionVariableDefinition::m_expression, Nektar::LibUtilities::FunctionVariableDefinition::m_filename, Nektar::LibUtilities::FunctionVariableDefinition::m_fileVariable, m_functions, m_interpreter, m_timeLevel, Nektar::LibUtilities::FunctionVariableDefinition::m_type, and NEKERROR.

Referenced by ParseDocument().

◆ ReadGlobalSysSolnInfo()

void Nektar::LibUtilities::SessionReader::ReadGlobalSysSolnInfo ( TiXmlElement *  conditions)
private

Reads the GLOBALSYSSOLNINFO section of the XML document.

Definition at line 1918 of file BasicUtils/SessionReader.cpp.

1919{
1920 GetGloSysSolnList().clear();
1921
1922 if (!conditions)
1923 {
1924 return;
1925 }
1926
1927 TiXmlElement *GlobalSys =
1928 conditions->FirstChildElement("GLOBALSYSSOLNINFO");
1930
1931 if (!GlobalSys)
1932 {
1933 return;
1934 }
1935
1936 TiXmlElement *VarInfo = GlobalSys->FirstChildElement("V");
1937
1938 while (VarInfo)
1939 {
1940 std::stringstream tagcontent;
1941 tagcontent << *VarInfo;
1942
1943 ASSERTL0(VarInfo->Attribute("VAR"),
1944 "Missing VAR attribute in GobalSysSolnInfo XML "
1945 "element: \n\t'" +
1946 tagcontent.str() + "'");
1947 std::string VarList = VarInfo->Attribute("VAR");
1948 ASSERTL0(!VarList.empty(),
1949 "VAR attribute must be non-empty in XML element:\n\t'" +
1950 tagcontent.str() + "'");
1951
1952 // generate a list of variables.
1953 std::vector<std::string> varStrings;
1954 bool valid = ParseUtils::GenerateVector(VarList, varStrings);
1955
1956 ASSERTL0(valid, "Unable to process list of variable in XML "
1957 "element \n\t'" +
1958 tagcontent.str() + "'");
1959
1960 if (varStrings.size())
1961 {
1962 TiXmlElement *SysSolnInfo = VarInfo->FirstChildElement("I");
1963
1964 while (SysSolnInfo)
1965 {
1966 tagcontent.clear();
1967 tagcontent << *SysSolnInfo;
1968 // read the property name
1969 ASSERTL0(SysSolnInfo->Attribute("PROPERTY"),
1970 "Missing PROPERTY attribute in "
1971 "GlobalSysSolnInfo for variable(s) '" +
1972 VarList + "' in XML element: \n\t'" +
1973 tagcontent.str() + "'");
1974 std::string SysSolnProperty =
1975 SysSolnInfo->Attribute("PROPERTY");
1976 ASSERTL0(!SysSolnProperty.empty(),
1977 "GlobalSysSolnIno properties must have a "
1978 "non-empty name for variable(s) : '" +
1979 VarList + "' in XML element: \n\t'" +
1980 tagcontent.str() + "'");
1981
1982 // make sure that solver property is capitalised
1983 std::string SysSolnPropertyUpper =
1984 boost::to_upper_copy(SysSolnProperty);
1985
1986 // read the value
1987 ASSERTL0(SysSolnInfo->Attribute("VALUE"),
1988 "Missing VALUE attribute in GlobalSysSolnInfo "
1989 "for variable(s) '" +
1990 VarList + "' in XML element: \n\t" +
1991 tagcontent.str() + "'");
1992 std::string SysSolnValue = SysSolnInfo->Attribute("VALUE");
1993 ASSERTL0(!SysSolnValue.empty(),
1994 "GlobalSysSolnInfo properties must have a "
1995 "non-empty value for variable(s) '" +
1996 VarList + "' in XML element: \n\t'" +
1997 tagcontent.str() + "'");
1998
1999 // Store values under variable map.
2000 for (int i = 0; i < varStrings.size(); ++i)
2001 {
2002 auto x = GetGloSysSolnList().find(varStrings[i]);
2003 if (x == GetGloSysSolnList().end())
2004 {
2005 (GetGloSysSolnList()[varStrings[i]])
2006 [SysSolnPropertyUpper] = SysSolnValue;
2007 }
2008 else
2009 {
2010 x->second[SysSolnPropertyUpper] = SysSolnValue;
2011 }
2012 }
2013 SysSolnInfo = SysSolnInfo->NextSiblingElement("I");
2014 }
2015 VarInfo = VarInfo->NextSiblingElement("V");
2016 }
2017 }
2018
2019 if (m_verbose && GetGloSysSolnList().size() > 0 && m_comm)
2020 {
2021 if (m_comm->GetRank() == 0)
2022 {
2023 cout << "GlobalSysSoln Info:" << endl;
2024
2025 for (auto &x : GetGloSysSolnList())
2026 {
2027 cout << "\t Variable: " << x.first << endl;
2028
2029 for (auto &y : x.second)
2030 {
2031 cout << "\t\t " << y.first << " = " << y.second << endl;
2032 }
2033 }
2034 cout << endl;
2035 }
2036 }
2037}

References ASSERTL0, Nektar::ParseUtils::GenerateVector(), GetGloSysSolnList(), GetXMLElementTimeLevel(), m_comm, m_timeLevel, and m_verbose.

Referenced by ParseDocument().

◆ ReadParameters()

void Nektar::LibUtilities::SessionReader::ReadParameters ( TiXmlElement *  conditions)
private

Reads the PARAMETERS section of the XML document.

Definition at line 1760 of file BasicUtils/SessionReader.cpp.

1761{
1762 m_parameters.clear();
1763
1764 if (!conditions)
1765 {
1766 return;
1767 }
1768
1769 TiXmlElement *parametersElement =
1770 conditions->FirstChildElement("PARAMETERS");
1771 GetXMLElementTimeLevel(parametersElement, m_timeLevel);
1772
1773 // See if we have parameters defined. They are optional so we go on
1774 // if not.
1775 if (parametersElement)
1776 {
1777 TiXmlElement *parameter = parametersElement->FirstChildElement("P");
1778
1779 ParameterMap caseSensitiveParameters;
1780
1781 // Multiple nodes will only occur if there is a comment in
1782 // between definitions.
1783 while (parameter)
1784 {
1785 stringstream tagcontent;
1786 tagcontent << *parameter;
1787 TiXmlNode *node = parameter->FirstChild();
1788
1789 while (node && node->Type() != TiXmlNode::TINYXML_TEXT)
1790 {
1791 node = node->NextSibling();
1792 }
1793
1794 if (node)
1795 {
1796 // Format is "paramName = value"
1797 std::string line = node->ToText()->Value(), lhs, rhs;
1798
1799 try
1800 {
1801 ParseEquals(line, lhs, rhs);
1802 }
1803 catch (...)
1804 {
1806 "Syntax error in parameter expression '" + line +
1807 "' in XML element: \n\t'" + tagcontent.str() +
1808 "'");
1809 }
1810
1811 // We want the list of parameters to have their RHS
1812 // evaluated, so we use the expression evaluator to do
1813 // the dirty work.
1814 if (!lhs.empty() && !rhs.empty())
1815 {
1816 NekDouble value = 0.0;
1817 try
1818 {
1819 LibUtilities::Equation expession(m_interpreter, rhs);
1820 value = expession.Evaluate();
1821 }
1822 catch (const std::runtime_error &)
1823 {
1825 "Error evaluating parameter expression"
1826 " '" +
1827 rhs + "' in XML element: \n\t'" +
1828 tagcontent.str() + "'");
1829 }
1830 m_interpreter->SetParameter(lhs, value);
1831 caseSensitiveParameters[lhs] = value;
1832 boost::to_upper(lhs);
1833 m_parameters[lhs] = value;
1834 }
1835 }
1836 parameter = parameter->NextSiblingElement();
1837 }
1838 }
1839}
std::map< std::string, NekDouble > ParameterMap

References Nektar::ErrorUtil::efatal, Nektar::LibUtilities::Equation::Evaluate(), GetXMLElementTimeLevel(), m_interpreter, m_parameters, m_timeLevel, NEKERROR, and ParseEquals().

Referenced by ParseDocument().

◆ ReadSolverInfo()

void Nektar::LibUtilities::SessionReader::ReadSolverInfo ( TiXmlElement *  conditions)
private

Reads the SOLVERINFO section of the XML document.

Definition at line 1844 of file BasicUtils/SessionReader.cpp.

1845{
1846 m_solverInfo.clear();
1848
1849 if (!conditions)
1850 {
1851 return;
1852 }
1853
1854 TiXmlElement *solverInfoElement =
1855 conditions->FirstChildElement("SOLVERINFO");
1856 GetXMLElementTimeLevel(solverInfoElement, m_timeLevel);
1857
1858 if (solverInfoElement)
1859 {
1860 TiXmlElement *solverInfo = solverInfoElement->FirstChildElement("I");
1861
1862 while (solverInfo)
1863 {
1864 std::stringstream tagcontent;
1865 tagcontent << *solverInfo;
1866 // read the property name
1867 ASSERTL0(solverInfo->Attribute("PROPERTY"),
1868 "Missing PROPERTY attribute in solver info "
1869 "XML element: \n\t'" +
1870 tagcontent.str() + "'");
1871 std::string solverProperty = solverInfo->Attribute("PROPERTY");
1872 ASSERTL0(!solverProperty.empty(),
1873 "PROPERTY attribute must be non-empty in XML "
1874 "element: \n\t'" +
1875 tagcontent.str() + "'");
1876
1877 // make sure that solver property is capitalised
1878 std::string solverPropertyUpper =
1879 boost::to_upper_copy(solverProperty);
1880
1881 // read the value
1882 ASSERTL0(solverInfo->Attribute("VALUE"),
1883 "Missing VALUE attribute in solver info "
1884 "XML element: \n\t'" +
1885 tagcontent.str() + "'");
1886 std::string solverValue = solverInfo->Attribute("VALUE");
1887 ASSERTL0(!solverValue.empty(),
1888 "VALUE attribute must be non-empty in XML "
1889 "element: \n\t'" +
1890 tagcontent.str() + "'");
1891
1892 // Set Variable
1893 m_solverInfo[solverPropertyUpper] = solverValue;
1894 solverInfo = solverInfo->NextSiblingElement("I");
1895 }
1896 }
1897
1898 if (m_comm && m_comm->GetRowComm()->GetSize() > 1)
1899 {
1900 ASSERTL0(
1901 m_solverInfo["GLOBALSYSSOLN"] == "IterativeFull" ||
1902 m_solverInfo["GLOBALSYSSOLN"] == "IterativeStaticCond" ||
1903 m_solverInfo["GLOBALSYSSOLN"] ==
1904 "IterativeMultiLevelStaticCond" ||
1905 m_solverInfo["GLOBALSYSSOLN"] == "XxtFull" ||
1906 m_solverInfo["GLOBALSYSSOLN"] == "XxtStaticCond" ||
1907 m_solverInfo["GLOBALSYSSOLN"] == "XxtMultiLevelStaticCond" ||
1908 m_solverInfo["GLOBALSYSSOLN"] == "PETScFull" ||
1909 m_solverInfo["GLOBALSYSSOLN"] == "PETScStaticCond" ||
1910 m_solverInfo["GLOBALSYSSOLN"] == "PETScMultiLevelStaticCond",
1911 "A parallel solver must be used when run in parallel.");
1912 }
1913}

References ASSERTL0, GetSolverInfoDefaults(), GetXMLElementTimeLevel(), m_comm, m_solverInfo, and m_timeLevel.

Referenced by ParseDocument().

◆ ReadTimeIntScheme()

void Nektar::LibUtilities::SessionReader::ReadTimeIntScheme ( TiXmlElement *  conditions)
private

Reads the TIMEINTEGRATIONSCHEME section of the XML document.

Read the time-integration scheme structure, if present.

Definition at line 2042 of file BasicUtils/SessionReader.cpp.

2043{
2044 if (!conditions)
2045 {
2046 return;
2047 }
2048
2049 TiXmlElement *timeInt =
2050 conditions->FirstChildElement("TIMEINTEGRATIONSCHEME");
2052
2053 if (!timeInt)
2054 {
2055 return;
2056 }
2057
2058 TiXmlElement *method = timeInt->FirstChildElement("METHOD");
2059 TiXmlElement *variant = timeInt->FirstChildElement("VARIANT");
2060 TiXmlElement *order = timeInt->FirstChildElement("ORDER");
2061 TiXmlElement *params = timeInt->FirstChildElement("FREEPARAMETERS");
2062
2063 // Only the method and order are required.
2064 ASSERTL0(method, "Missing METHOD tag inside "
2065 "TIMEINTEGRATIONSCHEME.");
2066 ASSERTL0(order, "Missing ORDER tag inside "
2067 "TIMEINTEGRATIONSCHEME.");
2068
2069 m_timeIntScheme.method = method->GetText();
2070
2071 std::string orderStr = order->GetText();
2072
2073 // Only the method and order are required.
2074 ASSERTL0(m_timeIntScheme.method.size() > 0,
2075 "Empty text inside METHOD tag in TIMEINTEGRATIONSCHEME.");
2076 ASSERTL0(orderStr.size() > 0,
2077 "Empty text inside ORDER tag in TIMEINTEGRATIONSCHEME.");
2078 try
2079 {
2080 m_timeIntScheme.order = std::stol(orderStr);
2081 }
2082 catch (...)
2083 {
2084 NEKERROR(ErrorUtil::efatal, "In ORDER tag, unable to convert "
2085 "string '" +
2086 orderStr + "' to an unsigned integer.");
2087 }
2088
2089 if (variant)
2090 {
2091 m_timeIntScheme.variant = variant->GetText();
2092 }
2093
2094 if (params)
2095 {
2096 std::string paramsStr = params->GetText();
2097 ASSERTL0(paramsStr.size() > 0,
2098 "Empty text inside FREEPARAMETERS tag in "
2099 "TIMEINTEGRATIONSCHEME.");
2100
2101 std::vector<std::string> pSplit;
2102 boost::split(pSplit, paramsStr, boost::is_any_of(" "));
2103
2104 m_timeIntScheme.freeParams.resize(pSplit.size());
2105 for (size_t i = 0; i < pSplit.size(); ++i)
2106 {
2107 try
2108 {
2109 m_timeIntScheme.freeParams[i] = std::stod(pSplit[i]);
2110 }
2111 catch (...)
2112 {
2113 NEKERROR(ErrorUtil::efatal, "In FREEPARAMETERS tag, "
2114 "unable to convert string '" +
2115 pSplit[i] +
2116 "' "
2117 "to a floating-point value.");
2118 }
2119 }
2120 }
2121
2122 if (m_verbose && m_comm)
2123 {
2124 if (m_comm->GetRank() == 0)
2125 {
2126 cout << "Trying to use time integration scheme:" << endl;
2127 cout << "\t Method : " << m_timeIntScheme.method << endl;
2128 cout << "\t Variant: " << m_timeIntScheme.variant << endl;
2129 cout << "\t Order : " << m_timeIntScheme.order << endl;
2130
2131 if (m_timeIntScheme.freeParams.size() > 0)
2132 {
2133 cout << "\t Params :";
2134 for (auto &x : m_timeIntScheme.freeParams)
2135 {
2136 cout << " " << x;
2137 }
2138 cout << endl;
2139 }
2140 }
2141 }
2142}
std::vector< NekDouble > freeParams

References ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::LibUtilities::TimeIntScheme::freeParams, GetXMLElementTimeLevel(), m_comm, m_timeIntScheme, m_timeLevel, m_verbose, Nektar::LibUtilities::TimeIntScheme::method, NEKERROR, Nektar::LibUtilities::TimeIntScheme::order, and Nektar::LibUtilities::TimeIntScheme::variant.

Referenced by ParseDocument().

◆ ReadVariables()

void Nektar::LibUtilities::SessionReader::ReadVariables ( TiXmlElement *  conditions)
private

Reads the VARIABLES section of the XML document.

All elements are of the form: "<V ID="#"> name = value </V>", with ? being the element type.

Definition at line 2147 of file BasicUtils/SessionReader.cpp.

2148{
2149 m_variables.clear();
2150
2151 if (!conditions)
2152 {
2153 return;
2154 }
2155
2156 TiXmlElement *variablesElement = conditions->FirstChildElement("VARIABLES");
2157 GetXMLElementTimeLevel(variablesElement, m_timeLevel);
2158
2159 // See if we have parameters defined. They are optional so we go on
2160 // if not.
2161 if (variablesElement)
2162 {
2163 TiXmlElement *varElement = variablesElement->FirstChildElement("V");
2164
2165 // Sequential counter for the composite numbers.
2166 int nextVariableNumber = -1;
2167
2168 while (varElement)
2169 {
2170 stringstream tagcontent;
2171 tagcontent << *varElement;
2172
2173 /// All elements are of the form: "<V ID="#"> name = value
2174 /// </V>", with ? being the element type.
2175 nextVariableNumber++;
2176
2177 int i;
2178 int err = varElement->QueryIntAttribute("ID", &i);
2179 ASSERTL0(err == TIXML_SUCCESS,
2180 "Variables must have a unique ID number attribute "
2181 "in XML element: \n\t'" +
2182 tagcontent.str() + "'");
2183 ASSERTL0(i == nextVariableNumber,
2184 "ID numbers for variables must begin with zero and"
2185 " be sequential in XML element: \n\t'" +
2186 tagcontent.str() + "'");
2187
2188 TiXmlNode *varChild = varElement->FirstChild();
2189 // This is primarily to skip comments that may be present.
2190 // Comments appear as nodes just like elements. We are
2191 // specifically looking for text in the body of the
2192 // definition.
2193 while (varChild && varChild->Type() != TiXmlNode::TINYXML_TEXT)
2194 {
2195 varChild = varChild->NextSibling();
2196 }
2197
2198 ASSERTL0(varChild, "Unable to read variable definition body for "
2199 "variable with ID " +
2200 std::to_string(i) +
2201 " in XML element: \n\t'" + tagcontent.str() +
2202 "'");
2203 std::string variableName = varChild->ToText()->ValueStr();
2204
2205 std::istringstream variableStrm(variableName);
2206 variableStrm >> variableName;
2207
2208 ASSERTL0(std::find(m_variables.begin(), m_variables.end(),
2209 variableName) == m_variables.end(),
2210 "Variable with ID " + std::to_string(i) +
2211 " in XML element \n\t'" + tagcontent.str() +
2212 "'\nhas already been defined.");
2213
2214 m_variables.push_back(variableName);
2215
2216 varElement = varElement->NextSiblingElement("V");
2217 }
2218
2219 ASSERTL0(nextVariableNumber > -1,
2220 "Number of variables must be greater than zero.");
2221 }
2222}

References ASSERTL0, GetXMLElementTimeLevel(), m_timeLevel, and m_variables.

Referenced by ParseDocument().

◆ RegisterCmdLineArgument()

std::string Nektar::LibUtilities::SessionReader::RegisterCmdLineArgument ( const std::string &  pName,
const std::string &  pShortName,
const std::string &  pDescription 
)
inlinestatic

Registers a command-line argument with the session reader.

Definition at line 621 of file SessionReader.h.

624{
625 ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
626 CmdLineArg x;
627 x.shortName = pShortName;
628 x.description = pDescription;
629 x.isFlag = false;
630 GetCmdLineArgMap()[pName] = x;
631 return pName;
632}

References ASSERTL0, Nektar::LibUtilities::CmdLineArg::description, GetCmdLineArgMap(), Nektar::LibUtilities::CmdLineArg::isFlag, and Nektar::LibUtilities::CmdLineArg::shortName.

◆ RegisterCmdLineFlag()

std::string Nektar::LibUtilities::SessionReader::RegisterCmdLineFlag ( const std::string &  pName,
const std::string &  pShortName,
const std::string &  pDescription 
)
inlinestatic

Registers a command-line flag with the session reader.

Definition at line 637 of file SessionReader.h.

640{
641 ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
642 CmdLineArg x;
643 x.shortName = pShortName;
644 x.description = pDescription;
645 x.isFlag = true;
646 GetCmdLineArgMap()[pName] = x;
647 return pName;
648}

References ASSERTL0, Nektar::LibUtilities::CmdLineArg::description, GetCmdLineArgMap(), Nektar::LibUtilities::CmdLineArg::isFlag, and Nektar::LibUtilities::CmdLineArg::shortName.

◆ RegisterDefaultSolverInfo()

std::string Nektar::LibUtilities::SessionReader::RegisterDefaultSolverInfo ( const std::string &  pName,
const std::string &  pValue 
)
inlinestatic

Registers the default string value of a solver info property.

A default value for a given solver info property may be registered using this function. The property will take this value until it is overwritten by a value specified in the XML document, or specified as a command-line argument. Usage has the form:

std::string GlobalLinSys::def
"GlobalSysSoln","DirectMultiLevelStaticCond");
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
Parameters
pNameThe name of the property.
pValueThe default value of the property.
Returns
The default value of the property provided by #pValue.

Definition at line 610 of file SessionReader.h.

612{
613 std::string vName = boost::to_upper_copy(pName);
614 GetSolverInfoDefaults()[vName] = pValue;
615 return pValue;
616}

References GetSolverInfoDefaults().

◆ RegisterEnumValue()

std::string Nektar::LibUtilities::SessionReader::RegisterEnumValue ( std::string  pEnum,
std::string  pString,
int  pEnumValue 
)
inlinestatic

Registers an enumeration value.

A set of valid values for a given solver info property may be registered using this function. It must be called statically during the initialisation of a static variable. For example:

std::string GlobalLinSys::lookupIds[2] = {
"GlobalSysSoln",
"DirectFull",
"GlobalSysSoln",
"DirectStaticCond",
}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
Parameters
pEnumThe name of the property.
pStringA valid value for the property.
pEnumValueAn enumeration value corresponding to this value.
Returns
The value for the property provided by #pString.

Definition at line 576 of file SessionReader.h.

579{
580 std::string vEnum = boost::to_upper_copy(pEnum);
581 auto x = GetSolverInfoEnums().find(vEnum);
582
583 if (x == GetSolverInfoEnums().end())
584 {
585 GetSolverInfoEnums()[vEnum] = EnumMap();
586 x = GetSolverInfoEnums().find(vEnum);
587 }
588
589 x->second[pString] = pEnumValue;
590 return pString;
591}
std::map< std::string, int > EnumMap

References GetSolverInfoEnums().

◆ SetParameter() [1/4]

void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
int &  var 
)

Set an integer parameter.

Definition at line 951 of file BasicUtils/SessionReader.cpp.

952{
953 std::string vName = boost::to_upper_copy(pName);
954 m_parameters[vName] = pVar;
955}

References m_parameters.

Referenced by export_SessionReader().

◆ SetParameter() [2/4]

void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
NekDouble var 
)

Set a double precision parameter.

Definition at line 978 of file BasicUtils/SessionReader.cpp.

979{
980 std::string vName = boost::to_upper_copy(pName);
981 m_parameters[vName] = pVar;
982}

References m_parameters.

◆ SetParameter() [3/4]

void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
size_t &  var 
)

Set an size_t parameter.

Definition at line 969 of file BasicUtils/SessionReader.cpp.

970{
971 std::string vName = boost::to_upper_copy(pName);
972 m_parameters[vName] = pVar;
973}

References m_parameters.

◆ SetParameter() [4/4]

void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
unsigned int &  var 
)

Set an unsigned integer parameter.

Definition at line 960 of file BasicUtils/SessionReader.cpp.

961{
962 std::string vName = boost::to_upper_copy(pName);
963 m_parameters[vName] = pVar;
964}

References m_parameters.

◆ SetSolverInfo()

void Nektar::LibUtilities::SessionReader::SetSolverInfo ( const std::string &  pProperty,
const std::string &  pValue 
)

Sets the value of the specified solver info property.

Definition at line 1012 of file BasicUtils/SessionReader.cpp.

1014{
1015 std::string vProperty = boost::to_upper_copy(pProperty);
1016 auto iter = m_solverInfo.find(vProperty);
1017
1018 ASSERTL1(iter != m_solverInfo.end(),
1019 "Unable to find requested property: " + pProperty);
1020
1021 iter->second = pValue;
1022}

References ASSERTL1, and m_solverInfo.

Referenced by export_SessionReader().

◆ SetTag()

void Nektar::LibUtilities::SessionReader::SetTag ( const std::string &  pName,
const std::string &  pValue 
)

Sets a specified tag.

Definition at line 1438 of file BasicUtils/SessionReader.cpp.

1439{
1440 std::string vName = boost::to_upper_copy(pName);
1441 m_tags[vName] = pValue;
1442}

References m_tags.

◆ SetUpdateOptFile()

void Nektar::LibUtilities::SessionReader::SetUpdateOptFile ( bool  flag)
inline

Set bool to update optimisation file.

Definition at line 384 of file SessionReader.h.

385 {
386 m_updateOptFile = flag;
387 }

References m_updateOptFile.

◆ SetVariable()

void Nektar::LibUtilities::SessionReader::SetVariable ( const unsigned int &  idx,
std::string  newname 
)

Definition at line 1189 of file BasicUtils/SessionReader.cpp.

1190{
1191 ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1192 m_variables[idx] = newname;
1193}

References ASSERTL0, and m_variables.

◆ TestSharedFilesystem()

void Nektar::LibUtilities::SessionReader::TestSharedFilesystem ( )
private

Definition at line 398 of file BasicUtils/SessionReader.cpp.

399{
400 m_sharedFilesystem = false;
401
402 if (m_comm->GetSize() > 1)
403 {
404 if (m_comm->GetRank() == 0)
405 {
406 std::ofstream testfile("shared-fs-testfile");
407 testfile << "" << std::endl;
408 ASSERTL1(!testfile.fail(), "Test file creation failed");
409 testfile.close();
410 }
411 m_comm->Block();
412
413 int exists = fs::exists("shared-fs-testfile");
414 m_comm->AllReduce(exists, LibUtilities::ReduceSum);
415
416 m_sharedFilesystem = (exists == m_comm->GetSize());
417
418 if ((m_sharedFilesystem && m_comm->GetRank() == 0) ||
420 {
421 std::remove("shared-fs-testfile");
422 }
423 }
424 else
425 {
426 m_sharedFilesystem = false;
427 }
428
429 if (m_verbose && m_comm->GetRank() == 0 && m_sharedFilesystem)
430 {
431 cout << "Shared filesystem detected" << endl;
432 }
433}

References ASSERTL1, m_comm, m_sharedFilesystem, m_verbose, and Nektar::LibUtilities::ReduceSum.

Referenced by SessionReader(), and SessionReader().

◆ VerifySolverInfo()

void Nektar::LibUtilities::SessionReader::VerifySolverInfo ( )
private

Check values of solver info options are valid.

Definition at line 2609 of file BasicUtils/SessionReader.cpp.

2610{
2611 for (auto &x : m_solverInfo)
2612 {
2613 std::string solverProperty = x.first;
2614 std::string solverValue = x.second;
2615
2616 auto propIt = GetSolverInfoEnums().find(solverProperty);
2617 if (propIt != GetSolverInfoEnums().end())
2618 {
2619 auto valIt = propIt->second.find(solverValue);
2620 ASSERTL0(valIt != propIt->second.end(), "Value '" + solverValue +
2621 "' is not valid for "
2622 "property '" +
2623 solverProperty + "'");
2624 }
2625 }
2626}

References ASSERTL0, GetSolverInfoEnums(), and m_solverInfo.

Referenced by InitSession().

Friends And Related Symbol Documentation

◆ MemoryManager< SessionReader >

friend class MemoryManager< SessionReader >
friend

Support creation through MemoryManager.

Definition at line 111 of file SessionReader.h.

Member Data Documentation

◆ m_backups

bool Nektar::LibUtilities::SessionReader::m_backups = true
private

Backups.

Definition at line 436 of file SessionReader.h.

Referenced by GetBackups(), InitSession(), and ParseCommandLineArguments().

◆ m_cmdLineOptions

boost::program_options::variables_map Nektar::LibUtilities::SessionReader::m_cmdLineOptions
private

◆ m_comm

CommSharedPtr Nektar::LibUtilities::SessionReader::m_comm
private

◆ m_filenames

std::vector<std::string> Nektar::LibUtilities::SessionReader::m_filenames
private

Filenames.

Definition at line 406 of file SessionReader.h.

Referenced by GetFilenames(), InitSession(), SessionReader(), and SessionReader().

◆ m_filters

FilterMap Nektar::LibUtilities::SessionReader::m_filters
private

Filters map.

Definition at line 426 of file SessionReader.h.

Referenced by GetFilters(), and ReadFilters().

◆ m_functions

FunctionMap Nektar::LibUtilities::SessionReader::m_functions
private

◆ m_geometricInfo

GeometricInfoMap Nektar::LibUtilities::SessionReader::m_geometricInfo
private

Geometric information properties.

Definition at line 416 of file SessionReader.h.

◆ m_interpreter

InterpreterSharedPtr Nektar::LibUtilities::SessionReader::m_interpreter
private

◆ m_parameters

ParameterMap Nektar::LibUtilities::SessionReader::m_parameters
private

◆ m_sessionName

std::string Nektar::LibUtilities::SessionReader::m_sessionName
private

Session name of the loaded XML document (filename minus ext).

Definition at line 408 of file SessionReader.h.

Referenced by GetSessionName(), InitSession(), SessionReader(), and SessionReader().

◆ m_sharedFilesystem

bool Nektar::LibUtilities::SessionReader::m_sharedFilesystem
private

Running on a shared filesystem.

Definition at line 434 of file SessionReader.h.

Referenced by GetSharedFilesystem(), and TestSharedFilesystem().

◆ m_solverInfo

SolverInfoMap Nektar::LibUtilities::SessionReader::m_solverInfo
private

◆ m_tags

TagMap Nektar::LibUtilities::SessionReader::m_tags
private

Custom tags.

Definition at line 424 of file SessionReader.h.

Referenced by DefinesTag(), GetTag(), and SetTag().

◆ m_timeIntScheme

TimeIntScheme Nektar::LibUtilities::SessionReader::m_timeIntScheme
private

Time integration scheme information.

Definition at line 428 of file SessionReader.h.

Referenced by DefinesTimeIntScheme(), GetTimeIntScheme(), and ReadTimeIntScheme().

◆ m_timeLevel

size_t Nektar::LibUtilities::SessionReader::m_timeLevel = 0
private

◆ m_updateOptFile

bool Nektar::LibUtilities::SessionReader::m_updateOptFile
private

Update optimisation file.

Definition at line 438 of file SessionReader.h.

Referenced by GetUpdateOptFile(), InitSession(), ParseCommandLineArguments(), and SetUpdateOptFile().

◆ m_variables

VariableList Nektar::LibUtilities::SessionReader::m_variables
private

◆ m_verbose

bool Nektar::LibUtilities::SessionReader::m_verbose
private

◆ m_xmlDoc

TiXmlDocument* Nektar::LibUtilities::SessionReader::m_xmlDoc
private

Pointer to the loaded XML document.

Definition at line 410 of file SessionReader.h.

Referenced by DefinesElement(), GetDocument(), GetElement(), GetGeometryType(), InitSession(), ParseDocument(), SessionReader(), SessionReader(), and ~SessionReader().