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

Static Public Member Functions

static SessionReaderSharedPtr CreateInstance (int argc, char *argv[])
 Creates an instance of the SessionReader class. More...
 
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. More...
 
static std::string RegisterEnumValue (std::string pEnum, std::string pString, int pEnumValue)
 Registers an enumeration value. More...
 
static std::string RegisterDefaultSolverInfo (const std::string &pName, const std::string &pValue)
 Registers the default string value of a solver info property. More...
 
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. More...
 
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. More...
 
static void GetXMLElementTimeLevel (TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
 Get XML elment time level (Parallel-in-Time) More...
 

Private Member Functions

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

Static Private Member Functions

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

Private Attributes

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

Friends

class MemoryManager< SessionReader >
 Support creation through MemoryManager. More...
 

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 116 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)
Definition: ErrorUtil.hpp:208
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 2474 of file BasicUtils/SessionReader.cpp.

2475{
2476 // Parse solver info overrides
2477 if (m_cmdLineOptions.count("solverinfo"))
2478 {
2479 std::vector<std::string> solverInfoList =
2480 m_cmdLineOptions["solverinfo"].as<std::vector<std::string>>();
2481
2482 for (size_t i = 0; i < solverInfoList.size(); ++i)
2483 {
2484 std::string lhs, rhs;
2485
2486 try
2487 {
2488 ParseEquals(solverInfoList[i], lhs, rhs);
2489 }
2490 catch (...)
2491 {
2492 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2493 "option: " +
2494 solverInfoList[i]);
2495 }
2496
2497 std::string lhsUpper = boost::to_upper_copy(lhs);
2498 m_solverInfo[lhsUpper] = rhs;
2499 }
2500 }
2501
2502 if (m_cmdLineOptions.count("parameter"))
2503 {
2504 std::vector<std::string> parametersList =
2505 m_cmdLineOptions["parameter"].as<std::vector<std::string>>();
2506
2507 for (size_t i = 0; i < parametersList.size(); ++i)
2508 {
2509 std::string lhs, rhs;
2510
2511 try
2512 {
2513 ParseEquals(parametersList[i], lhs, rhs);
2514 }
2515 catch (...)
2516 {
2517 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2518 "option: " +
2519 parametersList[i]);
2520 }
2521
2522 std::string lhsUpper = boost::to_upper_copy(lhs);
2523
2524 try
2525 {
2526 m_parameters[lhsUpper] = std::stod(rhs);
2527 }
2528 catch (...)
2529 {
2530 NEKERROR(ErrorUtil::efatal, "Unable to convert string: " + rhs +
2531 "to double value.");
2532 }
2533 }
2534 }
2535}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
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 1619 of file BasicUtils/SessionReader.cpp.

1620{
1621 if (argc == 0)
1622 {
1623 m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1624 }
1625 else
1626 {
1627 string vCommModule("Serial");
1628 if (GetCommFactory().ModuleExists("ParallelMPI"))
1629 {
1630 vCommModule = "ParallelMPI";
1631 }
1632 if (m_cmdLineOptions.count("cwipi") &&
1633 GetCommFactory().ModuleExists("CWIPI"))
1634 {
1635 vCommModule = "CWIPI";
1636 }
1637
1638 m_comm = GetCommFactory().CreateInstance(vCommModule, argc, argv);
1639 }
1640}
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().

◆ CreateInstance() [1/2]

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

◆ 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 151 of file SessionReader.h.

154 {
157 argc, argv, pFilenames, pComm, timelevel);
158 return p;
159 }

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ DefinesCmdLineArgument()

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

Checks if a specified cmdline argument has been given.

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

1425{
1426 return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1427}

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 1172 of file BasicUtils/SessionReader.cpp.

1173{
1174 std::string vName = boost::to_upper_copy(pName);
1175 return m_functions.find(vName) != m_functions.end();
1176}
FunctionMap m_functions
Functions.

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 1181 of file BasicUtils/SessionReader.cpp.

1184{
1185 std::string vName = boost::to_upper_copy(pName);
1186
1187 // Check function exists
1188 auto it1 = m_functions.find(vName);
1189 if (it1 != m_functions.end())
1190 {
1191 pair<std::string, int> key(pVariable, pDomain);
1192 pair<std::string, int> defkey("*", pDomain);
1193 bool varExists = it1->second.find(key) != it1->second.end() ||
1194 it1->second.find(defkey) != it1->second.end();
1195 return varExists;
1196 }
1197 return false;
1198}

References m_functions.

◆ DefinesGlobalSysSolnInfo()

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

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

1042{
1043 auto iter = GetGloSysSolnList().find(pVariable);
1044 if (iter == GetGloSysSolnList().end())
1045 {
1046 return false;
1047 }
1048
1049 std::string vProperty = boost::to_upper_copy(pProperty);
1050
1051 auto iter1 = iter->second.find(vProperty);
1052 if (iter1 == iter->second.end())
1053 {
1054 return false;
1055 }
1056
1057 return true;
1058}
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.

◆ DefinesSolverInfo()

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

Checks if a solver info property is specified.

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

946{
947 std::string vName = boost::to_upper_copy(pName);
948 auto infoIter = m_solverInfo.find(vName);
949 return (infoIter != m_solverInfo.end());
950}

References m_solverInfo.

Referenced by GetSolverInfoAsEnum(), and MatchSolverInfo().

◆ DefinesTag()

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

Checks if a specified tag is defined.

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

1388{
1389 std::string vName = boost::to_upper_copy(pName);
1390 return m_tags.find(vName) != m_tags.end();
1391}

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 1084 of file BasicUtils/SessionReader.cpp.

1085{
1086 return m_timeIntScheme.method != "";
1087}
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.

◆ GetBackups()

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

Returns the backups.

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

1165{
1166 return m_backups;
1167}

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
Definition: SessionReader.h:71

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 348 of file SessionReader.h.

349 {
350 return m_cmdLineOptions.find(pName)->second.as<T>();
351 }

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.

◆ 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....
Definition: ErrorUtil.hpp:242

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:

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 1416 of file BasicUtils/SessionReader.cpp.

1417{
1418 return m_filters;
1419}

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 1203 of file BasicUtils/SessionReader.cpp.

1206{
1207 std::string vName = boost::to_upper_copy(pName);
1208 auto it1 = m_functions.find(vName);
1209
1210 ASSERTL0(it1 != m_functions.end(),
1211 std::string("No such function '") + pName +
1212 std::string("' has been defined in the session file."));
1213
1214 // Check for specific and wildcard definitions
1215 pair<std::string, int> key(pVariable, pDomain);
1216 pair<std::string, int> defkey("*", pDomain);
1217
1218 auto it2 = it1->second.find(key);
1219 auto it3 = it1->second.find(defkey);
1220 bool specific = it2 != it1->second.end();
1221 bool wildcard = it3 != it1->second.end();
1222
1223 // Check function is defined somewhere
1224 ASSERTL0(specific || wildcard,
1225 "No such variable " + pVariable + " in domain " +
1226 boost::lexical_cast<string>(pDomain) +
1227 " defined for function " + pName + " in session file.");
1228
1229 // If not specific, must be wildcard
1230 if (!specific)
1231 {
1232 it2 = it3;
1233 }
1234
1235 ASSERTL0((it2->second.m_type == eFunctionTypeExpression),
1236 std::string("Function is defined by a file."));
1237 return it2->second.m_expression;
1238}

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

Referenced by 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 1243 of file BasicUtils/SessionReader.cpp.

1246{
1247 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1248 return GetFunction(pName, m_variables[pVar], pDomain);
1249}
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 1302 of file BasicUtils/SessionReader.cpp.

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

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 1339 of file BasicUtils/SessionReader.cpp.

1342{
1343 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1344 return GetFunctionFilename(pName, m_variables[pVar], pDomain);
1345}
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 1350 of file BasicUtils/SessionReader.cpp.

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

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 1254 of file BasicUtils/SessionReader.cpp.

1257{
1258 std::string vName = boost::to_upper_copy(pName);
1259 auto it1 = m_functions.find(vName);
1260
1261 ASSERTL0(it1 != m_functions.end(),
1262 std::string("Function '") + pName + std::string("' not found."));
1263
1264 // Check for specific and wildcard definitions
1265 pair<std::string, int> key(pVariable, pDomain);
1266 pair<std::string, int> defkey("*", pDomain);
1267
1268 auto it2 = it1->second.find(key);
1269 auto it3 = it1->second.find(defkey);
1270 bool specific = it2 != it1->second.end();
1271 bool wildcard = it3 != it1->second.end();
1272
1273 // Check function is defined somewhere
1274 ASSERTL0(specific || wildcard,
1275 "No such variable " + pVariable + " in domain " +
1276 boost::lexical_cast<string>(pDomain) +
1277 " defined for function " + pName + " in session file.");
1278
1279 // If not specific, must be wildcard
1280 if (!specific)
1281 {
1282 it2 = it3;
1283 }
1284
1285 return it2->second.m_type;
1286}

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 1291 of file BasicUtils/SessionReader.cpp.

1294{
1295 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1296 return GetFunctionType(pName, m_variables[pVar], pDomain);
1297}
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 1101 of file BasicUtils/SessionReader.cpp.

1102{
1103 TiXmlElement *xmlGeom =
1104 m_xmlDoc->FirstChildElement("NEKTAR")->FirstChildElement("GEOMETRY");
1105 ASSERTL1(xmlGeom, "Failed to find a GEOMETRY section in m_xmlDoc");
1106
1107 TiXmlAttribute *attr = xmlGeom->FirstAttribute();
1108 while (attr)
1109 {
1110 std::string attrName(attr->Name());
1111 if (attrName == "HDF5FILE")
1112 {
1113 // there is a file pointer, therefore is HDF5
1114 return "HDF5";
1115 }
1116 // Get the next attribute.
1117 attr = attr->Next();
1118 }
1119
1120 // Check the VERTEX block. If this is compressed, assume the file is
1121 // compressed, otherwise assume uncompressed.
1122 TiXmlElement *element = xmlGeom->FirstChildElement("VERTEX");
1123 string IsCompressed;
1124 element->QueryStringAttribute("COMPRESSED", &IsCompressed);
1125
1126 if (IsCompressed.size() > 0)
1127 {
1128 return "XmlCompressed";
1129 }
1130
1131 // no file pointer or compressed, just standard xml
1132 return "Xml";
1133}

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 1063 of file BasicUtils/SessionReader.cpp.

1065{
1066 auto iter = GetGloSysSolnList().find(pVariable);
1067 ASSERTL0(iter != GetGloSysSolnList().end(),
1068 "Failed to find variable in GlobalSysSolnInfoList");
1069
1070 std::string vProperty = boost::to_upper_copy(pProperty);
1071 auto iter1 = iter->second.find(vProperty);
1072
1073 ASSERTL0(iter1 != iter->second.end(),
1074 "Failed to find property: " + vProperty +
1075 " in GlobalSysSolnInfoList");
1076
1077 return iter1->second;
1078}

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
Definition: SessionReader.h:77

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 2562 of file BasicUtils/SessionReader.cpp.

2563{
2564 return m_interpreter;
2565}

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.

◆ 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.

◆ 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.

◆ 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.

◆ 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 955 of file BasicUtils/SessionReader.cpp.

957{
958 std::string vProperty = boost::to_upper_copy(pProperty);
959 auto iter = m_solverInfo.find(vProperty);
960
961 ASSERTL1(iter != m_solverInfo.end(),
962 "Unable to find requested property: " + pProperty);
963
964 return iter->second;
965}

References ASSERTL1, and m_solverInfo.

Referenced by 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 497 of file SessionReader.h.

499{
500 std::string vName = boost::to_upper_copy(pName);
502 "Solver info '" + pName + "' not defined.");
503
504 std::string vValue = GetSolverInfo(vName);
505 auto x = GetSolverInfoEnums().find(vName);
506 ASSERTL0(x != GetSolverInfoEnums().end(),
507 "Enum for SolverInfo property '" + pName + "' not found.");
508
509 auto y = x->second.find(vValue);
510 ASSERTL0(y != x->second.end(),
511 "Value of SolverInfo property '" + pName + "' is invalid.");
512
513 return T(y->second);
514}
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
Definition: SessionReader.h:56

Referenced by ReadSolverInfo(), RegisterDefaultSolverInfo(), 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
Definition: SessionReader.h:74

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 1405 of file BasicUtils/SessionReader.cpp.

1406{
1407 std::string vName = boost::to_upper_copy(pName);
1408 auto vTagIterator = m_tags.find(vName);
1409 ASSERTL0(vTagIterator != m_tags.end(), "Requested tag does not exist.");
1410 return vTagIterator->second;
1411}

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 1093 of file BasicUtils/SessionReader.cpp.

1094{
1095 return m_timeIntScheme;
1096}

References m_timeIntScheme.

◆ GetTimeLevel()

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

Get time level (Parallel-in-Time)

Definition at line 374 of file SessionReader.h.

375 {
376 return m_timeLevel;
377 }

References m_timeLevel.

◆ GetUpdateOptFile()

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

Get bool to update optimisation file.

Definition at line 362 of file SessionReader.h.

363 {
364 return m_updateOptFile;
365 }
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 520 of file SessionReader.h.

522{
523 std::string vName = boost::to_upper_copy(pName);
524
525 auto x = GetSolverInfoEnums().find(vName);
526 ASSERTL0(x != GetSolverInfoEnums().end(),
527 "Enum for property '" + pName + "' not found.");
528
529 auto y = x->second.find(pValue);
530 ASSERTL0(y != x->second.end(),
531 "Value of property '" + pValue + "' is invalid.");
532 return T(y->second);
533}

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 1138 of file BasicUtils/SessionReader.cpp.

1139{
1140 ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1141 return m_variables[idx];
1142}

References ASSERTL0, and m_variables.

◆ GetVariables()

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

Returns the names of all variables.

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

1157{
1158 return m_variables;
1159}

References m_variables.

◆ 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 1432 of file BasicUtils/SessionReader.cpp.

1435{
1436 if (Element && Element->FirstChildElement("TIMELEVEL"))
1437 {
1438 Element = Element->FirstChildElement("TIMELEVEL");
1439 std::string timeLevelStr;
1440 while (Element)
1441 {
1442 std::stringstream tagcontent;
1443 tagcontent << *Element;
1444 ASSERTL0(Element->Attribute("VALUE"),
1445 "Missing LEVEL attribute in solver info "
1446 "XML element: \n\t'" +
1447 tagcontent.str() + "'");
1448 timeLevelStr = Element->Attribute("VALUE");
1449 ASSERTL0(!timeLevelStr.empty(),
1450 "LEVEL attribute must be non-empty in XML "
1451 "element: \n\t'" +
1452 tagcontent.str() + "'");
1453 if (stoi(timeLevelStr) == timeLevel)
1454 {
1455 break;
1456 }
1457 Element = Element->NextSiblingElement("TIMELEVEL");
1458 }
1459 if (enableCheck)
1460 {
1461 ASSERTL0(stoi(timeLevelStr) == timeLevel,
1462 "TIMELEVEL value " + std::to_string(timeLevel) +
1463 " not found in solver info "
1464 "XML element: \n\t'");
1465 }
1466 }
1467}

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().

◆ 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 1472 of file BasicUtils/SessionReader.cpp.

1474{
1475 if (pFilename.size() > 3 &&
1476 pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1477 {
1478 ifstream file(pFilename.c_str(), ios_base::in | ios_base::binary);
1479 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1480 stringstream ss;
1481 io::filtering_streambuf<io::input> in;
1482 in.push(io::gzip_decompressor());
1483 in.push(file);
1484 try
1485 {
1486 io::copy(in, ss);
1487 ss >> (*pDoc);
1488 }
1489 catch (io::gzip_error &)
1490 {
1492 "Error: File '" + pFilename + "' is corrupt.");
1493 }
1494 }
1495 else if (pFilename.size() > 4 &&
1496 pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1497 {
1498 fs::path pdirname(pFilename);
1499 boost::format pad("P%1$07d.xml");
1500 pad % m_comm->GetSpaceComm()->GetRank();
1501 fs::path pRankFilename(pad.str());
1502 fs::path fullpath = pdirname / pRankFilename;
1503
1504 ifstream file(PortablePath(fullpath).c_str());
1505 ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1506 file >> (*pDoc);
1507 }
1508 else
1509 {
1510 ifstream file(pFilename.c_str());
1511 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1512 file >> (*pDoc);
1513 }
1514}
def copy(self)
Definition: pycml.py:2663
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56

References ASSERTL0, CellMLToNektar.pycml::copy(), Nektar::ErrorUtil::efatal, CellMLToNektar.pycml::format, m_comm, NEKERROR, and Nektar::LibUtilities::PortablePath().

Referenced by MergeDoc().

◆ LoadParameter() [1/6]

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}
double NekDouble

References ASSERTL0, and m_parameters.

◆ LoadParameter() [2/6]

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/6]

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

Load a double precision parameter.

Definition at line 887 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 pVar = paramIter->second;
895}

References ASSERTL0, and m_parameters.

◆ LoadParameter() [4/6]

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 900 of file BasicUtils/SessionReader.cpp.

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

References m_parameters.

◆ LoadParameter() [5/6]

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

Load an size_t parameter.

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

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

References ASSERTL0, and m_parameters.

◆ LoadParameter() [6/6]

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 868 of file BasicUtils/SessionReader.cpp.

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

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 985 of file BasicUtils/SessionReader.cpp.

987{
988 std::string vName = boost::to_upper_copy(pName);
989 auto infoIter = m_solverInfo.find(vName);
990 if (infoIter != m_solverInfo.end())
991 {
992 pVar = infoIter->second;
993 }
994 else
995 {
996 pVar = pDefault;
997 }
998}

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 1022 of file BasicUtils/SessionReader.cpp.

1024{
1025 if (DefinesSolverInfo(pName))
1026 {
1027 std::string vName = boost::to_upper_copy(pName);
1028 auto iter = m_solverInfo.find(vName);
1029 if (iter != m_solverInfo.end())
1030 {
1031 return boost::iequals(iter->second, pTrueVal);
1032 }
1033 }
1034 return false;
1035}

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 1003 of file BasicUtils/SessionReader.cpp.

1006{
1007 std::string vName = boost::to_upper_copy(pName);
1008 auto infoIter = m_solverInfo.find(vName);
1009 if (infoIter != m_solverInfo.end())
1010 {
1011 pVar = boost::iequals(infoIter->second, pTrueVal);
1012 }
1013 else
1014 {
1015 pVar = pDefault;
1016 }
1017}

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 487 of file SessionReader.h.

489{
490 return (GetSolverInfoAsEnum<T>(name) == trueval);
491}

References CellMLToNektar.pycml::name.

◆ 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 1519 of file BasicUtils/SessionReader.cpp.

1521{
1522 ASSERTL0(pFilenames.size() > 0, "No filenames for merging.");
1523
1524 // Read the first document
1525 TiXmlDocument *vMainDoc = new TiXmlDocument;
1526 LoadDoc(pFilenames[0], vMainDoc);
1527
1528 TiXmlHandle vMainHandle(vMainDoc);
1529 TiXmlElement *vMainNektar =
1530 GetChildElementOrThrow(pFilenames[0], "NEKTAR", vMainHandle);
1531
1532 // Read all subsequent XML documents.
1533 // For each element within the NEKTAR tag, use it to replace the
1534 // version already present in the loaded XML data.
1535 for (int i = 1; i < pFilenames.size(); ++i)
1536 {
1537 if ((pFilenames[i].compare(pFilenames[i].size() - 3, 3, "xml") == 0) ||
1538 (pFilenames[i].compare(pFilenames[i].size() - 6, 6, "xml.gz") ==
1539 0) ||
1540 (pFilenames[i].compare(pFilenames[i].size() - 3, 3, "opt") == 0))
1541 {
1542 TiXmlDocument *vTempDoc = new TiXmlDocument;
1543 LoadDoc(pFilenames[i], vTempDoc);
1544
1545 TiXmlHandle docHandle(vTempDoc);
1546 TiXmlElement *vTempNektar =
1547 GetChildElementOrThrow(pFilenames[i], "NEKTAR", docHandle);
1548 TiXmlElement *p = vTempNektar->FirstChildElement();
1549
1550 while (p)
1551 {
1552 TiXmlElement *vMainEntry =
1553 vMainNektar->FirstChildElement(p->Value());
1554
1555 // First check if the new item is in fact blank
1556 // replace if it is a COLLECTIONS section however.
1557 if (!p->FirstChild() && vMainEntry &&
1558 !boost::iequals(p->Value(), "COLLECTIONS"))
1559 {
1560 std::string warningmsg =
1561 "File " + pFilenames[i] + " contains " +
1562 "an empty XML element " + std::string(p->Value()) +
1563 " which will be ignored.";
1564 NEKERROR(ErrorUtil::ewarning, warningmsg.c_str());
1565 }
1566 else
1567 {
1568 if (vMainEntry)
1569 {
1570 vMainNektar->RemoveChild(vMainEntry);
1571 }
1572 TiXmlElement *q = new TiXmlElement(*p);
1573 vMainNektar->LinkEndChild(q);
1574 }
1575 p = p->NextSiblingElement();
1576 }
1577 delete vTempDoc;
1578 }
1579 }
1580 return vMainDoc;
1581}
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(), NEKERROR, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

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.
def description(self, force=False, cellml=False)
Definition: pycml.py:2670
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 Nektar::UnitTests::d(), CellMLToNektar.pycml::description(), GetCmdLineArgMap(), Nektar::NekConstants::kGitBranch, Nektar::NekConstants::kGitSha1, m_backups, m_cmdLineOptions, m_updateOptFile, m_verbose, NEKTAR_VERSION, and CellMLToNektar.cellml_metadata::p.

Referenced by SessionReader().

◆ ParseDocument()

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

Loads and parses the specified file.

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

1587{
1588 // Check we actually have a document loaded.
1589 ASSERTL0(m_xmlDoc, "No XML document loaded.");
1590
1591 TiXmlHandle docHandle(m_xmlDoc);
1592 TiXmlElement *e;
1593
1594 // Look for all data in CONDITIONS block.
1595 e = docHandle.FirstChildElement("NEKTAR")
1596 .FirstChildElement("CONDITIONS")
1597 .Element();
1598
1599 // Read the various sections of the CONDITIONS block
1600 ReadParameters(e);
1601 ReadSolverInfo(e);
1604 ReadVariables(e);
1605 ReadFunctions(e);
1606
1607 // Look for all data in FILTERS block.
1608 e = docHandle.FirstChildElement("NEKTAR")
1609 .FirstChildElement("FILTERS")
1610 .Element();
1611
1612 // Read the various sections of the FILTERS block
1613 ReadFilters(e);
1614}
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 2442 of file BasicUtils/SessionReader.cpp.

2444{
2445 /// Pull out lhs and rhs and eliminate any spaces.
2446 size_t beg = line.find_first_not_of(" ");
2447 size_t end = line.find_first_of("=");
2448 // Check for no parameter name
2449 if (beg == end)
2450 {
2451 throw 1;
2452 }
2453 // Check for no parameter value
2454 if (end != line.find_last_of("="))
2455 {
2456 throw 1;
2457 }
2458 // Check for no equals sign
2459 if (end == std::string::npos)
2460 {
2461 throw 1;
2462 }
2463
2464 lhs = line.substr(line.find_first_not_of(" "), end - beg);
2465 lhs = lhs.substr(0, lhs.find_last_not_of(" ") + 1);
2466 rhs = line.substr(line.find_last_of("=") + 1);
2467 rhs = rhs.substr(rhs.find_first_not_of(" "));
2468 rhs = rhs.substr(0, rhs.find_last_not_of(" ") + 1);
2469}

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().

◆ 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 1646 of file BasicUtils/SessionReader.cpp.

1647{
1648 if (m_comm->GetSize() > 1)
1649 {
1650 int nProcZ = 1;
1651 int nProcY = 1;
1652 int nProcX = 1;
1653 int nStripZ = 1;
1654 int nTime = 1;
1655 if (DefinesCmdLineArgument("npx"))
1656 {
1657 nProcX = GetCmdLineArgument<int>("npx");
1658 }
1659 if (DefinesCmdLineArgument("npy"))
1660 {
1661 nProcY = GetCmdLineArgument<int>("npy");
1662 }
1663 if (DefinesCmdLineArgument("npz"))
1664 {
1665 nProcZ = GetCmdLineArgument<int>("npz");
1666 }
1667 if (DefinesCmdLineArgument("nsz"))
1668 {
1669 nStripZ = GetCmdLineArgument<int>("nsz");
1670 }
1671 if (DefinesCmdLineArgument("npt"))
1672 {
1673 nTime = GetCmdLineArgument<int>("npt");
1674 }
1675 ASSERTL0(m_comm->GetSize() % nTime == 0,
1676 "Cannot exactly partition time using npt value.");
1677 ASSERTL0((m_comm->GetSize() / nTime) % (nProcZ * nProcY * nProcX) == 0,
1678 "Cannot exactly partition using PROC_Z value.");
1679 ASSERTL0(nProcZ % nProcY == 0,
1680 "Cannot exactly partition using PROC_Y value.");
1681 ASSERTL0(nProcY % nProcX == 0,
1682 "Cannot exactly partition using PROC_X value.");
1683
1684 // Number of processes associated with the spectral method
1685 int nProcSm = nProcZ * nProcY * nProcX;
1686
1687 // Number of processes associated with the spectral element
1688 // method.
1689 int nProcSem = m_comm->GetSize() / nTime / nProcSm;
1690
1691 m_comm->SplitComm(nProcSm, nProcSem, nTime);
1692 m_comm->GetColumnComm()->SplitComm(nProcZ / nStripZ, nStripZ);
1693 m_comm->GetColumnComm()->GetColumnComm()->SplitComm((nProcY * nProcX),
1694 nProcZ / nStripZ);
1695 m_comm->GetColumnComm()->GetColumnComm()->GetColumnComm()->SplitComm(
1696 nProcX, nProcY);
1697 }
1698}

References ASSERTL0, DefinesCmdLineArgument(), and m_comm.

Referenced by SessionReader().

◆ ReadFilters()

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

Reads the FILTERS section of the XML document.

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

2396{
2397 if (!filters)
2398 {
2399 return;
2400 }
2401
2402 m_filters.clear();
2403
2404 TiXmlElement *filter = filters->FirstChildElement("FILTER");
2405
2406 while (filter)
2407 {
2408 ASSERTL0(filter->Attribute("TYPE"),
2409 "Missing attribute 'TYPE' for filter.");
2410 std::string typeStr = filter->Attribute("TYPE");
2411
2412 std::map<std::string, std::string> vParams;
2413
2414 TiXmlElement *element = filter;
2415 GetXMLElementTimeLevel(element, m_timeLevel, false);
2416 TiXmlElement *param = element->FirstChildElement("PARAM");
2417 while (param)
2418 {
2419 ASSERTL0(param->Attribute("NAME"),
2420 "Missing attribute 'NAME' for parameter in filter " +
2421 typeStr + "'.");
2422 std::string nameStr = param->Attribute("NAME");
2423
2424 ASSERTL0(param->GetText(), "Empty value string for param.");
2425 std::string valueStr = param->GetText();
2426
2427 vParams[nameStr] = valueStr;
2428
2429 param = param->NextSiblingElement("PARAM");
2430 }
2431
2432 m_filters.push_back(
2433 std::pair<std::string, FilterParams>(typeStr, vParams));
2434
2435 filter = filter->NextSiblingElement("FILTER");
2436 }
2437}
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)

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

Referenced by ParseDocument().

◆ ReadFunctions()

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

Reads the FUNCTIONS section of the XML document.

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

2171{
2172 m_functions.clear();
2173
2174 if (!conditions)
2175 {
2176 return;
2177 }
2178
2179 // Scan through conditions section looking for functions.
2180 TiXmlElement *function = conditions->FirstChildElement("FUNCTION");
2181
2182 while (function)
2183 {
2184 stringstream tagcontent;
2185 tagcontent << *function;
2186
2187 // Every function must have a NAME attribute
2188 ASSERTL0(function->Attribute("NAME"),
2189 "Functions must have a NAME attribute defined in XML "
2190 "element: \n\t'" +
2191 tagcontent.str() + "'");
2192 std::string functionStr = function->Attribute("NAME");
2193 ASSERTL0(!functionStr.empty(),
2194 "Functions must have a non-empty name in XML "
2195 "element: \n\t'" +
2196 tagcontent.str() + "'");
2197
2198 // Store function names in uppercase to remain case-insensitive.
2199 boost::to_upper(functionStr);
2200
2201 // Retrieve first entry (variable, or file)
2202 TiXmlElement *element = function;
2203 GetXMLElementTimeLevel(element, m_timeLevel, false);
2204 TiXmlElement *variable = element->FirstChildElement();
2205
2206 // Create new function structure with default type of none.
2207 FunctionVariableMap functionVarMap;
2208
2209 // Process all entries in the function block
2210 while (variable)
2211 {
2212 FunctionVariableDefinition funcDef;
2213 std::string conditionType = variable->Value();
2214
2215 // If no var is specified, assume wildcard
2216 std::string variableStr;
2217 if (!variable->Attribute("VAR"))
2218 {
2219 variableStr = "*";
2220 }
2221 else
2222 {
2223 variableStr = variable->Attribute("VAR");
2224 }
2225
2226 // Parse list of variables
2227 std::vector<std::string> variableList;
2228 ParseUtils::GenerateVector(variableStr, variableList);
2229
2230 // If no domain is specified, put to 0
2231 std::string domainStr;
2232 if (!variable->Attribute("DOMAIN"))
2233 {
2234 domainStr = "0";
2235 }
2236 else
2237 {
2238 domainStr = variable->Attribute("DOMAIN");
2239 }
2240
2241 // Parse list of domains
2242 std::vector<std::string> varSplit;
2243 std::vector<unsigned int> domainList;
2244 ParseUtils::GenerateSeqVector(domainStr, domainList);
2245
2246 // if no evars is specified, put "x y z t"
2247 std::string evarsStr = "x y z t";
2248 if (variable->Attribute("EVARS"))
2249 {
2250 evarsStr =
2251 evarsStr + std::string(" ") + variable->Attribute("EVARS");
2252 }
2253
2254 // Expressions are denoted by E
2255 if (conditionType == "E")
2256 {
2257 funcDef.m_type = eFunctionTypeExpression;
2258
2259 // Expression must have a VALUE.
2260 ASSERTL0(variable->Attribute("VALUE"),
2261 "Attribute VALUE expected for function '" +
2262 functionStr + "'.");
2263 std::string fcnStr = variable->Attribute("VALUE");
2264 ASSERTL0(!fcnStr.empty(),
2265 (std::string("Expression for var: ") + variableStr +
2266 std::string(" must be specified."))
2267 .c_str());
2268
2269 // set expression
2270 funcDef.m_expression =
2272 m_interpreter, fcnStr, evarsStr);
2273 }
2274
2275 // Files are denoted by F
2276 else if (conditionType == "F")
2277 {
2278 // Check if transient or not
2279 if (variable->Attribute("TIMEDEPENDENT") &&
2280 boost::lexical_cast<bool>(
2281 variable->Attribute("TIMEDEPENDENT")))
2282 {
2283 funcDef.m_type = eFunctionTypeTransientFile;
2284 }
2285 else
2286 {
2287 funcDef.m_type = eFunctionTypeFile;
2288 }
2289
2290 // File must have a FILE.
2291 ASSERTL0(variable->Attribute("FILE"),
2292 "Attribute FILE expected for function '" +
2293 functionStr + "'.");
2294 std::string filenameStr = variable->Attribute("FILE");
2295 ASSERTL0(!filenameStr.empty(),
2296 "A filename must be specified for the FILE "
2297 "attribute of function '" +
2298 functionStr + "'.");
2299
2300 std::vector<std::string> fSplit;
2301 boost::split(fSplit, filenameStr, boost::is_any_of(":"));
2302 ASSERTL0(fSplit.size() == 1 || fSplit.size() == 2,
2303 "Incorrect filename specification in function " +
2304 functionStr +
2305 "'. "
2306 "Specify variables inside file as: "
2307 "filename:var1,var2");
2308
2309 // set the filename
2310 fs::path fullpath = fSplit[0];
2311 fs::path ftype = fullpath.extension();
2312 if (fullpath.parent_path().extension() == ".pit")
2313 {
2314 string filename = fullpath.stem().string();
2315 fullpath = fullpath.parent_path();
2316 size_t start = filename.find_last_of("_") + 1;
2317 int index =
2318 atoi(filename.substr(start, filename.size()).c_str());
2319 fullpath /= filename.substr(0, start) +
2320 std::to_string(
2321 index + m_comm->GetTimeComm()->GetRank()) +
2322 ftype.string();
2323 }
2324 funcDef.m_filename = fullpath.string();
2325
2326 if (fSplit.size() == 2)
2327 {
2328 ASSERTL0(variableList[0] != "*",
2329 "Filename variable mapping not valid "
2330 "when using * as a variable inside "
2331 "function '" +
2332 functionStr + "'.");
2333
2334 boost::split(varSplit, fSplit[1], boost::is_any_of(","));
2335 ASSERTL0(varSplit.size() == variableList.size(),
2336 "Filename variables should contain the "
2337 "same number of variables defined in "
2338 "VAR in function " +
2339 functionStr + "'.");
2340 }
2341 }
2342
2343 // Nothing else supported so throw an error
2344 else
2345 {
2346 stringstream tagcontent;
2347 tagcontent << *variable;
2348
2350 "Identifier " + conditionType + " in function " +
2351 std::string(function->Attribute("NAME")) +
2352 " is not recognised in XML element: \n\t'" +
2353 tagcontent.str() + "'");
2354 }
2355
2356 // Add variables to function
2357 for (unsigned int i = 0; i < variableList.size(); ++i)
2358 {
2359 for (unsigned int j = 0; j < domainList.size(); ++j)
2360 {
2361 // Check it has not already been defined
2362 pair<std::string, int> key(variableList[i], domainList[j]);
2363 auto fcnsIter = functionVarMap.find(key);
2364 ASSERTL0(fcnsIter == functionVarMap.end(),
2365 "Error setting expression '" + variableList[i] +
2366 " in domain " + std::to_string(domainList[j]) +
2367 "' in function '" + functionStr +
2368 "'. "
2369 "Expression has already been defined.");
2370
2371 if (varSplit.size() > 0)
2372 {
2373 FunctionVariableDefinition funcDef2 = funcDef;
2374 funcDef2.m_fileVariable = varSplit[i];
2375 functionVarMap[key] = funcDef2;
2376 }
2377 else
2378 {
2379 functionVarMap[key] = funcDef;
2380 }
2381 }
2382 }
2383 variable = variable->NextSiblingElement();
2384 }
2385
2386 // Add function definition to map
2387 m_functions[functionStr] = functionVarMap;
2388 function = function->NextSiblingElement("FUNCTION");
2389 }
2390}
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
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.
Definition: ParseUtils.cpp:104
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 1861 of file BasicUtils/SessionReader.cpp.

1862{
1863 GetGloSysSolnList().clear();
1864
1865 if (!conditions)
1866 {
1867 return;
1868 }
1869
1870 TiXmlElement *GlobalSys =
1871 conditions->FirstChildElement("GLOBALSYSSOLNINFO");
1873
1874 if (!GlobalSys)
1875 {
1876 return;
1877 }
1878
1879 TiXmlElement *VarInfo = GlobalSys->FirstChildElement("V");
1880
1881 while (VarInfo)
1882 {
1883 std::stringstream tagcontent;
1884 tagcontent << *VarInfo;
1885
1886 ASSERTL0(VarInfo->Attribute("VAR"),
1887 "Missing VAR attribute in GobalSysSolnInfo XML "
1888 "element: \n\t'" +
1889 tagcontent.str() + "'");
1890 std::string VarList = VarInfo->Attribute("VAR");
1891 ASSERTL0(!VarList.empty(),
1892 "VAR attribute must be non-empty in XML element:\n\t'" +
1893 tagcontent.str() + "'");
1894
1895 // generate a list of variables.
1896 std::vector<std::string> varStrings;
1897 bool valid = ParseUtils::GenerateVector(VarList, varStrings);
1898
1899 ASSERTL0(valid, "Unable to process list of variable in XML "
1900 "element \n\t'" +
1901 tagcontent.str() + "'");
1902
1903 if (varStrings.size())
1904 {
1905 TiXmlElement *SysSolnInfo = VarInfo->FirstChildElement("I");
1906
1907 while (SysSolnInfo)
1908 {
1909 tagcontent.clear();
1910 tagcontent << *SysSolnInfo;
1911 // read the property name
1912 ASSERTL0(SysSolnInfo->Attribute("PROPERTY"),
1913 "Missing PROPERTY attribute in "
1914 "GlobalSysSolnInfo for variable(s) '" +
1915 VarList + "' in XML element: \n\t'" +
1916 tagcontent.str() + "'");
1917 std::string SysSolnProperty =
1918 SysSolnInfo->Attribute("PROPERTY");
1919 ASSERTL0(!SysSolnProperty.empty(),
1920 "GlobalSysSolnIno properties must have a "
1921 "non-empty name for variable(s) : '" +
1922 VarList + "' in XML element: \n\t'" +
1923 tagcontent.str() + "'");
1924
1925 // make sure that solver property is capitalised
1926 std::string SysSolnPropertyUpper =
1927 boost::to_upper_copy(SysSolnProperty);
1928
1929 // read the value
1930 ASSERTL0(SysSolnInfo->Attribute("VALUE"),
1931 "Missing VALUE attribute in GlobalSysSolnInfo "
1932 "for variable(s) '" +
1933 VarList + "' in XML element: \n\t" +
1934 tagcontent.str() + "'");
1935 std::string SysSolnValue = SysSolnInfo->Attribute("VALUE");
1936 ASSERTL0(!SysSolnValue.empty(),
1937 "GlobalSysSolnInfo properties must have a "
1938 "non-empty value for variable(s) '" +
1939 VarList + "' in XML element: \n\t'" +
1940 tagcontent.str() + "'");
1941
1942 // Store values under variable map.
1943 for (int i = 0; i < varStrings.size(); ++i)
1944 {
1945 auto x = GetGloSysSolnList().find(varStrings[i]);
1946 if (x == GetGloSysSolnList().end())
1947 {
1948 (GetGloSysSolnList()[varStrings[i]])
1949 [SysSolnPropertyUpper] = SysSolnValue;
1950 }
1951 else
1952 {
1953 x->second[SysSolnPropertyUpper] = SysSolnValue;
1954 }
1955 }
1956 SysSolnInfo = SysSolnInfo->NextSiblingElement("I");
1957 }
1958 VarInfo = VarInfo->NextSiblingElement("V");
1959 }
1960 }
1961
1962 if (m_verbose && GetGloSysSolnList().size() > 0 && m_comm)
1963 {
1964 if (m_comm->GetRank() == 0)
1965 {
1966 cout << "GlobalSysSoln Info:" << endl;
1967
1968 for (auto &x : GetGloSysSolnList())
1969 {
1970 cout << "\t Variable: " << x.first << endl;
1971
1972 for (auto &y : x.second)
1973 {
1974 cout << "\t\t " << y.first << " = " << y.second << endl;
1975 }
1976 }
1977 cout << endl;
1978 }
1979 }
1980}

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 1703 of file BasicUtils/SessionReader.cpp.

1704{
1705 m_parameters.clear();
1706
1707 if (!conditions)
1708 {
1709 return;
1710 }
1711
1712 TiXmlElement *parametersElement =
1713 conditions->FirstChildElement("PARAMETERS");
1714 GetXMLElementTimeLevel(parametersElement, m_timeLevel);
1715
1716 // See if we have parameters defined. They are optional so we go on
1717 // if not.
1718 if (parametersElement)
1719 {
1720 TiXmlElement *parameter = parametersElement->FirstChildElement("P");
1721
1722 ParameterMap caseSensitiveParameters;
1723
1724 // Multiple nodes will only occur if there is a comment in
1725 // between definitions.
1726 while (parameter)
1727 {
1728 stringstream tagcontent;
1729 tagcontent << *parameter;
1730 TiXmlNode *node = parameter->FirstChild();
1731
1732 while (node && node->Type() != TiXmlNode::TINYXML_TEXT)
1733 {
1734 node = node->NextSibling();
1735 }
1736
1737 if (node)
1738 {
1739 // Format is "paramName = value"
1740 std::string line = node->ToText()->Value(), lhs, rhs;
1741
1742 try
1743 {
1744 ParseEquals(line, lhs, rhs);
1745 }
1746 catch (...)
1747 {
1749 "Syntax error in parameter expression '" + line +
1750 "' in XML element: \n\t'" + tagcontent.str() +
1751 "'");
1752 }
1753
1754 // We want the list of parameters to have their RHS
1755 // evaluated, so we use the expression evaluator to do
1756 // the dirty work.
1757 if (!lhs.empty() && !rhs.empty())
1758 {
1759 NekDouble value = 0.0;
1760 try
1761 {
1762 LibUtilities::Equation expession(m_interpreter, rhs);
1763 value = expession.Evaluate();
1764 }
1765 catch (const std::runtime_error &)
1766 {
1768 "Error evaluating parameter expression"
1769 " '" +
1770 rhs + "' in XML element: \n\t'" +
1771 tagcontent.str() + "'");
1772 }
1773 m_interpreter->SetParameter(lhs, value);
1774 caseSensitiveParameters[lhs] = value;
1775 boost::to_upper(lhs);
1776 m_parameters[lhs] = value;
1777 }
1778 }
1779 parameter = parameter->NextSiblingElement();
1780 }
1781 }
1782}
std::map< std::string, NekDouble > ParameterMap
Definition: SessionReader.h:57

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 1787 of file BasicUtils/SessionReader.cpp.

1788{
1789 m_solverInfo.clear();
1791
1792 if (!conditions)
1793 {
1794 return;
1795 }
1796
1797 TiXmlElement *solverInfoElement =
1798 conditions->FirstChildElement("SOLVERINFO");
1799 GetXMLElementTimeLevel(solverInfoElement, m_timeLevel);
1800
1801 if (solverInfoElement)
1802 {
1803 TiXmlElement *solverInfo = solverInfoElement->FirstChildElement("I");
1804
1805 while (solverInfo)
1806 {
1807 std::stringstream tagcontent;
1808 tagcontent << *solverInfo;
1809 // read the property name
1810 ASSERTL0(solverInfo->Attribute("PROPERTY"),
1811 "Missing PROPERTY attribute in solver info "
1812 "XML element: \n\t'" +
1813 tagcontent.str() + "'");
1814 std::string solverProperty = solverInfo->Attribute("PROPERTY");
1815 ASSERTL0(!solverProperty.empty(),
1816 "PROPERTY attribute must be non-empty in XML "
1817 "element: \n\t'" +
1818 tagcontent.str() + "'");
1819
1820 // make sure that solver property is capitalised
1821 std::string solverPropertyUpper =
1822 boost::to_upper_copy(solverProperty);
1823
1824 // read the value
1825 ASSERTL0(solverInfo->Attribute("VALUE"),
1826 "Missing VALUE attribute in solver info "
1827 "XML element: \n\t'" +
1828 tagcontent.str() + "'");
1829 std::string solverValue = solverInfo->Attribute("VALUE");
1830 ASSERTL0(!solverValue.empty(),
1831 "VALUE attribute must be non-empty in XML "
1832 "element: \n\t'" +
1833 tagcontent.str() + "'");
1834
1835 // Set Variable
1836 m_solverInfo[solverPropertyUpper] = solverValue;
1837 solverInfo = solverInfo->NextSiblingElement("I");
1838 }
1839 }
1840
1841 if (m_comm && m_comm->GetRowComm()->GetSize() > 1)
1842 {
1843 ASSERTL0(
1844 m_solverInfo["GLOBALSYSSOLN"] == "IterativeFull" ||
1845 m_solverInfo["GLOBALSYSSOLN"] == "IterativeStaticCond" ||
1846 m_solverInfo["GLOBALSYSSOLN"] ==
1847 "IterativeMultiLevelStaticCond" ||
1848 m_solverInfo["GLOBALSYSSOLN"] == "XxtFull" ||
1849 m_solverInfo["GLOBALSYSSOLN"] == "XxtStaticCond" ||
1850 m_solverInfo["GLOBALSYSSOLN"] == "XxtMultiLevelStaticCond" ||
1851 m_solverInfo["GLOBALSYSSOLN"] == "PETScFull" ||
1852 m_solverInfo["GLOBALSYSSOLN"] == "PETScStaticCond" ||
1853 m_solverInfo["GLOBALSYSSOLN"] == "PETScMultiLevelStaticCond",
1854 "A parallel solver must be used when run in parallel.");
1855 }
1856}

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 1985 of file BasicUtils/SessionReader.cpp.

1986{
1987 if (!conditions)
1988 {
1989 return;
1990 }
1991
1992 TiXmlElement *timeInt =
1993 conditions->FirstChildElement("TIMEINTEGRATIONSCHEME");
1995
1996 if (!timeInt)
1997 {
1998 return;
1999 }
2000
2001 TiXmlElement *method = timeInt->FirstChildElement("METHOD");
2002 TiXmlElement *variant = timeInt->FirstChildElement("VARIANT");
2003 TiXmlElement *order = timeInt->FirstChildElement("ORDER");
2004 TiXmlElement *params = timeInt->FirstChildElement("FREEPARAMETERS");
2005
2006 // Only the method and order are required.
2007 ASSERTL0(method, "Missing METHOD tag inside "
2008 "TIMEINTEGRATIONSCHEME.");
2009 ASSERTL0(order, "Missing ORDER tag inside "
2010 "TIMEINTEGRATIONSCHEME.");
2011
2012 m_timeIntScheme.method = method->GetText();
2013
2014 std::string orderStr = order->GetText();
2015
2016 // Only the method and order are required.
2017 ASSERTL0(m_timeIntScheme.method.size() > 0,
2018 "Empty text inside METHOD tag in TIMEINTEGRATIONSCHEME.");
2019 ASSERTL0(orderStr.size() > 0,
2020 "Empty text inside ORDER tag in TIMEINTEGRATIONSCHEME.");
2021 try
2022 {
2023 m_timeIntScheme.order = boost::lexical_cast<unsigned int>(orderStr);
2024 }
2025 catch (...)
2026 {
2027 NEKERROR(ErrorUtil::efatal, "In ORDER tag, unable to convert "
2028 "string '" +
2029 orderStr + "' to an unsigned integer.");
2030 }
2031
2032 if (variant)
2033 {
2034 m_timeIntScheme.variant = variant->GetText();
2035 }
2036
2037 if (params)
2038 {
2039 std::string paramsStr = params->GetText();
2040 ASSERTL0(paramsStr.size() > 0,
2041 "Empty text inside FREEPARAMETERS tag in "
2042 "TIMEINTEGRATIONSCHEME.");
2043
2044 std::vector<std::string> pSplit;
2045 boost::split(pSplit, paramsStr, boost::is_any_of(" "));
2046
2047 m_timeIntScheme.freeParams.resize(pSplit.size());
2048 for (size_t i = 0; i < pSplit.size(); ++i)
2049 {
2050 try
2051 {
2052 m_timeIntScheme.freeParams[i] = std::stod(pSplit[i]);
2053 }
2054 catch (...)
2055 {
2056 NEKERROR(ErrorUtil::efatal, "In FREEPARAMETERS tag, "
2057 "unable to convert string '" +
2058 pSplit[i] +
2059 "' "
2060 "to a floating-point value.");
2061 }
2062 }
2063 }
2064
2065 if (m_verbose && m_comm)
2066 {
2067 if (m_comm->GetRank() == 0)
2068 {
2069 cout << "Trying to use time integration scheme:" << endl;
2070 cout << "\t Method : " << m_timeIntScheme.method << endl;
2071 cout << "\t Variant: " << m_timeIntScheme.variant << endl;
2072 cout << "\t Order : " << m_timeIntScheme.order << endl;
2073
2074 if (m_timeIntScheme.freeParams.size() > 0)
2075 {
2076 cout << "\t Params :";
2077 for (auto &x : m_timeIntScheme.freeParams)
2078 {
2079 cout << " " << x;
2080 }
2081 cout << endl;
2082 }
2083 }
2084 }
2085}
std::vector< NekDouble > freeParams
Definition: SessionReader.h:87

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 2090 of file BasicUtils/SessionReader.cpp.

2091{
2092 m_variables.clear();
2093
2094 if (!conditions)
2095 {
2096 return;
2097 }
2098
2099 TiXmlElement *variablesElement = conditions->FirstChildElement("VARIABLES");
2100 GetXMLElementTimeLevel(variablesElement, m_timeLevel);
2101
2102 // See if we have parameters defined. They are optional so we go on
2103 // if not.
2104 if (variablesElement)
2105 {
2106 TiXmlElement *varElement = variablesElement->FirstChildElement("V");
2107
2108 // Sequential counter for the composite numbers.
2109 int nextVariableNumber = -1;
2110
2111 while (varElement)
2112 {
2113 stringstream tagcontent;
2114 tagcontent << *varElement;
2115
2116 /// All elements are of the form: "<V ID="#"> name = value
2117 /// </V>", with ? being the element type.
2118 nextVariableNumber++;
2119
2120 int i;
2121 int err = varElement->QueryIntAttribute("ID", &i);
2122 ASSERTL0(err == TIXML_SUCCESS,
2123 "Variables must have a unique ID number attribute "
2124 "in XML element: \n\t'" +
2125 tagcontent.str() + "'");
2126 ASSERTL0(i == nextVariableNumber,
2127 "ID numbers for variables must begin with zero and"
2128 " be sequential in XML element: \n\t'" +
2129 tagcontent.str() + "'");
2130
2131 TiXmlNode *varChild = varElement->FirstChild();
2132 // This is primarily to skip comments that may be present.
2133 // Comments appear as nodes just like elements. We are
2134 // specifically looking for text in the body of the
2135 // definition.
2136 while (varChild && varChild->Type() != TiXmlNode::TINYXML_TEXT)
2137 {
2138 varChild = varChild->NextSibling();
2139 }
2140
2141 ASSERTL0(varChild, "Unable to read variable definition body for "
2142 "variable with ID " +
2143 boost::lexical_cast<string>(i) +
2144 " in XML element: \n\t'" + tagcontent.str() +
2145 "'");
2146 std::string variableName = varChild->ToText()->ValueStr();
2147
2148 std::istringstream variableStrm(variableName);
2149 variableStrm >> variableName;
2150
2152 variableName) == m_variables.end(),
2153 "Variable with ID " + boost::lexical_cast<string>(i) +
2154 " in XML element \n\t'" + tagcontent.str() +
2155 "'\nhas already been defined.");
2156
2157 m_variables.push_back(variableName);
2158
2159 varElement = varElement->NextSiblingElement("V");
2160 }
2161
2162 ASSERTL0(nextVariableNumber > -1,
2163 "Number of variables must be greater than zero.");
2164 }
2165}
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:475

References ASSERTL0, Nektar::StdRegions::find(), 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 605 of file SessionReader.h.

608{
609 ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
610 CmdLineArg x;
611 x.shortName = pShortName;
612 x.description = pDescription;
613 x.isFlag = false;
614 GetCmdLineArgMap()[pName] = x;
615 return pName;
616}

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 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 = true;
630 GetCmdLineArgMap()[pName] = x;
631 return pName;
632}

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 594 of file SessionReader.h.

596{
597 std::string vName = boost::to_upper_copy(pName);
598 GetSolverInfoDefaults()[vName] = pValue;
599 return pValue;
600}

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.
static std::string lookupIds[]
Definition: GlobalLinSys.h:156
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 560 of file SessionReader.h.

563{
564 std::string vEnum = boost::to_upper_copy(pEnum);
565 auto x = GetSolverInfoEnums().find(vEnum);
566
567 if (x == GetSolverInfoEnums().end())
568 {
569 GetSolverInfoEnums()[vEnum] = EnumMap();
570 x = GetSolverInfoEnums().find(vEnum);
571 }
572
573 x->second[pString] = pEnumValue;
574 return pString;
575}
std::map< std::string, int > EnumMap
Definition: SessionReader.h:73

References GetSolverInfoEnums().

◆ SetParameter() [1/3]

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

Set an integer parameter.

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

919{
920 std::string vName = boost::to_upper_copy(pName);
921 m_parameters[vName] = pVar;
922}

References m_parameters.

◆ SetParameter() [2/3]

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

Set a double precision parameter.

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

937{
938 std::string vName = boost::to_upper_copy(pName);
939 m_parameters[vName] = pVar;
940}

References m_parameters.

◆ SetParameter() [3/3]

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

Set an size_t parameter.

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

928{
929 std::string vName = boost::to_upper_copy(pName);
930 m_parameters[vName] = pVar;
931}

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 970 of file BasicUtils/SessionReader.cpp.

972{
973 std::string vProperty = boost::to_upper_copy(pProperty);
974 auto iter = m_solverInfo.find(vProperty);
975
976 ASSERTL1(iter != m_solverInfo.end(),
977 "Unable to find requested property: " + pProperty);
978
979 iter->second = pValue;
980}

References ASSERTL1, and m_solverInfo.

◆ SetTag()

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

Sets a specified tag.

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

1397{
1398 std::string vName = boost::to_upper_copy(pName);
1399 m_tags[vName] = pValue;
1400}

References m_tags.

◆ SetUpdateOptFile()

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

Set bool to update optimisation file.

Definition at line 368 of file SessionReader.h.

369 {
370 m_updateOptFile = flag;
371 }

References m_updateOptFile.

◆ SetVariable()

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

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

1148{
1149 ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1150 m_variables[idx] = newname;
1151}

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().

◆ VerifySolverInfo()

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

Check values of solver info options are valid.

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

2541{
2542 for (auto &x : m_solverInfo)
2543 {
2544 std::string solverProperty = x.first;
2545 std::string solverValue = x.second;
2546
2547 auto propIt = GetSolverInfoEnums().find(solverProperty);
2548 if (propIt != GetSolverInfoEnums().end())
2549 {
2550 auto valIt = propIt->second.find(solverValue);
2551 ASSERTL0(valIt != propIt->second.end(), "Value '" + solverValue +
2552 "' is not valid for "
2553 "property '" +
2554 solverProperty + "'");
2555 }
2556 }
2557}

References ASSERTL0, GetSolverInfoEnums(), and m_solverInfo.

Referenced by InitSession().

Friends And Related Function Documentation

◆ MemoryManager< SessionReader >

friend class MemoryManager< SessionReader >
friend

Support creation through MemoryManager.

Definition at line 105 of file SessionReader.h.

Member Data Documentation

◆ m_backups

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

Backups.

Definition at line 420 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 390 of file SessionReader.h.

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

◆ m_filters

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

Filters map.

Definition at line 410 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 400 of file SessionReader.h.

◆ m_interpreter

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

Interpreter instance.

Definition at line 402 of file SessionReader.h.

Referenced by GetInterpreter(), ReadFunctions(), ReadParameters(), and SessionReader().

◆ 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 392 of file SessionReader.h.

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

◆ m_sharedFilesystem

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

Running on a shared filesystem.

Definition at line 418 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 408 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 412 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 422 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 394 of file SessionReader.h.

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