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

2466{
2467 // Parse solver info overrides
2468 if (m_cmdLineOptions.count("solverinfo"))
2469 {
2470 std::vector<std::string> solverInfoList =
2471 m_cmdLineOptions["solverinfo"].as<std::vector<std::string>>();
2472
2473 for (size_t i = 0; i < solverInfoList.size(); ++i)
2474 {
2475 std::string lhs, rhs;
2476
2477 try
2478 {
2479 ParseEquals(solverInfoList[i], lhs, rhs);
2480 }
2481 catch (...)
2482 {
2483 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2484 "option: " +
2485 solverInfoList[i]);
2486 }
2487
2488 std::string lhsUpper = boost::to_upper_copy(lhs);
2489 m_solverInfo[lhsUpper] = rhs;
2490 }
2491 }
2492
2493 if (m_cmdLineOptions.count("parameter"))
2494 {
2495 std::vector<std::string> parametersList =
2496 m_cmdLineOptions["parameter"].as<std::vector<std::string>>();
2497
2498 for (size_t i = 0; i < parametersList.size(); ++i)
2499 {
2500 std::string lhs, rhs;
2501
2502 try
2503 {
2504 ParseEquals(parametersList[i], lhs, rhs);
2505 }
2506 catch (...)
2507 {
2508 NEKERROR(ErrorUtil::efatal, "Parse error with command line "
2509 "option: " +
2510 parametersList[i]);
2511 }
2512
2513 std::string lhsUpper = boost::to_upper_copy(lhs);
2514
2515 try
2516 {
2517 m_parameters[lhsUpper] = boost::lexical_cast<NekDouble>(rhs);
2518 }
2519 catch (...)
2520 {
2521 NEKERROR(ErrorUtil::efatal, "Unable to convert string: " + rhs +
2522 "to double value.");
2523 }
2524 }
2525 }
2526}
#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 1609 of file BasicUtils/SessionReader.cpp.

1610{
1611 if (argc == 0)
1612 {
1613 m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1614 }
1615 else
1616 {
1617 string vCommModule("Serial");
1618 if (GetCommFactory().ModuleExists("ParallelMPI"))
1619 {
1620 vCommModule = "ParallelMPI";
1621 }
1622 if (m_cmdLineOptions.count("cwipi") &&
1623 GetCommFactory().ModuleExists("CWIPI"))
1624 {
1625 vCommModule = "CWIPI";
1626 }
1627
1628 m_comm = GetCommFactory().CreateInstance(vCommModule, argc, argv);
1629 }
1630}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
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 1414 of file BasicUtils/SessionReader.cpp.

1415{
1416 return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1417}

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

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

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

1163{
1164 std::string vName = boost::to_upper_copy(pName);
1165 return m_functions.find(vName) != m_functions.end();
1166}
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 1171 of file BasicUtils/SessionReader.cpp.

1174{
1175 std::string vName = boost::to_upper_copy(pName);
1176
1177 // Check function exists
1178 auto it1 = m_functions.find(vName);
1179 if (it1 != m_functions.end())
1180 {
1181 pair<std::string, int> key(pVariable, pDomain);
1182 pair<std::string, int> defkey("*", pDomain);
1183 bool varExists = it1->second.find(key) != it1->second.end() ||
1184 it1->second.find(defkey) != it1->second.end();
1185 return varExists;
1186 }
1187 return false;
1188}

References m_functions.

◆ DefinesGlobalSysSolnInfo()

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

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

1032{
1033 auto iter = GetGloSysSolnList().find(pVariable);
1034 if (iter == GetGloSysSolnList().end())
1035 {
1036 return false;
1037 }
1038
1039 std::string vProperty = boost::to_upper_copy(pProperty);
1040
1041 auto iter1 = iter->second.find(vProperty);
1042 if (iter1 == iter->second.end())
1043 {
1044 return false;
1045 }
1046
1047 return true;
1048}
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 785 of file BasicUtils/SessionReader.cpp.

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

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

936{
937 std::string vName = boost::to_upper_copy(pName);
938 auto infoIter = m_solverInfo.find(vName);
939 return (infoIter != m_solverInfo.end());
940}

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

1378{
1379 std::string vName = boost::to_upper_copy(pName);
1380 return m_tags.find(vName) != m_tags.end();
1381}

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

1075{
1076 return m_timeIntScheme.method != "";
1077}
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 777 of file BasicUtils/SessionReader.cpp.

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

References m_comm.

◆ GetBackups()

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

Returns the backups.

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

1155{
1156 return m_backups;
1157}

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

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

References m_cmdLineOptions.

◆ GetComm()

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

Returns the communication object.

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

762{
763 return m_comm;
764}

References m_comm.

◆ GetDocument()

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

Provides direct access to the TiXmlDocument object.

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

672{
673 ASSERTL1(m_xmlDoc, "XML Document not defined.");
674 return *m_xmlDoc;
675}
#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 699 of file BasicUtils/SessionReader.cpp.

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

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

746{
747 return m_filenames;
748}

References m_filenames.

◆ GetFilters()

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

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

1407{
1408 return m_filters;
1409}

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

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

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

1236{
1237 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1238 return GetFunction(pName, m_variables[pVar], pDomain);
1239}
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 1292 of file BasicUtils/SessionReader.cpp.

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

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

1332{
1333 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1334 return GetFunctionFilename(pName, m_variables[pVar], pDomain);
1335}
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 1340 of file BasicUtils/SessionReader.cpp.

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

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

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

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

1284{
1285 ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1286 return GetFunctionType(pName, m_variables[pVar], pDomain);
1287}
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 1091 of file BasicUtils/SessionReader.cpp.

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

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

1055{
1056 auto iter = GetGloSysSolnList().find(pVariable);
1057 ASSERTL0(iter != GetGloSysSolnList().end(),
1058 "Failed to find variable in GlobalSysSolnInfoList");
1059
1060 std::string vProperty = boost::to_upper_copy(pProperty);
1061 auto iter1 = iter->second.find(vProperty);
1062
1063 ASSERTL0(iter1 != iter->second.end(),
1064 "Failed to find property: " + vProperty +
1065 " in GlobalSysSolnInfoList");
1066
1067 return iter1->second;
1068}

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

2554{
2555 return m_interpreter;
2556}

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

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

References ASSERTL0, and m_parameters.

◆ GetSessionName()

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

Returns the session name of the loaded XML document.

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

754{
755 return m_sessionName;
756}

References m_sessionName.

◆ GetSharedFilesystem()

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

Returns if file system shared.

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

767{
768 return m_sharedFilesystem;
769}
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 945 of file BasicUtils/SessionReader.cpp.

947{
948 std::string vProperty = boost::to_upper_copy(pProperty);
949 auto iter = m_solverInfo.find(vProperty);
950
951 ASSERTL1(iter != m_solverInfo.end(),
952 "Unable to find requested property: " + pProperty);
953
954 return iter->second;
955}

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

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

1396{
1397 std::string vName = boost::to_upper_copy(pName);
1398 auto vTagIterator = m_tags.find(vName);
1399 ASSERTL0(vTagIterator != m_tags.end(), "Requested tag does not exist.");
1400 return vTagIterator->second;
1401}

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

1084{
1085 return m_timeIntScheme;
1086}

References m_timeIntScheme.

◆ GetTimeLevel()

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

Get time level (Parallel-in-Time)

Definition at line 373 of file SessionReader.h.

374 {
375 return m_timeLevel;
376 }

References m_timeLevel.

◆ GetUpdateOptFile()

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

Get bool to update optimisation file.

Definition at line 361 of file SessionReader.h.

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

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

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

1129{
1130 ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1131 return m_variables[idx];
1132}

References ASSERTL0, and m_variables.

◆ GetVariables()

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

Returns the names of all variables.

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

1147{
1148 return m_variables;
1149}

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

1425{
1426 if (Element && Element->FirstChildElement("TIMELEVEL"))
1427 {
1428 Element = Element->FirstChildElement("TIMELEVEL");
1429 std::string timeLevelStr;
1430 while (Element)
1431 {
1432 std::stringstream tagcontent;
1433 tagcontent << *Element;
1434 ASSERTL0(Element->Attribute("VALUE"),
1435 "Missing LEVEL attribute in solver info "
1436 "XML element: \n\t'" +
1437 tagcontent.str() + "'");
1438 timeLevelStr = Element->Attribute("VALUE");
1439 ASSERTL0(!timeLevelStr.empty(),
1440 "LEVEL attribute must be non-empty in XML "
1441 "element: \n\t'" +
1442 tagcontent.str() + "'");
1443 if (stoi(timeLevelStr) == timeLevel)
1444 {
1445 break;
1446 }
1447 Element = Element->NextSiblingElement("TIMELEVEL");
1448 }
1449 if (enableCheck)
1450 {
1451 ASSERTL0(stoi(timeLevelStr) == timeLevel,
1452 "TIMELEVEL value " + std::to_string(timeLevel) +
1453 " not found in solver info "
1454 "XML element: \n\t'");
1455 }
1456 }
1457}

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::MeshGraphHDF5::v_PartitionMesh(), and Nektar::SpatialDomains::MeshGraphXml::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 1462 of file BasicUtils/SessionReader.cpp.

1464{
1465 if (pFilename.size() > 3 &&
1466 pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1467 {
1468 ifstream file(pFilename.c_str(), ios_base::in | ios_base::binary);
1469 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1470 stringstream ss;
1471 io::filtering_streambuf<io::input> in;
1472 in.push(io::gzip_decompressor());
1473 in.push(file);
1474 try
1475 {
1476 io::copy(in, ss);
1477 ss >> (*pDoc);
1478 }
1479 catch (io::gzip_error &)
1480 {
1482 "Error: File '" + pFilename + "' is corrupt.");
1483 }
1484 }
1485 else if (pFilename.size() > 4 &&
1486 pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1487 {
1488 fs::path pdirname(pFilename);
1489 boost::format pad("P%1$07d.xml");
1490 pad % m_comm->GetSpaceComm()->GetRank();
1491 fs::path pRankFilename(pad.str());
1492 fs::path fullpath = pdirname / pRankFilename;
1493
1494 ifstream file(PortablePath(fullpath).c_str());
1495 ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1496 file >> (*pDoc);
1497 }
1498 else
1499 {
1500 ifstream file(pFilename.c_str());
1501 ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1502 file >> (*pDoc);
1503 }
1504}
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 813 of file BasicUtils/SessionReader.cpp.

814{
815 std::string vName = boost::to_upper_copy(pName);
816 auto paramIter = m_parameters.find(vName);
817 ASSERTL0(paramIter != m_parameters.end(),
818 "Required parameter '" + pName + "' not specified in session.");
819 NekDouble param = round(paramIter->second);
820 pVar = checked_cast<int>(param);
821}
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 826 of file BasicUtils/SessionReader.cpp.

828{
829 std::string vName = boost::to_upper_copy(pName);
830 auto paramIter = m_parameters.find(vName);
831 if (paramIter != m_parameters.end())
832 {
833 NekDouble param = round(paramIter->second);
834 pVar = checked_cast<int>(param);
835 }
836 else
837 {
838 pVar = pDefault;
839 }
840}

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

879{
880 std::string vName = boost::to_upper_copy(pName);
881 auto paramIter = m_parameters.find(vName);
882 ASSERTL0(paramIter != m_parameters.end(),
883 "Required parameter '" + pName + "' not specified in session.");
884 pVar = paramIter->second;
885}

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

892{
893 std::string vName = boost::to_upper_copy(pName);
894 auto paramIter = m_parameters.find(vName);
895 if (paramIter != m_parameters.end())
896 {
897 pVar = paramIter->second;
898 }
899 else
900 {
901 pVar = pDefault;
902 }
903}

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

846{
847 std::string vName = boost::to_upper_copy(pName);
848 auto paramIter = m_parameters.find(vName);
849 ASSERTL0(paramIter != m_parameters.end(),
850 "Required parameter '" + pName + "' not specified in session.");
851 NekDouble param = round(paramIter->second);
852 pVar = checked_cast<int>(param);
853}

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

860{
861 std::string vName = boost::to_upper_copy(pName);
862 auto paramIter = m_parameters.find(vName);
863 if (paramIter != m_parameters.end())
864 {
865 NekDouble param = round(paramIter->second);
866 pVar = checked_cast<int>(param);
867 }
868 else
869 {
870 pVar = pDefault;
871 }
872}

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

977{
978 std::string vName = boost::to_upper_copy(pName);
979 auto infoIter = m_solverInfo.find(vName);
980 if (infoIter != m_solverInfo.end())
981 {
982 pVar = infoIter->second;
983 }
984 else
985 {
986 pVar = pDefault;
987 }
988}

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

1014{
1015 if (DefinesSolverInfo(pName))
1016 {
1017 std::string vName = boost::to_upper_copy(pName);
1018 auto iter = m_solverInfo.find(vName);
1019 if (iter != m_solverInfo.end())
1020 {
1021 return boost::iequals(iter->second, pTrueVal);
1022 }
1023 }
1024 return false;
1025}

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

996{
997 std::string vName = boost::to_upper_copy(pName);
998 auto infoIter = m_solverInfo.find(vName);
999 if (infoIter != m_solverInfo.end())
1000 {
1001 pVar = boost::iequals(infoIter->second, pTrueVal);
1002 }
1003 else
1004 {
1005 pVar = pDefault;
1006 }
1007}

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

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

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

1511{
1512 ASSERTL0(pFilenames.size() > 0, "No filenames for merging.");
1513
1514 // Read the first document
1515 TiXmlDocument *vMainDoc = new TiXmlDocument;
1516 LoadDoc(pFilenames[0], vMainDoc);
1517
1518 TiXmlHandle vMainHandle(vMainDoc);
1519 TiXmlElement *vMainNektar =
1520 GetChildElementOrThrow(pFilenames[0], "NEKTAR", vMainHandle);
1521
1522 // Read all subsequent XML documents.
1523 // For each element within the NEKTAR tag, use it to replace the
1524 // version already present in the loaded XML data.
1525 for (int i = 1; i < pFilenames.size(); ++i)
1526 {
1527 if ((pFilenames[i].compare(pFilenames[i].size() - 3, 3, "xml") == 0) ||
1528 (pFilenames[i].compare(pFilenames[i].size() - 6, 6, "xml.gz") ==
1529 0) ||
1530 (pFilenames[i].compare(pFilenames[i].size() - 3, 3, "opt") == 0))
1531 {
1532 TiXmlDocument *vTempDoc = new TiXmlDocument;
1533 LoadDoc(pFilenames[i], vTempDoc);
1534
1535 TiXmlHandle docHandle(vTempDoc);
1536 TiXmlElement *vTempNektar =
1537 GetChildElementOrThrow(pFilenames[i], "NEKTAR", docHandle);
1538 TiXmlElement *p = vTempNektar->FirstChildElement();
1539
1540 while (p)
1541 {
1542 TiXmlElement *vMainEntry =
1543 vMainNektar->FirstChildElement(p->Value());
1544
1545 // First check if the new item is in fact blank
1546 // replace if it is a COLLECTIONS section however.
1547 if (!p->FirstChild() && vMainEntry &&
1548 !boost::iequals(p->Value(), "COLLECTIONS"))
1549 {
1550 std::string warningmsg =
1551 "File " + pFilenames[i] + " contains " +
1552 "an empty XML element " + std::string(p->Value()) +
1553 " which will be ignored.";
1554 NEKERROR(ErrorUtil::ewarning, warningmsg.c_str());
1555 }
1556 else
1557 {
1558 if (vMainEntry)
1559 {
1560 vMainNektar->RemoveChild(vMainEntry);
1561 }
1562 TiXmlElement *q = new TiXmlElement(*p);
1563 vMainNektar->LinkEndChild(q);
1564 }
1565 p = p->NextSiblingElement();
1566 }
1567 delete vTempDoc;
1568 }
1569 }
1570 return vMainDoc;
1571}
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 ("verbose,v", "be verbose")
449 ("version,V", "print version information")
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 ("npx", po::value<int>(),
456 "number of procs in X-dir")
457 ("npy", po::value<int>(),
458 "number of procs in Y-dir")
459 ("npz", po::value<int>(),
460 "number of procs in Z-dir")
461 ("nsz", po::value<int>(),
462 "number of slices in Z-dir")
463 ("npt", po::value<int>(),
464 "number of procs in T-dir (parareal)")
465 ("part-only", po::value<int>(),
466 "only partition mesh into N partitions.")
467 ("part-only-overlapping", po::value<int>(),
468 "only partition mesh into N overlapping partitions.")
469 ("part-info", "output partition information")
470 ("force-output,f","disables backups files and forces output to be "
471 "written without any checks")
472 ("write-opt-file","write an optimisation file")
473 ("use-opt-file", po::value<std::string>(),
474 "use an optimisation file");
475 // clang-format on
476
477#ifdef NEKTAR_USE_CWIPI
478 desc.add_options()("cwipi", po::value<std::string>(), "set CWIPI name");
479#endif
480
481 for (auto &cmdIt : GetCmdLineArgMap())
482 {
483 std::string names = cmdIt.first;
484 if (cmdIt.second.shortName != "")
485 {
486 names += "," + cmdIt.second.shortName;
487 }
488 if (cmdIt.second.isFlag)
489 {
490 desc.add_options()(names.c_str(), cmdIt.second.description.c_str());
491 }
492 else
493 {
494 desc.add_options()(names.c_str(), po::value<std::string>(),
495 cmdIt.second.description.c_str());
496 }
497 }
498
499 // Deprecated options: introduced in 5.4.0 to homogenise command-line
500 // options to use '-' instead of camelCase or no spaces.
501 std::map<std::string, std::string> deprecated = {
502 {"forceoutput", "force-output"},
503 {"writeoptfile", "write-opt-file"},
504 {"useoptfile", "use-opt-file"}};
505
506 for (auto &d : deprecated)
507 {
508 std::string description = "Deprecated: use --" + d.second;
509 dep.add_options()(d.first.c_str(), description.c_str());
510 }
511
512 // List hidden options (e.g. session file arguments are not actually
513 // specified using the input-file option by the user).
514 po::options_description hidden("Hidden options");
515
516 hidden.add_options()("input-file", po::value<vector<string>>(),
517 "input filename");
518
519 // Combine all options for the parser
520 po::options_description all("All options");
521 all.add(desc).add(dep).add(hidden);
522
523 // Session file is a positional option
524 po::positional_options_description p;
525 p.add("input-file", -1);
526
527 // Parse the command-line options
528 po::parsed_options parsed = po::command_line_parser(argc, argv)
529 .options(all)
530 .positional(p)
531 .allow_unregistered()
532 .run();
533
534 // Extract known options to map and update
535 po::store(parsed, m_cmdLineOptions);
536 po::notify(m_cmdLineOptions);
537
538 // Help message
539 if (m_cmdLineOptions.count("help"))
540 {
541 cout << desc;
542 exit(0);
543 }
544
545 // Version information
546 if (m_cmdLineOptions.count("version"))
547 {
548 cout << "Nektar++ version " << NEKTAR_VERSION;
549
550 if (NekConstants::kGitSha1 != "GITDIR-NOTFOUND")
551 {
553 string branch(NekConstants::kGitBranch);
554 boost::replace_all(branch, "refs/heads/", "");
555
556 cout << " (git changeset " << sha1.substr(0, 8) << ", ";
557
558 if (branch == "")
559 {
560 cout << "detached head";
561 }
562 else
563 {
564 cout << "head " << branch;
565 }
566
567 cout << ")";
568 }
569
570 cout << endl;
571 exit(0);
572 }
573
574 // Deal with deprecated options.
575 for (auto &d : deprecated)
576 {
577 if (m_cmdLineOptions.count(d.first))
578 {
579 std::cerr << "Warning: --" << d.first << " deprecated: use --"
580 << d.second << std::endl;
581 m_cmdLineOptions.emplace(
582 d.second, po::variable_value(m_cmdLineOptions[d.first]));
583 }
584 }
585
586 // Enable verbose mode
587 if (m_cmdLineOptions.count("verbose"))
588 {
589 m_verbose = true;
590 }
591 else
592 {
593 m_verbose = false;
594 }
595
596 // Disable backups
597 if (m_cmdLineOptions.count("force-output"))
598 {
599 m_backups = false;
600 }
601 else
602 {
603 m_backups = true;
604 }
605
606 // Enable update optimisation file
607 if (m_cmdLineOptions.count("write-opt-file"))
608 {
609 m_updateOptFile = true;
610 }
611 else
612 {
613 m_updateOptFile = false;
614 }
615
616 // Print a warning for unknown options
617 for (auto &x : parsed.options)
618 {
619 if (x.unregistered)
620 {
621 cout << "Warning: Unknown option: " << x.string_key << endl;
622 }
623 }
624
625 // Return the vector of filename(s) given as positional options
626 if (m_cmdLineOptions.count("input-file"))
627 {
628 return m_cmdLineOptions["input-file"].as<std::vector<std::string>>();
629 }
630 else
631 {
632 return std::vector<std::string>();
633 }
634}
#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 1576 of file BasicUtils/SessionReader.cpp.

1577{
1578 // Check we actually have a document loaded.
1579 ASSERTL0(m_xmlDoc, "No XML document loaded.");
1580
1581 TiXmlHandle docHandle(m_xmlDoc);
1582 TiXmlElement *e;
1583
1584 // Look for all data in CONDITIONS block.
1585 e = docHandle.FirstChildElement("NEKTAR")
1586 .FirstChildElement("CONDITIONS")
1587 .Element();
1588
1589 // Read the various sections of the CONDITIONS block
1590 ReadParameters(e);
1591 ReadSolverInfo(e);
1594 ReadVariables(e);
1595 ReadFunctions(e);
1596
1597 // Look for all data in FILTERS block.
1598 e = docHandle.FirstChildElement("NEKTAR")
1599 .FirstChildElement("FILTERS")
1600 .Element();
1601
1602 // Read the various sections of the FILTERS block
1603 ReadFilters(e);
1604}
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 2433 of file BasicUtils/SessionReader.cpp.

2435{
2436 /// Pull out lhs and rhs and eliminate any spaces.
2437 size_t beg = line.find_first_not_of(" ");
2438 size_t end = line.find_first_of("=");
2439 // Check for no parameter name
2440 if (beg == end)
2441 {
2442 throw 1;
2443 }
2444 // Check for no parameter value
2445 if (end != line.find_last_of("="))
2446 {
2447 throw 1;
2448 }
2449 // Check for no equals sign
2450 if (end == std::string::npos)
2451 {
2452 throw 1;
2453 }
2454
2455 lhs = line.substr(line.find_first_not_of(" "), end - beg);
2456 lhs = lhs.substr(0, lhs.find_last_not_of(" ") + 1);
2457 rhs = line.substr(line.find_last_of("=") + 1);
2458 rhs = rhs.substr(rhs.find_first_not_of(" "));
2459 rhs = rhs.substr(0, rhs.find_last_not_of(" ") + 1);
2460}

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

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

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

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

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

2387{
2388 if (!filters)
2389 {
2390 return;
2391 }
2392
2393 m_filters.clear();
2394
2395 TiXmlElement *filter = filters->FirstChildElement("FILTER");
2396
2397 while (filter)
2398 {
2399 ASSERTL0(filter->Attribute("TYPE"),
2400 "Missing attribute 'TYPE' for filter.");
2401 std::string typeStr = filter->Attribute("TYPE");
2402
2403 std::map<std::string, std::string> vParams;
2404
2405 TiXmlElement *element = filter;
2406 GetXMLElementTimeLevel(element, m_timeLevel, false);
2407 TiXmlElement *param = element->FirstChildElement("PARAM");
2408 while (param)
2409 {
2410 ASSERTL0(param->Attribute("NAME"),
2411 "Missing attribute 'NAME' for parameter in filter " +
2412 typeStr + "'.");
2413 std::string nameStr = param->Attribute("NAME");
2414
2415 ASSERTL0(param->GetText(), "Empty value string for param.");
2416 std::string valueStr = param->GetText();
2417
2418 vParams[nameStr] = valueStr;
2419
2420 param = param->NextSiblingElement("PARAM");
2421 }
2422
2423 m_filters.push_back(
2424 std::pair<std::string, FilterParams>(typeStr, vParams));
2425
2426 filter = filter->NextSiblingElement("FILTER");
2427 }
2428}
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 2161 of file BasicUtils/SessionReader.cpp.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

562{
563 std::string vEnum = boost::to_upper_copy(pEnum);
564 auto x = GetSolverInfoEnums().find(vEnum);
565
566 if (x == GetSolverInfoEnums().end())
567 {
568 GetSolverInfoEnums()[vEnum] = EnumMap();
569 x = GetSolverInfoEnums().find(vEnum);
570 }
571
572 x->second[pString] = pEnumValue;
573 return pString;
574}
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 908 of file BasicUtils/SessionReader.cpp.

909{
910 std::string vName = boost::to_upper_copy(pName);
911 m_parameters[vName] = pVar;
912}

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

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

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

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

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

962{
963 std::string vProperty = boost::to_upper_copy(pProperty);
964 auto iter = m_solverInfo.find(vProperty);
965
966 ASSERTL1(iter != m_solverInfo.end(),
967 "Unable to find requested property: " + pProperty);
968
969 iter->second = pValue;
970}

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

1387{
1388 std::string vName = boost::to_upper_copy(pName);
1389 m_tags[vName] = pValue;
1390}

References m_tags.

◆ SetUpdateOptFile()

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

Set bool to update optimisation file.

Definition at line 367 of file SessionReader.h.

368 {
369 m_updateOptFile = flag;
370 }

References m_updateOptFile.

◆ SetVariable()

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

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

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

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

2532{
2533 for (auto &x : m_solverInfo)
2534 {
2535 std::string solverProperty = x.first;
2536 std::string solverValue = x.second;
2537
2538 auto propIt = GetSolverInfoEnums().find(solverProperty);
2539 if (propIt != GetSolverInfoEnums().end())
2540 {
2541 auto valIt = propIt->second.find(solverValue);
2542 ASSERTL0(valIt != propIt->second.end(), "Value '" + solverValue +
2543 "' is not valid for "
2544 "property '" +
2545 solverProperty + "'");
2546 }
2547 }
2548}

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

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

◆ m_filters

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

Filters map.

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

◆ m_interpreter

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

Interpreter instance.

Definition at line 401 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 391 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 417 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 407 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 411 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 421 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 393 of file SessionReader.h.

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