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>

Inheritance diagram for Nektar::LibUtilities::SessionReader:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LibUtilities::SessionReader:
Collaboration graph
[legend]

Public Member Functions

 SessionReader (int argc, char *argv[], const std::vector< std::string > &pFilenames, const CommSharedPtr &pComm)
 
 ~SessionReader ()
 Destructor. More...
 
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...
 
const std::string GetSessionNameRank () const
 Returns the session name with process rank. More...
 
CommSharedPtrGetComm ()
 Returns the communication object. 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, int &var, const int &def) const
 Check for and load an integer 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, 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 DefinesGlobalSysSolnInfo (const std::string &variable, const std::string &property) const
 
const std::string & GetGlobalSysSolnInfo (const std::string &variable, const std::string &property) const
 
bool DefinesGeometricInfo (const std::string &name) const
 Checks if a geometric info property is defined. More...
 
void LoadGeometricInfo (const std::string &name, std::string &var, const std::string &def="") const
 Checks for and load a geometric info string property. More...
 
void LoadGeometricInfo (const std::string &name, bool &var, const bool &def=false) const
 Checks for and loads a geometric info boolean property. More...
 
void LoadGeometricInfo (const std::string &name, NekDouble &var, const NekDouble &def=0.0) const
 Checks for and loads a geometric info double-precision property. More...
 
void MatchGeometricInfo (const std::string &name, const std::string &trueval, bool &var, const bool &def=false) const
 Check if the value of a geometric info string property matches. More...
 
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...
 
AnalyticExpressionEvaluatorGetExpressionEvaluator ()
 Returns the instance of AnalyticExpressionEvaluator 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...
 
void SubstituteExpressions (std::string &expr)
 Substitutes expressions defined in the XML document. More...
 
CompositeOrdering GetCompositeOrdering () const
 
BndRegionOrdering GetBndRegionOrdering () const
 
void SetUpXmlDoc ()
 

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

Private Member Functions

 SessionReader (int argc, char *argv[])
 Main constructor. More...
 
void InitSession ()
 
SessionReaderSharedPtr GetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
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 PartitionMesh ()
 Partitions the mesh when running in parallel. 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 ReadExpressions (TiXmlElement *conditions)
 Reads the EXPRESSIONS 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 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...
 
ExpressionMap m_expressions
 Expressions. More...
 
AnalyticExpressionEvaluator m_exprEvaluator
 Analytic expression evaluator instance. More...
 
FunctionMap m_functions
 Functions. More...
 
VariableList m_variables
 Variables. More...
 
TagMap m_tags
 Custom tags. More...
 
FilterMap m_filters
 Filters map. More...
 
bool m_verbose
 Be verbose. More...
 
CompositeOrdering m_compOrder
 Map of original composite ordering for parallel periodic bcs. More...
 
BndRegionOrdering m_bndRegOrder
 Map of original boundary region ordering for parallel periodic bcs. 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:

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

Constructor & Destructor Documentation

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

Definition at line 210 of file SessionReader.cpp.

References ASSERTL0.

215  {
216  ASSERTL0(pFilenames.size() > 0, "No filenames specified.");
217 
218  ParseCommandLineArguments(argc, argv);
219  m_xmlDoc = 0;
220  m_filenames = pFilenames;
221 
223 
224  // Create communicator
225  if (!pComm.get())
226  {
227  CreateComm(argc, argv);
228  }
229  else
230  {
231  m_comm = pComm;
232 
233  if (m_comm->GetSize() > 1)
234  {
235  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
236  "IterativeStaticCond";
237  }
238  }
239 
240  // If running in parallel change the default global sys solution
241  // type.
242  if (m_comm->GetSize() > 1)
243  {
244  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
245  "IterativeStaticCond";
246  }
247  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
std::vector< std::string > ParseCommandLineArguments(int argc, char *argv[])
Parse the program arguments and fill m_cmdLineOptions.
CommSharedPtr m_comm
Communication object.
std::vector< std::string > m_filenames
Filenames.
void CreateComm(int &argc, char *argv[])
Loads the given XML document and instantiates an appropriate communication object.
std::string ParseSessionName(std::vector< std::string > &filenames)
Parse the session name.
static SolverInfoMap & GetSolverInfoDefaults()
Default solver info options.
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
std::string m_sessionName
Session name of the loaded XML document (filename minus ext).
Nektar::LibUtilities::SessionReader::~SessionReader ( )

Destructor.

Definition at line 253 of file SessionReader.cpp.

254  {
255  delete m_xmlDoc;
256  }
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
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 185 of file SessionReader.cpp.

References ASSERTL0.

186  {
187  m_xmlDoc = 0;
189 
190  ASSERTL0(m_filenames.size() > 0, "No session file(s) given.");
191 
193 
194  // Create communicator
195  CreateComm(argc, argv);
196 
197  // If running in parallel change the default global sys solution
198  // type.
199  if (m_comm->GetSize() > 1)
200  {
201  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
202  "IterativeStaticCond";
203  }
204  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
std::vector< std::string > ParseCommandLineArguments(int argc, char *argv[])
Parse the program arguments and fill m_cmdLineOptions.
CommSharedPtr m_comm
Communication object.
std::vector< std::string > m_filenames
Filenames.
void CreateComm(int &argc, char *argv[])
Loads the given XML document and instantiates an appropriate communication object.
std::string ParseSessionName(std::vector< std::string > &filenames)
Parse the session name.
static SolverInfoMap & GetSolverInfoDefaults()
Default solver info options.
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
std::string m_sessionName
Session name of the loaded XML document (filename minus ext).

Member Function Documentation

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

Enforce parameters from command line arguments.

Definition at line 2570 of file SessionReader.cpp.

References ASSERTL0, and Nektar::lhs.

2571  {
2572  // Parse solver info overrides
2573  if (m_cmdLineOptions.count("solverinfo"))
2574  {
2575  std::vector<std::string> solverInfoList =
2576  m_cmdLineOptions["solverinfo"].as<
2577  std::vector<std::string> >();
2578 
2579  for (int i = 0; i < solverInfoList.size(); ++i)
2580  {
2581  std::string lhs, rhs;
2582 
2583  try
2584  {
2585  ParseEquals(solverInfoList[i], lhs, rhs);
2586  }
2587  catch (...)
2588  {
2589  ASSERTL0(false, "Parse error with command line "
2590  "option: "+solverInfoList[i]);
2591  }
2592 
2593  std::string lhsUpper = boost::to_upper_copy(lhs);
2594  m_solverInfo[lhsUpper] = rhs;
2595  }
2596  }
2597 
2598  if (m_cmdLineOptions.count("parameter"))
2599  {
2600  std::vector<std::string> parametersList =
2601  m_cmdLineOptions["parameter"].as<
2602  std::vector<std::string> >();
2603 
2604  for (int i = 0; i < parametersList.size(); ++i)
2605  {
2606  std::string lhs, rhs;
2607 
2608  try
2609  {
2610  ParseEquals(parametersList[i], lhs, rhs);
2611  }
2612  catch (...)
2613  {
2614  ASSERTL0(false, "Parse error with command line "
2615  "option: "+parametersList[i]);
2616  }
2617 
2618  std::string lhsUpper = boost::to_upper_copy(lhs);
2619 
2620  try
2621  {
2622  m_parameters[lhsUpper] =
2623  boost::lexical_cast<NekDouble>(rhs);
2624  }
2625  catch (...)
2626  {
2627  ASSERTL0(false, "Unable to convert string: "+rhs+
2628  "to double value.");
2629  }
2630  }
2631  }
2632  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
SolverInfoMap m_solverInfo
Solver information properties.
StandardMatrixTag & lhs
void ParseEquals(const std::string &line, std::string &lhs, std::string &rhs)
Parse a string in the form lhs = rhs.
boost::program_options::variables_map m_cmdLineOptions
double NekDouble
void Nektar::LibUtilities::SessionReader::CreateComm ( int &  argc,
char *  argv[] 
)
private

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

Definition at line 1509 of file SessionReader.cpp.

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

1512  {
1513  if (argc == 0)
1514  {
1515  m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1516  }
1517  else
1518  {
1519  string vCommModule("Serial");
1520  if (GetCommFactory().ModuleExists("ParallelMPI"))
1521  {
1522  vCommModule = "ParallelMPI";
1523  }
1524 
1525  m_comm = GetCommFactory().CreateInstance(vCommModule,argc,argv);
1526  }
1527  }
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
CommSharedPtr m_comm
Communication object.
CommFactory & GetCommFactory()
Definition: Comm.cpp:64
static SessionReaderSharedPtr Nektar::LibUtilities::SessionReader::CreateInstance ( int  argc,
char *  argv[] 
)
inlinestatic

Creates an instance of the SessionReader class.

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

Definition at line 140 of file SessionReader.h.

Referenced by Diffusion::Diffusion(), main(), Nektar::Utilities::InputNekpp::Process(), Nektar::Utilities::InputXml::Process(), Nektar::Utilities::ProcessDisplacement::Process(), Nektar::Utilities::OutputNekpp::Process(), Nektar::Utilities::ProcessInterpPoints::Process(), Nektar::Utilities::ProcessInterpField::Process(), Nektar::SolverUtils::Driver::v_InitObject(), and Nektar::VortexWaveInteraction::VortexWaveInteraction().

142  {
143  SessionReaderSharedPtr p = MemoryManager<
144  LibUtilities::SessionReader>::AllocateSharedPtr(argc, argv);
145  p->InitSession();
146  return p;
147  }
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
static SessionReaderSharedPtr Nektar::LibUtilities::SessionReader::CreateInstance ( int  argc,
char *  argv[],
std::vector< std::string > &  pFilenames,
const CommSharedPtr pComm = CommSharedPtr() 
)
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 159 of file SessionReader.h.

164  {
165  SessionReaderSharedPtr p = MemoryManager<
166  LibUtilities::SessionReader>
167  ::AllocateSharedPtr(argc, argv, pFilenames, pComm);
168  p->InitSession();
169  return p;
170  }
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
bool Nektar::LibUtilities::SessionReader::DefinesCmdLineArgument ( const std::string &  pName) const

Checks if a specified cmdline argument has been given.

Definition at line 1345 of file SessionReader.cpp.

1347  {
1348  return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1349  }
boost::program_options::variables_map m_cmdLineOptions
bool Nektar::LibUtilities::SessionReader::DefinesElement ( const std::string &  pPath) const

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

Definition at line 555 of file SessionReader.cpp.

References ASSERTL0.

556  {
557  std::string vPath = boost::to_upper_copy(pPath);
558  std::vector<std::string> st;
559  boost::split(st, vPath, boost::is_any_of("\\/ "));
560  ASSERTL0(st.size() > 0, "No path given in XML element request.");
561 
562  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
563  ASSERTL0(vReturn, std::string("Cannot find element '")
564  + st[0] + std::string("'."));
565  for (int i = 1; i < st.size(); ++i)
566  {
567  vReturn = vReturn->FirstChildElement(st[i].c_str());
568  if (!vReturn) return false;
569  }
570  return true;
571  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
bool Nektar::LibUtilities::SessionReader::DefinesFunction ( const std::string &  name) const

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

Definition at line 1041 of file SessionReader.cpp.

1042  {
1043  FunctionMap::const_iterator it1;
1044  std::string vName = boost::to_upper_copy(pName);
1045 
1046  if ((it1 = m_functions.find(vName)) != m_functions.end())
1047  {
1048  return true;
1049  }
1050  return false;
1051  }
FunctionMap m_functions
Functions.
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 1057 of file SessionReader.cpp.

1061  {
1062  FunctionMap::const_iterator it1;
1063  FunctionVariableMap::const_iterator it2;
1064  std::string vName = boost::to_upper_copy(pName);
1065 
1066  // Check function exists
1067  if ((it1 = m_functions.find(vName)) != m_functions.end())
1068  {
1069  pair<std::string, int> key(pVariable,pDomain);
1070  pair<std::string, int> defkey("*",pDomain);
1071  bool varExists =
1072  (it2 = it1->second.find(key)) != it1->second.end() ||
1073  (it2 = it1->second.find(defkey)) != it1->second.end();
1074  return varExists;
1075  }
1076  return false;
1077  }
FunctionMap m_functions
Functions.
bool Nektar::LibUtilities::SessionReader::DefinesGeometricInfo ( const std::string &  name) const

Checks if a geometric info property is defined.

Definition at line 906 of file SessionReader.cpp.

907  {
908  std::string vName = boost::to_upper_copy(pName);
909  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
910  return (iter != m_geometricInfo.end());
911  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
bool Nektar::LibUtilities::SessionReader::DefinesGlobalSysSolnInfo ( const std::string &  variable,
const std::string &  property 
) const

Definition at line 860 of file SessionReader.cpp.

862  {
863 
864  GloSysSolnInfoList::const_iterator iter =
865  GetGloSysSolnList().find(pVariable);
866  if(iter == GetGloSysSolnList().end())
867  {
868  return false;
869  }
870 
871  std::string vProperty = boost::to_upper_copy(pProperty);
872 
873  GloSysInfoMap::const_iterator iter1 = iter->second.find(vProperty);
874  if(iter1 == iter->second.end())
875  {
876  return false;
877  }
878 
879  return true;
880  }
static GloSysSolnInfoList & GetGloSysSolnList()
GlobalSysSoln Info map.
bool Nektar::LibUtilities::SessionReader::DefinesParameter ( const std::string &  name) const

Checks if a parameter is specified in the XML document.

Definition at line 633 of file SessionReader.cpp.

634  {
635  std::string vName = boost::to_upper_copy(pName);
636  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
637  return (paramIter != m_parameters.end());
638  }
ParameterMap m_parameters
Parameters.
bool Nektar::LibUtilities::SessionReader::DefinesSolverInfo ( const std::string &  name) const

Checks if a solver info property is specified.

Definition at line 756 of file SessionReader.cpp.

Referenced by GetSolverInfoAsEnum().

757  {
758  std::string vName = boost::to_upper_copy(pName);
759  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
760  return (infoIter != m_solverInfo.end());
761  }
SolverInfoMap m_solverInfo
Solver information properties.
bool Nektar::LibUtilities::SessionReader::DefinesTag ( const std::string &  pName) const

Checks if a specified tag is defined.

Definition at line 1300 of file SessionReader.cpp.

1301  {
1302  std::string vName = boost::to_upper_copy(pName);
1303  TagMap::const_iterator vTagIterator = m_tags.find(vName);
1304  return (vTagIterator != m_tags.end());
1305  }
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 624 of file SessionReader.cpp.

625  {
626  m_comm->Finalise();
627  }
CommSharedPtr m_comm
Communication object.
BndRegionOrdering Nektar::LibUtilities::SessionReader::GetBndRegionOrdering ( ) const

Definition at line 1370 of file SessionReader.cpp.

1371  {
1372  return m_bndRegOrder;
1373  }
BndRegionOrdering m_bndRegOrder
Map of original boundary region ordering for parallel periodic bcs.
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 170 of file SessionReader.cpp.

Referenced by RegisterCmdLineArgument(), and RegisterCmdLineFlag().

171  {
172  static CmdLineArgMap cmdLineArguments;
173  return cmdLineArguments;
174  }
std::map< std::string, CmdLineArg > CmdLineArgMap
Definition: SessionReader.h:75
template<typename T >
T Nektar::LibUtilities::SessionReader::GetCmdLineArgument ( const std::string &  pName) const
inline

Retrieves a command-line argument value.

Definition at line 404 of file SessionReader.h.

References m_cmdLineOptions.

406  {
407  return m_cmdLineOptions.find(pName)->second.as<T>();
408  }
boost::program_options::variables_map m_cmdLineOptions
CommSharedPtr & Nektar::LibUtilities::SessionReader::GetComm ( )

Returns the communication object.

Definition at line 612 of file SessionReader.cpp.

613  {
614  return m_comm;
615  }
CommSharedPtr m_comm
Communication object.
CompositeOrdering Nektar::LibUtilities::SessionReader::GetCompositeOrdering ( ) const

Definition at line 1365 of file SessionReader.cpp.

1366  {
1367  return m_compOrder;
1368  }
CompositeOrdering m_compOrder
Map of original composite ordering for parallel periodic bcs.
TiXmlDocument & Nektar::LibUtilities::SessionReader::GetDocument ( )

Provides direct access to the TiXmlDocument object.

Definition at line 503 of file SessionReader.cpp.

References ASSERTL1.

504  {
505  ASSERTL1(m_xmlDoc, "XML Document not defined.");
506  return *m_xmlDoc;
507  }
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:165
TiXmlElement * Nektar::LibUtilities::SessionReader::GetElement ( const std::string &  pPath)

Provides direct access to the TiXmlElement specified.

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

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

the PARAMETERS element would be retrieved by requesting the path:

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

Definition at line 532 of file SessionReader.cpp.

References ASSERTL0.

533  {
534  std::string vPath = boost::to_upper_copy(pPath);
535  std::vector<std::string> st;
536  boost::split(st, vPath, boost::is_any_of("\\/ "));
537  ASSERTL0(st.size() > 0, "No path given in XML element request.");
538 
539  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
540  ASSERTL0(vReturn, std::string("Cannot find element '")
541  + st[0] + std::string("'."));
542  for (int i = 1; i < st.size(); ++i)
543  {
544  vReturn = vReturn->FirstChildElement(st[i].c_str());
545  ASSERTL0(vReturn, std::string("Cannot find element '")
546  + st[i] + std::string("'."));
547  }
548  return vReturn;
549  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
AnalyticExpressionEvaluator & Nektar::LibUtilities::SessionReader::GetExpressionEvaluator ( )

Returns the instance of AnalyticExpressionEvaluator specific to this session.

Definition at line 1291 of file SessionReader.cpp.

1292  {
1293  return m_exprEvaluator;
1294  }
AnalyticExpressionEvaluator m_exprEvaluator
Analytic expression evaluator instance.
const std::vector< std::string > & Nektar::LibUtilities::SessionReader::GetFilenames ( ) const

Returns the filename of the loaded XML document.

Definition at line 577 of file SessionReader.cpp.

578  {
579  return m_filenames;
580  }
std::vector< std::string > m_filenames
Filenames.
const FilterMap & Nektar::LibUtilities::SessionReader::GetFilters ( ) const

Definition at line 1336 of file SessionReader.cpp.

1337  {
1338  return m_filters;
1339  }
FilterMap m_filters
Filters map.
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 1083 of file SessionReader.cpp.

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

1087  {
1088  FunctionMap::const_iterator it1;
1089  FunctionVariableMap::const_iterator it2, it3;
1090  std::string vName = boost::to_upper_copy(pName);
1091 
1092  ASSERTL0((it1 = m_functions.find(vName)) != m_functions.end(),
1093  std::string("No such function '") + pName
1094  + std::string("' has been defined in the session file."));
1095 
1096  // Check for specific and wildcard definitions
1097  pair<std::string,int> key(pVariable,pDomain);
1098  pair<std::string,int> defkey("*",pDomain);
1099  bool specific = (it2 = it1->second.find(key)) !=
1100  it1->second.end();
1101  bool wildcard = (it3 = it1->second.find(defkey)) !=
1102  it1->second.end();
1103 
1104  // Check function is defined somewhere
1105  ASSERTL0(specific || wildcard,
1106  "No such variable " + pVariable
1107  + " in domain " + boost::lexical_cast<string>(pDomain)
1108  + " defined for function " + pName
1109  + " in session file.");
1110 
1111  // If not specific, must be wildcard
1112  if (!specific)
1113  {
1114  it2 = it3;
1115  }
1116 
1117  ASSERTL0((it2->second.m_type == eFunctionTypeExpression),
1118  std::string("Function is defined by a file."));
1119  return it2->second.m_expression;
1120  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
FunctionMap m_functions
Functions.
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 1126 of file SessionReader.cpp.

References ASSERTL0.

1130  {
1131  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1132  return GetFunction(pName, m_variables[pVar],pDomain);
1133  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
EquationSharedPtr GetFunction(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns an EquationSharedPtr to a given function variable.
VariableList m_variables
Variables.
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 1194 of file SessionReader.cpp.

References ASSERTL0.

1198  {
1199  FunctionMap::const_iterator it1;
1200  FunctionVariableMap::const_iterator it2, it3;
1201  std::string vName = boost::to_upper_copy(pName);
1202 
1203  it1 = m_functions.find(vName);
1204  ASSERTL0 (it1 != m_functions.end(),
1205  std::string("Function '") + pName
1206  + std::string("' not found."));
1207 
1208  // Check for specific and wildcard definitions
1209  pair<std::string,int> key(pVariable,pDomain);
1210  pair<std::string,int> defkey("*",pDomain);
1211  bool specific = (it2 = it1->second.find(key)) !=
1212  it1->second.end();
1213  bool wildcard = (it3 = it1->second.find(defkey)) !=
1214  it1->second.end();
1215 
1216  // Check function is defined somewhere
1217  ASSERTL0(specific || wildcard,
1218  "No such variable " + pVariable
1219  + " in domain " + boost::lexical_cast<string>(pDomain)
1220  + " defined for function " + pName
1221  + " in session file.");
1222 
1223  // If not specific, must be wildcard
1224  if (!specific)
1225  {
1226  it2 = it3;
1227  }
1228 
1229  return it2->second.m_filename;
1230  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
FunctionMap m_functions
Functions.
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 1236 of file SessionReader.cpp.

References ASSERTL0.

1240  {
1241  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1242  return GetFunctionFilename(pName, m_variables[pVar],pDomain);
1243  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
VariableList m_variables
Variables.
std::string GetFunctionFilename(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns the filename to be loaded for a given variable.
std::string 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 1249 of file SessionReader.cpp.

References ASSERTL0.

1253  {
1254  FunctionMap::const_iterator it1;
1255  FunctionVariableMap::const_iterator it2, it3;
1256  std::string vName = boost::to_upper_copy(pName);
1257 
1258  it1 = m_functions.find(vName);
1259  ASSERTL0 (it1 != m_functions.end(),
1260  std::string("Function '") + pName
1261  + std::string("' not found."));
1262 
1263  // Check for specific and wildcard definitions
1264  pair<std::string,int> key(pVariable,pDomain);
1265  pair<std::string,int> defkey("*",pDomain);
1266  bool specific = (it2 = it1->second.find(key)) !=
1267  it1->second.end();
1268  bool wildcard = (it3 = it1->second.find(defkey)) !=
1269  it1->second.end();
1270 
1271  // Check function is defined somewhere
1272  ASSERTL0(specific || wildcard,
1273  "No such variable " + pVariable
1274  + " in domain " + boost::lexical_cast<string>(pDomain)
1275  + " defined for function " + pName
1276  + " in session file.");
1277 
1278  // If not specific, must be wildcard
1279  if (!specific)
1280  {
1281  it2 = it3;
1282  }
1283 
1284  return it2->second.m_fileVariable;
1285  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
FunctionMap m_functions
Functions.
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 1139 of file SessionReader.cpp.

References ASSERTL0.

1143  {
1144  FunctionMap::const_iterator it1;
1145  FunctionVariableMap::const_iterator it2, it3;
1146  std::string vName = boost::to_upper_copy(pName);
1147 
1148  it1 = m_functions.find(vName);
1149  ASSERTL0 (it1 != m_functions.end(),
1150  std::string("Function '") + pName
1151  + std::string("' not found."));
1152 
1153  // Check for specific and wildcard definitions
1154  pair<std::string,int> key(pVariable,pDomain);
1155  pair<std::string,int> defkey("*",pDomain);
1156  bool specific = (it2 = it1->second.find(key)) !=
1157  it1->second.end();
1158  bool wildcard = (it3 = it1->second.find(defkey)) !=
1159  it1->second.end();
1160 
1161  // Check function is defined somewhere
1162  ASSERTL0(specific || wildcard,
1163  "No such variable " + pVariable
1164  + " in domain " + boost::lexical_cast<string>(pDomain)
1165  + " defined for function " + pName
1166  + " in session file.");
1167 
1168  // If not specific, must be wildcard
1169  if (!specific)
1170  {
1171  it2 = it3;
1172  }
1173 
1174  return it2->second.m_type;
1175  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
FunctionMap m_functions
Functions.
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 1181 of file SessionReader.cpp.

References ASSERTL0.

1185  {
1186  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1187  return GetFunctionType(pName, m_variables[pVar],pDomain);
1188  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
VariableList m_variables
Variables.
enum FunctionType GetFunctionType(const std::string &name, const std::string &variable, const int pDomain=0) const
Returns the type of a given function variable.
const std::string & Nektar::LibUtilities::SessionReader::GetGlobalSysSolnInfo ( const std::string &  variable,
const std::string &  property 
) const

Definition at line 886 of file SessionReader.cpp.

References ASSERTL0, and Nektar::StdRegions::find().

887  {
888  GloSysSolnInfoList::const_iterator iter;
889 
890  ASSERTL0( (iter = GetGloSysSolnList().find(pVariable)) !=
891  GetGloSysSolnList().end(),
892  "Failed to find variable in GlobalSysSolnInfoList");
893 
894  std::string vProperty = boost::to_upper_copy(pProperty);
895  GloSysInfoMap::const_iterator iter1;
896 
897  ASSERTL0( (iter1 = iter->second.find(vProperty)) != iter->second.end(),
898  "Failed to find property: " + vProperty + " in GlobalSysSolnInfoList");
899 
900  return iter1->second;
901  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
static GloSysSolnInfoList & GetGloSysSolnList()
GlobalSysSoln Info map.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
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 ReadGlobalSysSoln 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 156 of file SessionReader.cpp.

157  {
158  static GloSysSolnInfoList gloSysSolnInfoList;
159  return gloSysSolnInfoList;
160  }
std::map< std::string, GloSysInfoMap > GloSysSolnInfoList
Definition: SessionReader.h:81
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 649 of file SessionReader.cpp.

References ASSERTL0.

651  {
652  std::string vName = boost::to_upper_copy(pName);
653  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
654 
655  ASSERTL0(paramIter != m_parameters.end(),
656  "Unable to find requested parameter: " + pName);
657 
658  return paramIter->second;
659  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
const std::string & Nektar::LibUtilities::SessionReader::GetSessionName ( ) const

Returns the session name of the loaded XML document.

Definition at line 586 of file SessionReader.cpp.

587  {
588  return m_sessionName;
589  }
std::string m_sessionName
Session name of the loaded XML document (filename minus ext).
const std::string Nektar::LibUtilities::SessionReader::GetSessionNameRank ( ) const

Returns the session name with process rank.

Output is of the form [sessionName]_P[idx] where idx is the rank of the process.

Definition at line 596 of file SessionReader.cpp.

References Nektar::LibUtilities::PortablePath().

597  {
598  std::string dirname = m_sessionName + "_xml";
599  fs::path pdirname(dirname);
600 
601  std::string vFilename = "P" + boost::lexical_cast<std::string>(m_comm->GetRowComm()->GetRank());
602  fs::path pFilename(vFilename);
603 
604  fs::path fullpath = pdirname / pFilename;
605 
606  return PortablePath(fullpath);
607  }
CommSharedPtr m_comm
Communication object.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
std::string m_sessionName
Session name of the loaded XML document (filename minus ext).
SessionReaderSharedPtr Nektar::LibUtilities::SessionReader::GetSharedThisPtr ( )
inlineprivate

Returns a shared pointer to the current object.

This allows a member function to pass a shared pointer to itself during a call to another function.

Definition at line 706 of file SessionReader.h.

707  {
708  return shared_from_this();
709  }
const std::string & Nektar::LibUtilities::SessionReader::GetSolverInfo ( const std::string &  pProperty) const

Returns the value of the specified solver info property.

Definition at line 767 of file SessionReader.cpp.

References ASSERTL1.

Referenced by GetSolverInfoAsEnum().

769  {
770  std::string vProperty = boost::to_upper_copy(pProperty);
771  SolverInfoMap::const_iterator iter = m_solverInfo.find(vProperty);
772 
773  ASSERTL1(iter != m_solverInfo.end(),
774  "Unable to find requested property: " + pProperty);
775 
776  return iter->second;
777  }
SolverInfoMap m_solverInfo
Solver information properties.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:165
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 552 of file SessionReader.h.

References ASSERTL0, DefinesSolverInfo(), Nektar::StdRegions::find(), GetSolverInfo(), GetSolverInfoEnums(), and Nektar::iterator.

554  {
555  std::string vName = boost::to_upper_copy(pName);
557  "Solver info '" + pName + "' not defined.");
558 
559  std::string vValue = GetSolverInfo(vName);
561  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
562  GetSolverInfoEnums().end(),
563  "Enum for SolverInfo property '" + pName + "' not found.");
564 
566  ASSERTL0((y = x->second.find(vValue)) != x->second.end(),
567  "Value of SolverInfo property '" + pName +
568  "' is invalid.");
569 
570  return T(y->second);
571  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
const std::string & GetSolverInfo(const std::string &pProperty) const
Returns the value of the specified solver info property.
bool DefinesSolverInfo(const std::string &name) const
Checks if a solver info property is specified.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
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 SessionReader.cpp.

Referenced by RegisterDefaultSolverInfo().

141  {
142  static SolverInfoMap solverInfoMap;
143  return solverInfoMap;
144  }
std::map< std::string, std::string > SolverInfoMap
Definition: SessionReader.h:58
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 125 of file SessionReader.cpp.

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

126  {
127  static EnumMapList solverInfoEnums;
128  return solverInfoEnums;
129  }
std::map< std::string, EnumMap > EnumMapList
Definition: SessionReader.h:78
const std::string & Nektar::LibUtilities::SessionReader::GetTag ( const std::string &  pName) const

Returns the value of a specified tag.

Definition at line 1323 of file SessionReader.cpp.

References ASSERTL0.

1324  {
1325  std::string vName = boost::to_upper_copy(pName);
1326  TagMap::const_iterator vTagIterator = m_tags.find(vName);
1327  ASSERTL0(vTagIterator != m_tags.end(),
1328  "Requested tag does not exist.");
1329  return vTagIterator->second;
1330  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 579 of file SessionReader.h.

References ASSERTL0, Nektar::StdRegions::find(), GetSolverInfoEnums(), and Nektar::iterator.

582  {
583  std::string vName = boost::to_upper_copy(pName);
584 
586  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
587  GetSolverInfoEnums().end(),
588  "Enum for property '" + pName + "' not found.");
589 
591  ASSERTL0((y = x->second.find(pValue)) != x->second.end(),
592  "Value of property '" + pValue + "' is invalid.");
593  return T(y->second);
594  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
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 1009 of file SessionReader.cpp.

References ASSERTL0.

1011  {
1012  ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1013  return m_variables[idx];
1014  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
VariableList m_variables
Variables.
std::vector< std::string > Nektar::LibUtilities::SessionReader::GetVariables ( ) const

Returns the names of all variables.

Definition at line 1032 of file SessionReader.cpp.

1033  {
1034  return m_variables;
1035  }
VariableList m_variables
Variables.
void Nektar::LibUtilities::SessionReader::InitSession ( )
private

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 265 of file SessionReader.cpp.

References Nektar::iterator.

266  {
267  m_exprEvaluator.SetRandomSeed((m_comm->GetRank() + 1) * time(NULL));
268 
269  // Split up the communicator
270  PartitionComm();
271 
272  // Partition mesh
273  PartitionMesh();
274 
275  // Parse the XML data in #m_xmlDoc
276  ParseDocument();
277 
278  // Override SOLVERINFO and parameters with any specified on the
279  // command line.
280  CmdLineOverride();
281 
282  // In verbose mode, print out parameters and solver info sections
283  if (m_verbose && m_comm)
284  {
285  if (m_comm->GetRank() == 0 && m_parameters.size() > 0)
286  {
287  cout << "Parameters:" << endl;
289  for (x = m_parameters.begin(); x != m_parameters.end(); ++x)
290  {
291  cout << "\t" << x->first << " = " << x->second << endl;
292  }
293  cout << endl;
294  }
295 
296  if (m_comm->GetRank() == 0 && m_solverInfo.size() > 0)
297  {
298  cout << "Solver Info:" << endl;
300  for (x = m_solverInfo.begin(); x != m_solverInfo.end(); ++x)
301  {
302  cout << "\t" << x->first << " = " << x->second << endl;
303  }
304  cout << endl;
305  }
306  }
307  }
ParameterMap m_parameters
Parameters.
void CmdLineOverride()
Enforce parameters from command line arguments.
SolverInfoMap m_solverInfo
Solver information properties.
AnalyticExpressionEvaluator m_exprEvaluator
Analytic expression evaluator instance.
void PartitionComm()
Partitions the comm object based on session parameters.
CommSharedPtr m_comm
Communication object.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void ParseDocument()
Loads and parses the specified file.
void PartitionMesh()
Partitions the mesh when running in parallel.
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 1378 of file SessionReader.cpp.

References ASSERTL0, and Nektar::LibUtilities::PortablePath().

1381  {
1382  if (pFilename.size() > 3 &&
1383  pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1384  {
1385  ifstream file(pFilename.c_str(),
1386  ios_base::in | ios_base::binary);
1387  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1388  stringstream ss;
1389  io::filtering_streambuf<io::input> in;
1390  in.push(io::gzip_decompressor());
1391  in.push(file);
1392  try
1393  {
1394  io::copy(in, ss);
1395  ss >> (*pDoc);
1396  }
1397  catch (io::gzip_error& e)
1398  {
1399  ASSERTL0(false,
1400  "Error: File '" + pFilename + "' is corrupt.");
1401  }
1402  }
1403  else if (pFilename.size() > 4 &&
1404  pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1405  {
1406  fs::path pdirname(pFilename);
1407  boost::format pad("P%1$07d.xml");
1408  pad % m_comm->GetRank();
1409  fs::path pRankFilename(pad.str());
1410  fs::path fullpath = pdirname / pRankFilename;
1411 
1412  ifstream file(PortablePath(fullpath).c_str());
1413  ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1414  file >> (*pDoc);
1415  }
1416  else
1417  {
1418  ifstream file(pFilename.c_str());
1419  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1420  file >> (*pDoc);
1421  }
1422  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
CommSharedPtr m_comm
Communication object.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
void Nektar::LibUtilities::SessionReader::LoadGeometricInfo ( const std::string &  name,
std::string &  var,
const std::string &  def = "" 
) const

Checks for and load a geometric info string property.

Definition at line 917 of file SessionReader.cpp.

921  {
922  std::string vName = boost::to_upper_copy(pName);
923  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
924  if(iter != m_geometricInfo.end())
925  {
926  pVar = iter->second;
927  }
928  else
929  {
930  pVar = pDefault;
931  }
932  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
void Nektar::LibUtilities::SessionReader::LoadGeometricInfo ( const std::string &  name,
bool &  var,
const bool &  def = false 
) const

Checks for and loads a geometric info boolean property.

Definition at line 938 of file SessionReader.cpp.

942  {
943  std::string vName = boost::to_upper_copy(pName);
944  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
945  if(iter != m_geometricInfo.end())
946  {
947  if (iter->second == "TRUE")
948  {
949  pVar = true;
950  }
951  else
952  {
953  pVar = false;
954  }
955  }
956  else
957  {
958  pVar = pDefault;
959  }
960  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
void Nektar::LibUtilities::SessionReader::LoadGeometricInfo ( const std::string &  name,
NekDouble var,
const NekDouble def = 0.0 
) const

Checks for and loads a geometric info double-precision property.

Definition at line 966 of file SessionReader.cpp.

970  {
971  std::string vName = boost::to_upper_copy(pName);
972  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
973  if(iter != m_geometricInfo.end())
974  {
975  pVar = std::atoi(iter->second.c_str());
976  }
977  else
978  {
979  pVar = pDefault;
980  }
981  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
int &  var 
) const

Load an integer parameter.

Definition at line 665 of file SessionReader.cpp.

References ASSERTL0.

667  {
668  std::string vName = boost::to_upper_copy(pName);
669  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
670  ASSERTL0(paramIter != m_parameters.end(), "Required parameter '" +
671  pName + "' not specified in session.");
672  pVar = (int)floor(paramIter->second);
673  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 679 of file SessionReader.cpp.

681  {
682  std::string vName = boost::to_upper_copy(pName);
683  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
684  if(paramIter != m_parameters.end())
685  {
686  pVar = (int)floor(paramIter->second);
687  }
688  else
689  {
690  pVar = pDefault;
691  }
692  }
ParameterMap m_parameters
Parameters.
void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
NekDouble var 
) const

Load a double precision parameter.

Definition at line 698 of file SessionReader.cpp.

References ASSERTL0.

700  {
701  std::string vName = boost::to_upper_copy(pName);
702  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
703  ASSERTL0(paramIter != m_parameters.end(), "Required parameter '" +
704  pName + "' not specified in session.");
705  pVar = paramIter->second;
706  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 712 of file SessionReader.cpp.

716  {
717  std::string vName = boost::to_upper_copy(pName);
718  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
719  if(paramIter != m_parameters.end())
720  {
721  pVar = paramIter->second;
722  }
723  else
724  {
725  pVar = pDefault;
726  }
727  }
ParameterMap m_parameters
Parameters.
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 797 of file SessionReader.cpp.

801  {
802  std::string vName = boost::to_upper_copy(pName);
803  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
804  if(infoIter != m_solverInfo.end())
805  {
806  pVar = infoIter->second;
807  }
808  else
809  {
810  pVar = pDefault;
811  }
812  }
SolverInfoMap m_solverInfo
Solver information properties.
void Nektar::LibUtilities::SessionReader::MatchGeometricInfo ( const std::string &  name,
const std::string &  trueval,
bool &  var,
const bool &  def = false 
) const

Check if the value of a geometric info string property matches.

Definition at line 987 of file SessionReader.cpp.

992  {
993  std::string vName = boost::to_upper_copy(pName);
994  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
995  if(iter != m_geometricInfo.end())
996  {
997  pVar = boost::iequals(iter->second, pTrueVal);
998  }
999  else
1000  {
1001  pVar = pDefault;
1002  }
1003  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
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 818 of file SessionReader.cpp.

823  {
824  std::string vName = boost::to_upper_copy(pName);
825  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
826  if(infoIter != m_solverInfo.end())
827  {
828  pVar = boost::iequals(infoIter->second, pTrueVal);
829  }
830  else
831  {
832  pVar = pDefault;
833  }
834  }
SolverInfoMap m_solverInfo
Solver information properties.
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 840 of file SessionReader.cpp.

843  {
844  if (DefinesSolverInfo(pName))
845  {
846  std::string vName = boost::to_upper_copy(pName);
847  SolverInfoMap::const_iterator iter = m_solverInfo.find(vName);
848  if(iter != m_solverInfo.end())
849  {
850  return true;
851  }
852  }
853  return false;
854  }
SolverInfoMap m_solverInfo
Solver information properties.
bool DefinesSolverInfo(const std::string &name) const
Checks if a solver info property is specified.
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 541 of file SessionReader.h.

543  {
544  return (GetSolverInfoAsEnum<T>(name) == trueval);
545  }
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 1427 of file SessionReader.cpp.

References ASSERTL0.

1429  {
1430  ASSERTL0(pFilenames.size() > 0, "No filenames for merging.");
1431 
1432  // Read the first document
1433  TiXmlDocument *vMainDoc = new TiXmlDocument;
1434  LoadDoc(pFilenames[0], vMainDoc);
1435 
1436  TiXmlHandle vMainHandle(vMainDoc);
1437  TiXmlElement* vMainNektar =
1438  vMainHandle.FirstChildElement("NEKTAR").Element();
1439 
1440  // Read all subsequent XML documents.
1441  // For each element within the NEKTAR tag, use it to replace the
1442  // version already present in the loaded XML data.
1443  for (int i = 1; i < pFilenames.size(); ++i)
1444  {
1445  if((pFilenames[i].compare(pFilenames[i].size()-3,3,"xml") == 0)
1446  ||(pFilenames[i].compare(pFilenames[i].size()-6,6,"xml.gz") == 0))
1447  {
1448  TiXmlDocument* vTempDoc = new TiXmlDocument;
1449  LoadDoc(pFilenames[i], vTempDoc);
1450 
1451  TiXmlHandle docHandle(vTempDoc);
1452  TiXmlElement* vTempNektar;
1453  vTempNektar = docHandle.FirstChildElement("NEKTAR").Element();
1454  ASSERTL0(vTempNektar, "Unable to find NEKTAR tag in file.");
1455  TiXmlElement* p = vTempNektar->FirstChildElement();
1456 
1457  while (p)
1458  {
1459  TiXmlElement *vMainEntry =
1460  vMainNektar->FirstChildElement(p->Value());
1461  TiXmlElement *q = new TiXmlElement(*p);
1462  if (vMainEntry)
1463  {
1464  vMainNektar->RemoveChild(vMainEntry);
1465  }
1466  vMainNektar->LinkEndChild(q);
1467  p = p->NextSiblingElement();
1468  }
1469 
1470  delete vTempDoc;
1471  }
1472  }
1473  return vMainDoc;
1474  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
void LoadDoc(const std::string &pFilename, TiXmlDocument *pDoc) const
Loads an xml file into a tinyxml doc and decompresses if needed.
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 314 of file SessionReader.cpp.

References Nektar::iterator, Nektar::NekConstants::kGitBranch, Nektar::NekConstants::kGitSha1, and NEKTAR_VERSION.

316  {
317  // List the publically visible options (listed using --help)
318  po::options_description desc("Allowed options");
319  desc.add_options()
320  ("verbose,v", "be verbose")
321  ("version,V", "print version information")
322  ("help,h", "print this help message")
323  ("solverinfo,I", po::value<vector<std::string> >(),
324  "override a SOLVERINFO property")
325  ("parameter,P", po::value<vector<std::string> >(),
326  "override a parameter")
327  ("shared-filesystem,s", "Using shared filesystem.")
328  ("npx", po::value<int>(),
329  "number of procs in X-dir")
330  ("npy", po::value<int>(),
331  "number of procs in Y-dir")
332  ("npz", po::value<int>(),
333  "number of procs in Z-dir")
334  ("nsz", po::value<int>(),
335  "number of slices in Z-dir")
336  ("part-only", po::value<int>(),
337  "only partition mesh into N partitions.")
338  ("part-info", "Output partition information")
339  ;
340 
341  CmdLineArgMap::const_iterator cmdIt;
342  for (cmdIt = GetCmdLineArgMap().begin();
343  cmdIt != GetCmdLineArgMap().end(); ++cmdIt)
344  {
345  std::string names = cmdIt->first;
346  if (cmdIt->second.shortName != "")
347  {
348  names += "," + cmdIt->second.shortName;
349  }
350  if (cmdIt->second.isFlag)
351  {
352  desc.add_options()
353  (names.c_str(), cmdIt->second.description.c_str())
354  ;
355  }
356  else
357  {
358  desc.add_options()
359  (names.c_str(), po::value<std::string>(),
360  cmdIt->second.description.c_str())
361  ;
362  }
363  }
364 
365  // List hidden options (e.g. session file arguments are not actually
366  // specified using the input-file option by the user).
367  po::options_description hidden("Hidden options");
368  hidden.add_options()
369  ("input-file", po::value< vector<string> >(),
370  "input filename")
371  ;
372 
373  // Combine all options for the parser
374  po::options_description all("All options");
375  all.add(desc).add(hidden);
376 
377  // Session file is a positional option
378  po::positional_options_description p;
379  p.add("input-file", -1);
380 
381  // Parse the command-line options
382  po::parsed_options parsed = po::command_line_parser(argc, argv).
383  options(all).
384  positional(p).
385  allow_unregistered().
386  run();
387 
388  // Extract known options to map and update
389  po::store(parsed, m_cmdLineOptions);
390  po::notify(m_cmdLineOptions);
391 
392  // Help message
393  if (m_cmdLineOptions.count("help"))
394  {
395  cout << desc;
396  exit(0);
397  }
398 
399  // Version information
400  if (m_cmdLineOptions.count("version"))
401  {
402  cout << "Nektar++ version " << NEKTAR_VERSION;
403 
404  if (NekConstants::kGitSha1 != "GITDIR-NOTFOUND")
405  {
406  string sha1(NekConstants::kGitSha1);
407  string branch(NekConstants::kGitBranch);
408  boost::replace_all(branch, "refs/heads/", "");
409 
410  cout << " (git changeset " << sha1.substr(0, 8) << ", ";
411 
412  if (branch == "")
413  {
414  cout << "detached head";
415  }
416  else
417  {
418  cout << "head " << branch;
419  }
420 
421  cout << ")";
422  }
423 
424  cout << endl;
425  exit(0);
426  }
427 
428  // Enable verbose mode
429  if (m_cmdLineOptions.count("verbose"))
430  {
431  m_verbose = true;
432  }
433  else
434  {
435  m_verbose = false;
436  }
437 
438  // Print a warning for unknown options
439  std::vector< po::basic_option<char> >::iterator x;
440  for (x = parsed.options.begin(); x != parsed.options.end(); ++x)
441  {
442  if (x->unregistered)
443  {
444  cout << "Warning: Unknown option: " << x->string_key
445  << endl;
446  }
447  }
448 
449  // Return the vector of filename(s) given as positional options
450  if (m_cmdLineOptions.count("input-file"))
451  {
452  return m_cmdLineOptions["input-file"].as<
453  std::vector<std::string> >();
454  }
455  else
456  {
457  return std::vector<std::string>();
458  }
459  }
#define NEKTAR_VERSION
const std::string kGitSha1
boost::program_options::variables_map m_cmdLineOptions
const std::string kGitBranch
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
static CmdLineArgMap & GetCmdLineArgMap()
CmdLine argument map.
void Nektar::LibUtilities::SessionReader::ParseDocument ( )
private

Loads and parses the specified file.

Definition at line 1480 of file SessionReader.cpp.

References ASSERTL0.

1481  {
1482  // Check we actually have a document loaded.
1483  ASSERTL0(m_xmlDoc, "No XML document loaded.");
1484 
1485  // Look for all data in CONDITIONS block.
1486  TiXmlHandle docHandle(m_xmlDoc);
1487  TiXmlElement* e;
1488  e = docHandle.FirstChildElement("NEKTAR").
1489  FirstChildElement("CONDITIONS").Element();
1490 
1491  // Read the various sections of the CONDITIONS block
1492  ReadParameters (e);
1493  ReadSolverInfo (e);
1495  ReadExpressions (e);
1496  ReadVariables (e);
1497  ReadFunctions (e);
1498 
1499  e = docHandle.FirstChildElement("NEKTAR").
1500  FirstChildElement("FILTERS").Element();
1501 
1502  ReadFilters(e);
1503  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
void ReadExpressions(TiXmlElement *conditions)
Reads the EXPRESSIONS section of the XML document.
void ReadSolverInfo(TiXmlElement *conditions)
Reads the SOLVERINFO section of the XML document.
void ReadGlobalSysSolnInfo(TiXmlElement *conditions)
Reads the GLOBALSYSSOLNINFO section of the XML document.
void ReadFunctions(TiXmlElement *conditions)
Reads the FUNCTIONS section of the XML document.
void ReadVariables(TiXmlElement *conditions)
Reads the VARIABLES section of the XML document.
void ReadFilters(TiXmlElement *filters)
Reads the FILTERS section of the XML document.
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
void ReadParameters(TiXmlElement *conditions)
Reads the PARAMETERS section of the XML document.
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 2544 of file SessionReader.cpp.

2548  {
2549  /// Pull out lhs and rhs and eliminate any spaces.
2550  int beg = line.find_first_not_of(" ");
2551  int end = line.find_first_of("=");
2552  // Check for no parameter name
2553  if (beg == end) throw 1;
2554  // Check for no parameter value
2555  if (end != line.find_last_of("=")) throw 1;
2556  // Check for no equals sign
2557  if (end == std::string::npos) throw 1;
2558 
2559  lhs = line.substr(line.find_first_not_of(" "),
2560  end-beg);
2561  lhs = lhs .substr(0, lhs.find_last_not_of(" ")+1);
2562  rhs = line.substr(line.find_last_of("=")+1);
2563  rhs = rhs .substr(rhs.find_first_not_of(" "));
2564  rhs = rhs .substr(0, rhs.find_last_not_of(" ")+1);
2565  }
StandardMatrixTag & lhs
std::string Nektar::LibUtilities::SessionReader::ParseSessionName ( std::vector< std::string > &  filenames)
private

Parse the session name.

Definition at line 465 of file SessionReader.cpp.

References ASSERTL0.

467  {
468  ASSERTL0(!filenames.empty(),
469  "At least one filename expected.");
470 
471  std::string retval = "";
472 
473  // First input file defines the session name
474  std::string fname = filenames[0];
475 
476  // If loading a pre-partitioned mesh, remove _xml extension
477  if (fname.size() > 4 &&
478  fname.substr(fname.size() - 4, 4) == "_xml")
479  {
480  retval = fname.substr(0, fname.find_last_of("_"));
481  }
482  // otherwise remove the .xml extension
483  else if (fname.size() > 4 &&
484  fname.substr(fname.size() - 4, 4) == ".xml")
485  {
486  retval = fname.substr(0, fname.find_last_of("."));
487  }
488  // If compressed .xml.gz, remove both extensions
489  else if (fname.size() > 7 &&
490  fname.substr(fname.size() - 7, 7) == ".xml.gz")
491  {
492  retval = fname.substr(0, fname.find_last_of("."));
493  retval = retval.substr(0, retval.find_last_of("."));
494  }
495 
496  return retval;
497  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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. The grid is defined by the PROC_X parameter which, if specified, gives the number of processes spanned by the Fourier direction. PROC_X must exactly divide the total number of processes or an error is thrown.

Definition at line 1813 of file SessionReader.cpp.

References ASSERTL0.

1814  {
1815  if (m_comm->GetSize() > 1)
1816  {
1817  int nProcZ = 1;
1818  int nProcY = 1;
1819  int nProcX = 1;
1820  int nStripZ = 1;
1821  if (DefinesCmdLineArgument("npx")) {
1822  nProcX = GetCmdLineArgument<int>("npx");
1823  }
1824  if (DefinesCmdLineArgument("npy")) {
1825  nProcY = GetCmdLineArgument<int>("npy");
1826  }
1827  if (DefinesCmdLineArgument("npz")) {
1828  nProcZ = GetCmdLineArgument<int>("npz");
1829  }
1830  if (DefinesCmdLineArgument("nsz")) {
1831  nStripZ = GetCmdLineArgument<int>("nsz");
1832  }
1833  ASSERTL0(m_comm->GetSize() % (nProcZ*nProcY*nProcX) == 0,
1834  "Cannot exactly partition using PROC_Z value.");
1835  ASSERTL0(nProcZ % nProcY == 0,
1836  "Cannot exactly partition using PROC_Y value.");
1837  ASSERTL0(nProcY % nProcX == 0,
1838  "Cannot exactly partition using PROC_X value.");
1839 
1840  // Number of processes associated with the spectral method
1841  int nProcSm = nProcZ * nProcY * nProcX;
1842 
1843  // Number of processes associated with the spectral element
1844  // method.
1845  int nProcSem = m_comm->GetSize() / nProcSm;
1846 
1847  m_comm->SplitComm(nProcSm,nProcSem);
1848  m_comm->GetColumnComm()->SplitComm(nProcZ/nStripZ,nStripZ);
1849  m_comm->GetColumnComm()->GetColumnComm()->SplitComm(
1850  (nProcY*nProcX),nProcZ/nStripZ);
1851  m_comm->GetColumnComm()->GetColumnComm()->GetColumnComm()
1852  ->SplitComm(nProcX,nProcY);
1853  }
1854  }
bool DefinesCmdLineArgument(const std::string &pName) const
Checks if a specified cmdline argument has been given.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
CommSharedPtr m_comm
Communication object.
void Nektar::LibUtilities::SessionReader::PartitionMesh ( )
private

Partitions the mesh when running in parallel.

Definition at line 1533 of file SessionReader.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Gs::Finalise(), Nektar::LibUtilities::GetMeshPartitionFactory(), Nektar::LibUtilities::GetSize(), Nektar::iterator, Nektar::LibUtilities::PortablePath(), and Nektar::LibUtilities::ReduceMax.

1534  {
1535  ASSERTL0(m_comm.get(), "Communication not initialised.");
1536 
1537  // Get row of comm, or the whole comm if not split
1538  CommSharedPtr vCommMesh = m_comm->GetRowComm();
1539  const bool isRoot = (m_comm->GetRank() == 0);
1540 
1541  // Delete any existing loaded mesh
1542  if (m_xmlDoc)
1543  {
1544  delete m_xmlDoc;
1545  }
1546 
1547  // Load file for root process only (since this is always needed)
1548  // and determine if the provided geometry has already been
1549  // partitioned. This will be the case if the user provides the
1550  // directory of mesh partitions as an input. Partitioned geometries
1551  // have the attribute
1552  // PARTITION=X
1553  // where X is the number of the partition (and should match the
1554  // process rank). The result is shared with all other processes.
1555  int isPartitioned = 0;
1556  if (isRoot)
1557  {
1559  if (DefinesElement("Nektar/Geometry"))
1560  {
1561  if (GetElement("Nektar/Geometry")->Attribute("PARTITION"))
1562  {
1563  cout << "Using pre-partitioned mesh." << endl;
1564  isPartitioned = 1;
1565  }
1566  }
1567  }
1568  GetComm()->AllReduce(isPartitioned, LibUtilities::ReduceMax);
1569 
1570  // If the mesh is already partitioned, we are done. Remaining
1571  // processes must load their partitions.
1572  if (isPartitioned) {
1573  if (!isRoot)
1574  {
1576  }
1577  return;
1578  }
1579 
1580  // Default partitioner to use is Metis. Use Scotch as default
1581  // if it is installed. Override default with command-line flags
1582  // if they are set.
1583  string vPartitionerName = "Metis";
1584  if (GetMeshPartitionFactory().ModuleExists("Scotch"))
1585  {
1586  vPartitionerName = "Scotch";
1587  }
1588  if (DefinesCmdLineArgument("use-metis"))
1589  {
1590  vPartitionerName = "Metis";
1591  }
1592  if (DefinesCmdLineArgument("use-scotch"))
1593  {
1594  vPartitionerName = "Scotch";
1595  }
1596 
1597  // Mesh has not been partitioned so do partitioning if required.
1598  // Note in the serial case nothing is done as we have already loaded
1599  // the mesh.
1600  if (DefinesCmdLineArgument("part-only"))
1601  {
1602  // Perform partitioning of the mesh only. For this we insist
1603  // the code is run in serial (parallel execution is pointless).
1604  ASSERTL0(GetComm()->GetSize() == 1,
1605  "The 'part-only' option should be used in serial.");
1606 
1607  // Number of partitions is specified by the parameter.
1608  int nParts = GetCmdLineArgument<int>("part-only");
1610  MeshPartitionSharedPtr vPartitioner =
1612  vPartitionerName, vSession);
1613  vPartitioner->PartitionMesh(nParts, true);
1614  vPartitioner->WriteAllPartitions(vSession);
1615  vPartitioner->GetCompositeOrdering(m_compOrder);
1616  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1617 
1618  if (isRoot && DefinesCmdLineArgument("part-info"))
1619  {
1620  vPartitioner->PrintPartInfo(std::cout);
1621  }
1622 
1623  Finalise();
1624  exit(0);
1625  }
1626  else if (vCommMesh->GetSize() > 1)
1627  {
1629  int nParts = vCommMesh->GetSize();
1630  if (DefinesCmdLineArgument("shared-filesystem"))
1631  {
1632  CommSharedPtr vComm = GetComm();
1633  vector<unsigned int> keys, vals;
1634  int i;
1635 
1636  if (vComm->GetRank() == 0)
1637  {
1639 
1640  MeshPartitionSharedPtr vPartitioner =
1642  vPartitionerName, vSession);
1643  vPartitioner->PartitionMesh(nParts, true);
1644  vPartitioner->WriteAllPartitions(vSession);
1645  vPartitioner->GetCompositeOrdering(m_compOrder);
1646  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1647 
1648  // Communicate orderings to the other processors.
1649 
1650  // First send sizes of the orderings and boundary
1651  // regions to allocate storage on the remote end.
1652  keys.resize(2);
1653  keys[0] = m_compOrder.size();
1654  keys[1] = m_bndRegOrder.size();
1655 
1656  for (i = 1; i < vComm->GetSize(); ++i)
1657  {
1658  vComm->Send(i, keys);
1659  }
1660 
1661  // Construct the keys and sizes of values for composite
1662  // ordering
1664  keys.resize(m_compOrder.size());
1665  vals.resize(m_compOrder.size());
1666 
1667  for (cIt = m_compOrder.begin(), i = 0;
1668  cIt != m_compOrder.end(); ++cIt, ++i)
1669  {
1670  keys[i] = cIt->first;
1671  vals[i] = cIt->second.size();
1672  }
1673 
1674  // Send across data.
1675  for (i = 1; i < vComm->GetSize(); ++i)
1676  {
1677  vComm->Send(i, keys);
1678  vComm->Send(i, vals);
1679 
1680  for (cIt = m_compOrder.begin();
1681  cIt != m_compOrder.end(); ++cIt)
1682  {
1683  vComm->Send(i, cIt->second);
1684  }
1685  }
1686 
1687  // Construct the keys and sizes of values for composite
1688  // ordering
1690  keys.resize(m_bndRegOrder.size());
1691  vals.resize(m_bndRegOrder.size());
1692 
1693  for (bIt = m_bndRegOrder.begin(), i = 0;
1694  bIt != m_bndRegOrder.end(); ++bIt, ++i)
1695  {
1696  keys[i] = bIt->first;
1697  vals[i] = bIt->second.size();
1698  }
1699 
1700  // Send across data.
1701  for (i = 1; i < vComm->GetSize(); ++i)
1702  {
1703  vComm->Send(i, keys);
1704  vComm->Send(i, vals);
1705 
1706  for (bIt = m_bndRegOrder.begin();
1707  bIt != m_bndRegOrder.end(); ++bIt)
1708  {
1709  vComm->Send(i, bIt->second);
1710  }
1711  }
1712 
1713  if (DefinesCmdLineArgument("part-info"))
1714  {
1715  vPartitioner->PrintPartInfo(std::cout);
1716  }
1717  }
1718  else
1719  {
1720  keys.resize(2);
1721  vComm->Recv(0, keys);
1722 
1723  int cmpSize = keys[0];
1724  int bndSize = keys[1];
1725 
1726  keys.resize(cmpSize);
1727  vals.resize(cmpSize);
1728  vComm->Recv(0, keys);
1729  vComm->Recv(0, vals);
1730 
1731  for (int i = 0; i < keys.size(); ++i)
1732  {
1733  vector<unsigned int> tmp(vals[i]);
1734  vComm->Recv(0, tmp);
1735  m_compOrder[keys[i]] = tmp;
1736  }
1737 
1738  keys.resize(bndSize);
1739  vals.resize(bndSize);
1740  vComm->Recv(0, keys);
1741  vComm->Recv(0, vals);
1742 
1743  for (int i = 0; i < keys.size(); ++i)
1744  {
1745  vector<unsigned int> tmp(vals[i]);
1746  vComm->Recv(0, tmp);
1747  m_bndRegOrder[keys[i]] = tmp;
1748  }
1749  }
1750  }
1751  else
1752  {
1753  // Need to load mesh on non-root processes.
1754  if (!isRoot)
1755  {
1757  }
1758 
1759  // Partitioner now operates in parallel
1760  // Each process receives partitioning over interconnect
1761  // and writes its own session file to the working directory.
1762  MeshPartitionSharedPtr vPartitioner =
1764  vPartitionerName, vSession);
1765  vPartitioner->PartitionMesh(nParts, false);
1766  vPartitioner->WriteLocalPartition(vSession);
1767  vPartitioner->GetCompositeOrdering(m_compOrder);
1768  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1769 
1770  if (DefinesCmdLineArgument("part-info") && isRoot)
1771  {
1772  vPartitioner->PrintPartInfo(std::cout);
1773  }
1774  }
1775  m_comm->Block();
1776 
1777  std::string dirname = GetSessionName() + "_xml";
1778  fs::path pdirname(dirname);
1779  boost::format pad("P%1$07d.xml");
1780  pad % m_comm->GetRowComm()->GetRank();
1781  fs::path pFilename(pad.str());
1782  fs::path fullpath = pdirname / pFilename;
1783 
1784  std::string vFilename = PortablePath(fullpath);
1785 
1786  if (m_xmlDoc)
1787  {
1788  delete m_xmlDoc;
1789  }
1790  m_xmlDoc = new TiXmlDocument(vFilename);
1791 
1792  ASSERTL0(m_xmlDoc, "Failed to create XML document object.");
1793 
1794  bool loadOkay = m_xmlDoc->LoadFile(vFilename);
1795  ASSERTL0(loadOkay, "Unable to load file: " + vFilename +
1796  ". Check XML standards compliance. Error on line: " +
1797  boost::lexical_cast<std::string>(m_xmlDoc->Row()));
1798  }
1799  else
1800  {
1802  }
1803  }
bool DefinesCmdLineArgument(const std::string &pName) const
Checks if a specified cmdline argument has been given.
TiXmlElement * GetElement(const std::string &pPath)
Provides direct access to the TiXmlElement specified.
TiXmlDocument * MergeDoc(const std::vector< std::string > &pFilenames) const
Creates an XML document from a list of input files.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
CommSharedPtr & GetComm()
Returns the communication object.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
boost::shared_ptr< MeshPartition > MeshPartitionSharedPtr
SessionReaderSharedPtr GetSharedThisPtr()
Returns a shared pointer to the current object.
CommSharedPtr m_comm
Communication object.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
void Finalise()
Finalises the session.
const std::string & GetSessionName() const
Returns the session name of the loaded XML document.
std::vector< std::string > m_filenames
Filenames.
int GetSize(const Array< OneD, const NekDouble > &x)
Definition: NodalUtil.cpp:111
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
bool DefinesElement(const std::string &pPath) const
Tests if a specified element is defined in the XML document.
CompositeOrdering m_compOrder
Map of original composite ordering for parallel periodic bcs.
MeshPartitionFactory & GetMeshPartitionFactory()
BndRegionOrdering m_bndRegOrder
Map of original boundary region ordering for parallel periodic bcs.
void Nektar::LibUtilities::SessionReader::ReadExpressions ( TiXmlElement *  conditions)
private

Reads the EXPRESSIONS section of the XML document.

Definition at line 2165 of file SessionReader.cpp.

References ASSERTL0, and Nektar::iterator.

2166  {
2167  m_expressions.clear();
2168 
2169  if (!conditions)
2170  {
2171  return;
2172  }
2173 
2174  TiXmlElement *expressionsElement =
2175  conditions->FirstChildElement("EXPRESSIONS");
2176 
2177  if (expressionsElement)
2178  {
2179  TiXmlElement *expr = expressionsElement->FirstChildElement("E");
2180 
2181  while (expr)
2182  {
2183  stringstream tagcontent;
2184  tagcontent << *expr;
2185  ASSERTL0(expr->Attribute("NAME"),
2186  "Missing NAME attribute in expression "
2187  "definition: \n\t'" + tagcontent.str() + "'");
2188  std::string nameString = expr->Attribute("NAME");
2189  ASSERTL0(!nameString.empty(),
2190  "Expressions must have a non-empty name: \n\t'"
2191  + tagcontent.str() + "'");
2192 
2193  ASSERTL0(expr->Attribute("VALUE"),
2194  "Missing VALUE attribute in expression "
2195  "definition: \n\t'" + tagcontent.str() + "'");
2196  std::string valString = expr->Attribute("VALUE");
2197  ASSERTL0(!valString.empty(),
2198  "Expressions must have a non-empty value: \n\t'"
2199  + tagcontent.str() + "'");
2200 
2201  ExpressionMap::iterator exprIter
2202  = m_expressions.find(nameString);
2203  ASSERTL0(exprIter == m_expressions.end(),
2204  std::string("Expression '") + nameString
2205  + std::string("' already specified."));
2206 
2207  m_expressions[nameString] = valString;
2208  expr = expr->NextSiblingElement("E");
2209  }
2210  }
2211  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
ExpressionMap m_expressions
Expressions.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void Nektar::LibUtilities::SessionReader::ReadFilters ( TiXmlElement *  filters)
private

Reads the FILTERS section of the XML document.

Definition at line 2503 of file SessionReader.cpp.

References ASSERTL0.

2504  {
2505  if (!filters)
2506  {
2507  return;
2508  }
2509 
2510  m_filters.clear();
2511 
2512  TiXmlElement *filter = filters->FirstChildElement("FILTER");
2513  while (filter)
2514  {
2515  ASSERTL0(filter->Attribute("TYPE"),
2516  "Missing attribute 'TYPE' for filter.");
2517  std::string typeStr = filter->Attribute("TYPE");
2518 
2519  std::map<std::string, std::string> vParams;
2520 
2521  TiXmlElement *param = filter->FirstChildElement("PARAM");
2522  while (param)
2523  {
2524  ASSERTL0(param->Attribute("NAME"),
2525  "Missing attribute 'NAME' for parameter in filter "
2526  + typeStr + "'.");
2527  std::string nameStr = param->Attribute("NAME");
2528 
2529  ASSERTL0(param->GetText(), "Empty value string for param.");
2530  std::string valueStr = param->GetText();
2531 
2532  vParams[nameStr] = valueStr;
2533 
2534  param = param->NextSiblingElement("PARAM");
2535  }
2536 
2537  m_filters.push_back(
2538  std::pair<std::string, FilterParams>(typeStr, vParams));
2539 
2540  filter = filter->NextSiblingElement("FILTER");
2541  }
2542  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
FilterMap m_filters
Filters map.
void Nektar::LibUtilities::SessionReader::ReadFunctions ( TiXmlElement *  conditions)
private

Reads the FUNCTIONS section of the XML document.

Definition at line 2300 of file SessionReader.cpp.

References ASSERTL0, Nektar::LibUtilities::eFunctionTypeExpression, Nektar::LibUtilities::eFunctionTypeFile, Nektar::LibUtilities::eFunctionTypeTransientFile, Nektar::iterator, Nektar::LibUtilities::FunctionVariableDefinition::m_expression, Nektar::LibUtilities::FunctionVariableDefinition::m_filename, Nektar::LibUtilities::FunctionVariableDefinition::m_fileVariable, and Nektar::LibUtilities::FunctionVariableDefinition::m_type.

2301  {
2302  m_functions.clear();
2303 
2304  if (!conditions)
2305  {
2306  return;
2307  }
2308 
2309  // Scan through conditions section looking for functions.
2310  TiXmlElement *function = conditions->FirstChildElement("FUNCTION");
2311  while (function)
2312  {
2313  stringstream tagcontent;
2314  tagcontent << *function;
2315 
2316  // Every function must have a NAME attribute
2317  ASSERTL0(function->Attribute("NAME"),
2318  "Functions must have a NAME attribute defined in XML "
2319  "element: \n\t'" + tagcontent.str() + "'");
2320  std::string functionStr = function->Attribute("NAME");
2321  ASSERTL0(!functionStr.empty(),
2322  "Functions must have a non-empty name in XML "
2323  "element: \n\t'" + tagcontent.str() + "'");
2324 
2325  // Store function names in uppercase to remain case-insensitive.
2326  boost::to_upper(functionStr);
2327 
2328  // Retrieve first entry (variable, or file)
2329  TiXmlElement *variable = function->FirstChildElement();
2330 
2331  // Create new function structure with default type of none.
2332  FunctionVariableMap functionVarMap;
2333 
2334  // Process all entries in the function block
2335  while (variable)
2336  {
2337  FunctionVariableDefinition funcDef;
2338  std::string conditionType = variable->Value();
2339 
2340  // If no var is specified, assume wildcard
2341  std::string variableStr;
2342  if (!variable->Attribute("VAR"))
2343  {
2344  variableStr = "*";
2345  }
2346  else
2347  {
2348  variableStr = variable->Attribute("VAR");
2349  }
2350 
2351  // Parse list of variables
2352  std::vector<std::string> variableList;
2353  ParseUtils::GenerateOrderedStringVector(variableStr.c_str(),
2354  variableList);
2355 
2356  // If no domain string put to 0
2357  std::string domainStr;
2358  if (!variable->Attribute("DOMAIN"))
2359  {
2360  domainStr = "0";
2361  }
2362  else
2363  {
2364  domainStr = variable->Attribute("DOMAIN");
2365  }
2366 
2367  // Parse list of variables
2368  std::vector<std::string> varSplit;
2369  std::vector<unsigned int> domainList;
2370  ParseUtils::GenerateSeqVector(domainStr.c_str(), domainList);
2371 
2372  // Expressions are denoted by E
2373  if (conditionType == "E")
2374  {
2375  funcDef.m_type = eFunctionTypeExpression;
2376 
2377  // Expression must have a VALUE.
2378  ASSERTL0(variable->Attribute("VALUE"),
2379  "Attribute VALUE expected for function '"
2380  + functionStr + "'.");
2381  std::string fcnStr = variable->Attribute("VALUE");
2382 
2383  ASSERTL0(!fcnStr.empty(),
2384  (std::string("Expression for var: ")
2385  + variableStr
2386  + std::string(" must be specified.")).c_str());
2387 
2388  SubstituteExpressions(fcnStr);
2389 
2390  // set expression
2391  funcDef.m_expression = MemoryManager<Equation>
2393  }
2394 
2395  // Files are denoted by F
2396  else if (conditionType == "F")
2397  {
2398  if (variable->Attribute("TIMEDEPENDENT") &&
2399  boost::lexical_cast<bool>(variable->Attribute("TIMEDEPENDENT")))
2400  {
2401  funcDef.m_type = eFunctionTypeTransientFile;
2402  }
2403  else
2404  {
2405  funcDef.m_type = eFunctionTypeFile;
2406  }
2407 
2408  // File must have a FILE.
2409  ASSERTL0(variable->Attribute("FILE"),
2410  "Attribute FILE expected for function '"
2411  + functionStr + "'.");
2412  std::string filenameStr = variable->Attribute("FILE");
2413 
2414  ASSERTL0(!filenameStr.empty(),
2415  "A filename must be specified for the FILE "
2416  "attribute of function '" + functionStr
2417  + "'.");
2418 
2419  std::vector<std::string> fSplit;
2420  boost::split(fSplit, filenameStr, boost::is_any_of(":"));
2421 
2422  ASSERTL0(fSplit.size() == 1 || fSplit.size() == 2,
2423  "Incorrect filename specification in function "
2424  + functionStr + "'. "
2425  "Specify variables inside file as: "
2426  "filename:var1,var2");
2427 
2428  // set the filename
2429  funcDef.m_filename = fSplit[0];
2430 
2431  if (fSplit.size() == 2)
2432  {
2433  ASSERTL0(variableList[0] != "*",
2434  "Filename variable mapping not valid "
2435  "when using * as a variable inside "
2436  "function '" + functionStr + "'.");
2437 
2438  boost::split(
2439  varSplit, fSplit[1], boost::is_any_of(","));
2440  ASSERTL0(varSplit.size() == variableList.size(),
2441  "Filename variables should contain the "
2442  "same number of variables defined in "
2443  "VAR in function " + functionStr + "'.");
2444  }
2445  }
2446 
2447  // Nothing else supported so throw an error
2448  else
2449  {
2450  stringstream tagcontent;
2451  tagcontent << *variable;
2452 
2453  ASSERTL0(false,
2454  "Identifier " + conditionType + " in function "
2455  + std::string(function->Attribute("NAME"))
2456  + " is not recognised in XML element: \n\t'"
2457  + tagcontent.str() + "'");
2458  }
2459 
2460 
2461 
2462  // Add variables to function
2463  for (unsigned int i = 0; i < variableList.size(); ++i)
2464  {
2465  for(unsigned int j = 0; j < domainList.size(); ++j)
2466  {
2467  // Check it has not already been defined
2468  pair<std::string,int> key(variableList[i],domainList[j]);
2470  = functionVarMap.find(key);
2471  ASSERTL0(fcnsIter == functionVarMap.end(),
2472  "Error setting expression '" + variableList[i]
2473  + " in domain "
2474  + boost::lexical_cast<std::string>(domainList[j])
2475  + "' in function '" + functionStr + "'. "
2476  "Expression has already been defined.");
2477 
2478  if (varSplit.size() > 0)
2479  {
2480  FunctionVariableDefinition funcDef2 = funcDef;
2481  funcDef2.m_fileVariable = varSplit[i];
2482  functionVarMap[key] = funcDef2;
2483  }
2484  else
2485  {
2486  functionVarMap[key] = funcDef;
2487  }
2488  }
2489  }
2490 
2491  variable = variable->NextSiblingElement();
2492  }
2493  // Add function definition to map
2494  m_functions[functionStr] = functionVarMap;
2495  function = function->NextSiblingElement("FUNCTION");
2496  }
2497  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:142
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::map< std::pair< std::string, int >, FunctionVariableDefinition > FunctionVariableMap
void SubstituteExpressions(std::string &expr)
Substitutes expressions defined in the XML document.
SessionReaderSharedPtr GetSharedThisPtr()
Returns a shared pointer to the current object.
static bool GenerateSeqVector(const char *const str, std::vector< unsigned int > &vec)
Definition: ParseUtils.hpp:78
FunctionMap m_functions
Functions.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void Nektar::LibUtilities::SessionReader::ReadGlobalSysSolnInfo ( TiXmlElement *  conditions)
private

Reads the GLOBALSYSSOLNINFO section of the XML document.

Definition at line 2031 of file SessionReader.cpp.

References ASSERTL0, Nektar::StdRegions::find(), and Nektar::iterator.

2032  {
2033  GetGloSysSolnList().clear();
2034 
2035  if (!conditions)
2036  {
2037  return;
2038  }
2039 
2040  TiXmlElement *GlobalSys =
2041  conditions->FirstChildElement("GLOBALSYSSOLNINFO");
2042 
2043  if(!GlobalSys)
2044  {
2045  return;
2046  }
2047 
2048  TiXmlElement *VarInfo = GlobalSys->FirstChildElement("V");
2049 
2050  while (VarInfo)
2051  {
2052  std::stringstream tagcontent;
2053  tagcontent << *VarInfo;
2054  ASSERTL0(VarInfo->Attribute("VAR"),
2055  "Missing VAR attribute in GobalSysSolnInfo XML "
2056  "element: \n\t'" + tagcontent.str() + "'");
2057 
2058  std::string VarList = VarInfo->Attribute("VAR");
2059  ASSERTL0(!VarList.empty(),
2060  "VAR attribute must be non-empty in XML element:\n\t'"
2061  + tagcontent.str() + "'");
2062 
2063  // generate a list of variables.
2064  std::vector<std::string> varStrings;
2066  VarList.c_str(),varStrings);
2067 
2068  ASSERTL0(valid,"Unable to process list of variable in XML "
2069  "element \n\t'" + tagcontent.str() + "'");
2070 
2071  if(varStrings.size())
2072  {
2073  TiXmlElement *SysSolnInfo = VarInfo->FirstChildElement("I");
2074 
2075  while (SysSolnInfo)
2076  {
2077  tagcontent.clear();
2078  tagcontent << *SysSolnInfo;
2079  // read the property name
2080  ASSERTL0(SysSolnInfo->Attribute("PROPERTY"),
2081  "Missing PROPERTY attribute in "
2082  "GlobalSysSolnInfo for variable(s) '"
2083  + VarList + "' in XML element: \n\t'"
2084  + tagcontent.str() + "'");
2085 
2086  std::string SysSolnProperty =
2087  SysSolnInfo->Attribute("PROPERTY");
2088 
2089  ASSERTL0(!SysSolnProperty.empty(),
2090  "GlobalSysSolnIno properties must have a "
2091  "non-empty name for variable(s) : '"
2092  + VarList + "' in XML element: \n\t'"
2093  + tagcontent.str() + "'");
2094 
2095  // make sure that solver property is capitalised
2096  std::string SysSolnPropertyUpper =
2097  boost::to_upper_copy(SysSolnProperty);
2098 
2099  // read the value
2100  ASSERTL0(SysSolnInfo->Attribute("VALUE"),
2101  "Missing VALUE attribute in GlobalSysSolnInfo "
2102  "for variable(s) '" + VarList
2103  + "' in XML element: \n\t"
2104  + tagcontent.str() + "'");
2105 
2106  std::string SysSolnValue =
2107  SysSolnInfo->Attribute("VALUE");
2108  ASSERTL0(!SysSolnValue.empty(),
2109  "GlobalSysSolnInfo properties must have a "
2110  "non-empty value for variable(s) '"
2111  + VarList + "' in XML element: \n\t'"
2112  + tagcontent.str() + "'");
2113 
2114  // Store values under variable map.
2115  for(int i = 0; i < varStrings.size(); ++i)
2116  {
2118  if ((x = GetGloSysSolnList().find(varStrings[i])) ==
2119  GetGloSysSolnList().end())
2120  {
2121  (GetGloSysSolnList()[varStrings[i]])[
2122  SysSolnPropertyUpper] = SysSolnValue;
2123  }
2124  else
2125  {
2126  x->second[SysSolnPropertyUpper] = SysSolnValue;
2127  }
2128  }
2129 
2130  SysSolnInfo = SysSolnInfo->NextSiblingElement("I");
2131  }
2132  VarInfo = VarInfo->NextSiblingElement("V");
2133  }
2134  }
2135 
2136  if (m_verbose && GetGloSysSolnList().size() > 0 && m_comm)
2137  {
2138  if(m_comm->GetRank() == 0)
2139  {
2140  cout << "GlobalSysSoln Info:" << endl;
2141 
2143  for (x = GetGloSysSolnList().begin();
2144  x != GetGloSysSolnList().end();
2145  ++x)
2146  {
2147  cout << "\t Variable: " << x->first << endl;
2148 
2150  for (y = x->second.begin(); y != x->second.end(); ++y)
2151  {
2152  cout << "\t\t " << y->first << " = " << y->second
2153  << endl;
2154  }
2155  }
2156  cout << endl;
2157  }
2158  }
2159  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:142
static GloSysSolnInfoList & GetGloSysSolnList()
GlobalSysSoln Info map.
CommSharedPtr m_comm
Communication object.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
void Nektar::LibUtilities::SessionReader::ReadParameters ( TiXmlElement *  conditions)
private

Reads the PARAMETERS section of the XML document.

Definition at line 1860 of file SessionReader.cpp.

References ASSERTL0, Nektar::LibUtilities::Equation::Evaluate(), and Nektar::lhs.

1861  {
1862  m_parameters.clear();
1863 
1864  if (!conditions)
1865  {
1866  return;
1867  }
1868 
1869  TiXmlElement *parametersElement = conditions->FirstChildElement(
1870  "PARAMETERS");
1871 
1872  // See if we have parameters defined. They are optional so we go on
1873  // if not.
1874  if (parametersElement)
1875  {
1876  TiXmlElement *parameter =
1877  parametersElement->FirstChildElement("P");
1878 
1879  ParameterMap caseSensitiveParameters;
1880 
1881  // Multiple nodes will only occur if there is a comment in
1882  // between definitions.
1883  while (parameter)
1884  {
1885  stringstream tagcontent;
1886  tagcontent << *parameter;
1887  TiXmlNode *node = parameter->FirstChild();
1888 
1889  while (node && node->Type() != TiXmlNode::TINYXML_TEXT)
1890  {
1891  node = node->NextSibling();
1892  }
1893 
1894  if (node)
1895  {
1896  // Format is "paramName = value"
1897  std::string line = node->ToText()->Value(), lhs, rhs;
1898 
1899  try {
1900  ParseEquals(line, lhs, rhs);
1901  }
1902  catch (...)
1903  {
1904  ASSERTL0(false, "Syntax error in parameter "
1905  "expression '" + line
1906  + "' in XML element: \n\t'"
1907  + tagcontent.str() + "'");
1908  }
1909 
1910  // We want the list of parameters to have their RHS
1911  // evaluated, so we use the expression evaluator to do
1912  // the dirty work.
1913  if (!lhs.empty() && !rhs.empty())
1914  {
1915  NekDouble value=0.0;
1916  try
1917  {
1918  LibUtilities::Equation expession(
1919  GetSharedThisPtr(), rhs);
1920  value = expession.Evaluate();
1921  }
1922  catch (const std::runtime_error &)
1923  {
1924  ASSERTL0(false,
1925  "Error evaluating parameter expression"
1926  " '" + rhs + "' in XML element: \n\t'"
1927  + tagcontent.str() + "'");
1928  }
1930  caseSensitiveParameters[lhs] = value;
1931  boost::to_upper(lhs);
1932  m_parameters[lhs] = value;
1933  }
1934  }
1935 
1936  parameter = parameter->NextSiblingElement();
1937  }
1938  }
1939  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
void SetParameter(std::string const &name, NekDouble value)
This function behaves in the same way as SetParameters, but it only adds one parameter and it does no...
AnalyticExpressionEvaluator m_exprEvaluator
Analytic expression evaluator instance.
SessionReaderSharedPtr GetSharedThisPtr()
Returns a shared pointer to the current object.
StandardMatrixTag & lhs
void ParseEquals(const std::string &line, std::string &lhs, std::string &rhs)
Parse a string in the form lhs = rhs.
double NekDouble
std::map< std::string, NekDouble > ParameterMap
Definition: SessionReader.h:59
void Nektar::LibUtilities::SessionReader::ReadSolverInfo ( TiXmlElement *  conditions)
private

Reads the SOLVERINFO section of the XML document.

Definition at line 1945 of file SessionReader.cpp.

References ASSERTL0.

1946  {
1947  m_solverInfo.clear();
1949 
1950  if (!conditions)
1951  {
1952  return;
1953  }
1954 
1955  TiXmlElement *solverInfoElement =
1956  conditions->FirstChildElement("SOLVERINFO");
1957 
1958  if (solverInfoElement)
1959  {
1960  TiXmlElement *solverInfo =
1961  solverInfoElement->FirstChildElement("I");
1962 
1963  while (solverInfo)
1964  {
1965  std::stringstream tagcontent;
1966  tagcontent << *solverInfo;
1967  // read the property name
1968  ASSERTL0(solverInfo->Attribute("PROPERTY"),
1969  "Missing PROPERTY attribute in solver info "
1970  "XML element: \n\t'" + tagcontent.str() + "'");
1971  std::string solverProperty =
1972  solverInfo->Attribute("PROPERTY");
1973  ASSERTL0(!solverProperty.empty(),
1974  "PROPERTY attribute must be non-empty in XML "
1975  "element: \n\t'" + tagcontent.str() + "'");
1976 
1977  // make sure that solver property is capitalised
1978  std::string solverPropertyUpper =
1979  boost::to_upper_copy(solverProperty);
1980 
1981  // read the value
1982  ASSERTL0(solverInfo->Attribute("VALUE"),
1983  "Missing VALUE attribute in solver info "
1984  "XML element: \n\t'" + tagcontent.str() + "'");
1985  std::string solverValue = solverInfo->Attribute("VALUE");
1986  ASSERTL0(!solverValue.empty(),
1987  "VALUE attribute must be non-empty in XML "
1988  "element: \n\t'" + tagcontent.str() + "'");
1989 
1990  EnumMapList::const_iterator propIt =
1991  GetSolverInfoEnums().find(solverPropertyUpper);
1992  if (propIt != GetSolverInfoEnums().end())
1993  {
1994  EnumMap::const_iterator valIt =
1995  propIt->second.find(solverValue);
1996  ASSERTL0(valIt != propIt->second.end(),
1997  "Value '" + solverValue + "' is not valid for "
1998  "property '" + solverProperty + "'");
1999  }
2000 
2001  // Set Variable
2002  m_solverInfo[solverPropertyUpper] = solverValue;
2003  solverInfo = solverInfo->NextSiblingElement("I");
2004  }
2005  }
2006 
2007  if (m_comm && m_comm->GetRowComm()->GetSize() > 1)
2008  {
2009  ASSERTL0(
2010  m_solverInfo["GLOBALSYSSOLN"] == "IterativeFull" ||
2011  m_solverInfo["GLOBALSYSSOLN"] == "IterativeStaticCond" ||
2012  m_solverInfo["GLOBALSYSSOLN"] ==
2013  "IterativeMultiLevelStaticCond" ||
2014  m_solverInfo["GLOBALSYSSOLN"] == "XxtFull" ||
2015  m_solverInfo["GLOBALSYSSOLN"] == "XxtStaticCond" ||
2016  m_solverInfo["GLOBALSYSSOLN"] ==
2017  "XxtMultiLevelStaticCond" ||
2018  m_solverInfo["GLOBALSYSSOLN"] == "PETScFull" ||
2019  m_solverInfo["GLOBALSYSSOLN"] == "PETScStaticCond" ||
2020  m_solverInfo["GLOBALSYSSOLN"] ==
2021  "PETScMultiLevelStaticCond",
2022  "A parallel solver must be used when run in parallel.");
2023  }
2024  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
SolverInfoMap m_solverInfo
Solver information properties.
CommSharedPtr m_comm
Communication object.
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.
static SolverInfoMap & GetSolverInfoDefaults()
Default solver info options.
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 2217 of file SessionReader.cpp.

References ASSERTL0, and Nektar::StdRegions::find().

2218  {
2219  m_variables.clear();
2220 
2221  if (!conditions)
2222  {
2223  return;
2224  }
2225 
2226  TiXmlElement *variablesElement =
2227  conditions->FirstChildElement("VARIABLES");
2228 
2229  // See if we have parameters defined. They are optional so we go on
2230  // if not.
2231  if (variablesElement)
2232  {
2233  TiXmlElement *varElement =
2234  variablesElement->FirstChildElement("V");
2235 
2236  // Sequential counter for the composite numbers.
2237  int nextVariableNumber = -1;
2238 
2239  while (varElement)
2240  {
2241  stringstream tagcontent;
2242  tagcontent << *varElement;
2243 
2244  /// All elements are of the form: "<V ID="#"> name = value
2245  /// </V>", with ? being the element type.
2246  nextVariableNumber++;
2247 
2248  int i;
2249  int err = varElement->QueryIntAttribute("ID", &i);
2250  ASSERTL0(err == TIXML_SUCCESS,
2251  "Variables must have a unique ID number attribute "
2252  "in XML element: \n\t'" + tagcontent.str() + "'");
2253  ASSERTL0(i == nextVariableNumber,
2254  "ID numbers for variables must begin with zero and"
2255  " be sequential in XML element: \n\t'"
2256  + tagcontent.str() + "'");
2257 
2258  TiXmlNode* varChild = varElement->FirstChild();
2259  // This is primarily to skip comments that may be present.
2260  // Comments appear as nodes just like elements. We are
2261  // specifically looking for text in the body of the
2262  // definition.
2263  while(varChild && varChild->Type() != TiXmlNode::TINYXML_TEXT)
2264  {
2265  varChild = varChild->NextSibling();
2266  }
2267 
2268  ASSERTL0(varChild,
2269  "Unable to read variable definition body for "
2270  "variable with ID "
2271  + boost::lexical_cast<string>(i)
2272  + " in XML element: \n\t'"
2273  + tagcontent.str() + "'");
2274  std::string variableName = varChild->ToText()->ValueStr();
2275 
2276  std::istringstream variableStrm(variableName);
2277  variableStrm >> variableName;
2278 
2279  ASSERTL0(std::find(m_variables.begin(), m_variables.end(),
2280  variableName) == m_variables.end(),
2281  "Variable with ID "
2282  + boost::lexical_cast<string>(i)
2283  + " in XML element \n\t'" + tagcontent.str()
2284  + "'\nhas already been defined.");
2285 
2286  m_variables.push_back(variableName);
2287 
2288  varElement = varElement->NextSiblingElement("V");
2289  }
2290 
2291  ASSERTL0(nextVariableNumber > -1,
2292  "Number of variables must be greater than zero.");
2293  }
2294  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
VariableList m_variables
Variables.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
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 669 of file SessionReader.h.

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

673  {
674  ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
675  CmdLineArg x;
676  x.shortName = pShortName;
677  x.description = pDescription;
678  x.isFlag = false;
679  GetCmdLineArgMap()[pName] = x;
680  return pName;
681  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
static CmdLineArgMap & GetCmdLineArgMap()
CmdLine argument map.
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 687 of file SessionReader.h.

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

Referenced by main().

691  {
692  ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
693  CmdLineArg x;
694  x.shortName = pShortName;
695  x.description = pDescription;
696  x.isFlag = true;
697  GetCmdLineArgMap()[pName] = x;
698  return pName;
699  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
static CmdLineArgMap & GetCmdLineArgMap()
CmdLine argument map.
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");
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 656 of file SessionReader.h.

References GetSolverInfoDefaults().

659  {
660  std::string vName = boost::to_upper_copy(pName);
661  GetSolverInfoDefaults()[vName] = pValue;
662  return pValue;
663  }
static SolverInfoMap & GetSolverInfoDefaults()
Default solver info options.
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",
}
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 623 of file SessionReader.h.

References Nektar::StdRegions::find(), GetSolverInfoEnums(), and Nektar::iterator.

625  {
626  std::string vEnum = boost::to_upper_copy(pEnum);
628  if ((x = GetSolverInfoEnums().find(vEnum)) ==
629  GetSolverInfoEnums().end())
630  {
631  GetSolverInfoEnums()[vEnum] = EnumMap();
632  x = GetSolverInfoEnums().find(vEnum);
633  }
634  x->second[pString] = pEnumValue;
635  return pString;
636  }
std::map< std::string, int > EnumMap
Definition: SessionReader.h:77
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:312
void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
int &  var 
)

Set an integer parameter.

Definition at line 734 of file SessionReader.cpp.

735  {
736  std::string vName = boost::to_upper_copy(pName);
737  m_parameters[vName] = pVar;
738  }
ParameterMap m_parameters
Parameters.
void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
NekDouble var 
)

Set a double precision parameter.

Definition at line 744 of file SessionReader.cpp.

746  {
747  std::string vName = boost::to_upper_copy(pName);
748  m_parameters[vName] = pVar;
749  }
ParameterMap m_parameters
Parameters.
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 782 of file SessionReader.cpp.

References ASSERTL1, and Nektar::iterator.

784  {
785  std::string vProperty = boost::to_upper_copy(pProperty);
786  SolverInfoMap::iterator iter = m_solverInfo.find(vProperty);
787 
788  ASSERTL1(iter != m_solverInfo.end(),
789  "Unable to find requested property: " + pProperty);
790 
791  iter->second = pValue;
792  }
SolverInfoMap m_solverInfo
Solver information properties.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:165
void Nektar::LibUtilities::SessionReader::SetTag ( const std::string &  pName,
const std::string &  pValue 
)

Sets a specified tag.

Definition at line 1311 of file SessionReader.cpp.

1314  {
1315  std::string vName = boost::to_upper_copy(pName);
1316  m_tags[vName] = pValue;
1317  }
void Nektar::LibUtilities::SessionReader::SetUpXmlDoc ( void  )

Definition at line 2634 of file SessionReader.cpp.

2635  {
2637  }
TiXmlDocument * MergeDoc(const std::vector< std::string > &pFilenames) const
Creates an XML document from a list of input files.
std::vector< std::string > m_filenames
Filenames.
TiXmlDocument * m_xmlDoc
Pointer to the loaded XML document.
void Nektar::LibUtilities::SessionReader::SetVariable ( const unsigned int &  idx,
std::string  newname 
)

Definition at line 1021 of file SessionReader.cpp.

References ASSERTL0.

1023  {
1024  ASSERTL0(idx < m_variables.size(), "Variable index out of range.");
1025  m_variables[idx] = newname;
1026  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
VariableList m_variables
Variables.
void Nektar::LibUtilities::SessionReader::SubstituteExpressions ( std::string &  expr)

Substitutes expressions defined in the XML document.

Definition at line 1355 of file SessionReader.cpp.

References Nektar::iterator.

1356  {
1357  ExpressionMap::iterator exprIter;
1358  for (exprIter = m_expressions.begin();
1359  exprIter != m_expressions.end(); ++exprIter)
1360  {
1361  boost::replace_all(pExpr, exprIter->first, exprIter->second);
1362  }
1363  }
ExpressionMap m_expressions
Expressions.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator

Friends And Related Function Documentation

friend class MemoryManager< SessionReader >
friend

Support creation through MemoryManager.

Definition at line 128 of file SessionReader.h.

Member Data Documentation

BndRegionOrdering Nektar::LibUtilities::SessionReader::m_bndRegOrder
private

Map of original boundary region ordering for parallel periodic bcs.

Definition at line 464 of file SessionReader.h.

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

Definition at line 430 of file SessionReader.h.

Referenced by GetCmdLineArgument().

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

Communication object.

Definition at line 433 of file SessionReader.h.

CompositeOrdering Nektar::LibUtilities::SessionReader::m_compOrder
private

Map of original composite ordering for parallel periodic bcs.

Definition at line 461 of file SessionReader.h.

ExpressionMap Nektar::LibUtilities::SessionReader::m_expressions
private

Expressions.

Definition at line 447 of file SessionReader.h.

AnalyticExpressionEvaluator Nektar::LibUtilities::SessionReader::m_exprEvaluator
private

Analytic expression evaluator instance.

Definition at line 449 of file SessionReader.h.

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

Filenames.

Definition at line 435 of file SessionReader.h.

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

Filters map.

Definition at line 457 of file SessionReader.h.

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

Functions.

Definition at line 451 of file SessionReader.h.

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

Geometric information properties.

Definition at line 445 of file SessionReader.h.

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

Parameters.

Definition at line 441 of file SessionReader.h.

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

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

Definition at line 437 of file SessionReader.h.

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

Solver information properties.

Definition at line 443 of file SessionReader.h.

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

Custom tags.

Definition at line 455 of file SessionReader.h.

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

Variables.

Definition at line 453 of file SessionReader.h.

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

Be verbose.

Definition at line 459 of file SessionReader.h.

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

Pointer to the loaded XML document.

Definition at line 439 of file SessionReader.h.