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 VerifySolverInfo ()
 Check values of solver info options are valid. More...
 
void ParseEquals (const std::string &line, std::string &lhs, std::string &rhs)
 Parse a string in the form lhs = rhs. More...
 

Static Private Member Functions

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

Private Attributes

boost::program_options::variables_map m_cmdLineOptions
 
CommSharedPtr m_comm
 Communication object. More...
 
std::vector< std::string > m_filenames
 Filenames. More...
 
std::string m_sessionName
 Session name of the loaded XML document (filename minus ext). More...
 
TiXmlDocument * m_xmlDoc
 Pointer to the loaded XML document. More...
 
ParameterMap m_parameters
 Parameters. More...
 
SolverInfoMap m_solverInfo
 Solver information properties. More...
 
GeometricInfoMap m_geometricInfo
 Geometric information properties. More...
 
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:161
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:161
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 2562 of file SessionReader.cpp.

References ASSERTL0, and Nektar::lhs.

2563  {
2564  // Parse solver info overrides
2565  if (m_cmdLineOptions.count("solverinfo"))
2566  {
2567  std::vector<std::string> solverInfoList =
2568  m_cmdLineOptions["solverinfo"].as<
2569  std::vector<std::string> >();
2570 
2571  for (int i = 0; i < solverInfoList.size(); ++i)
2572  {
2573  std::string lhs, rhs;
2574 
2575  try
2576  {
2577  ParseEquals(solverInfoList[i], lhs, rhs);
2578  }
2579  catch (...)
2580  {
2581  ASSERTL0(false, "Parse error with command line "
2582  "option: "+solverInfoList[i]);
2583  }
2584 
2585  std::string lhsUpper = boost::to_upper_copy(lhs);
2586  m_solverInfo[lhsUpper] = rhs;
2587  }
2588  }
2589 
2590  if (m_cmdLineOptions.count("parameter"))
2591  {
2592  std::vector<std::string> parametersList =
2593  m_cmdLineOptions["parameter"].as<
2594  std::vector<std::string> >();
2595 
2596  for (int i = 0; i < parametersList.size(); ++i)
2597  {
2598  std::string lhs, rhs;
2599 
2600  try
2601  {
2602  ParseEquals(parametersList[i], lhs, rhs);
2603  }
2604  catch (...)
2605  {
2606  ASSERTL0(false, "Parse error with command line "
2607  "option: "+parametersList[i]);
2608  }
2609 
2610  std::string lhsUpper = boost::to_upper_copy(lhs);
2611 
2612  try
2613  {
2614  m_parameters[lhsUpper] =
2615  boost::lexical_cast<NekDouble>(rhs);
2616  }
2617  catch (...)
2618  {
2619  ASSERTL0(false, "Unable to convert string: "+rhs+
2620  "to double value.");
2621  }
2622  }
2623  }
2624  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1512 of file SessionReader.cpp.

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

1515  {
1516  if (argc == 0)
1517  {
1518  m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1519  }
1520  else
1521  {
1522  string vCommModule("Serial");
1523  if (GetCommFactory().ModuleExists("ParallelMPI"))
1524  {
1525  vCommModule = "ParallelMPI";
1526  }
1527 
1528  m_comm = GetCommFactory().CreateInstance(vCommModule,argc,argv);
1529  }
1530  }
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 1348 of file SessionReader.cpp.

1350  {
1351  return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1352  }
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 558 of file SessionReader.cpp.

References ASSERTL0.

559  {
560  std::string vPath = boost::to_upper_copy(pPath);
561  std::vector<std::string> st;
562  boost::split(st, vPath, boost::is_any_of("\\/ "));
563  ASSERTL0(st.size() > 0, "No path given in XML element request.");
564 
565  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
566  ASSERTL0(vReturn, std::string("Cannot find element '")
567  + st[0] + std::string("'."));
568  for (int i = 1; i < st.size(); ++i)
569  {
570  vReturn = vReturn->FirstChildElement(st[i].c_str());
571  if (!vReturn) return false;
572  }
573  return true;
574  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1044 of file SessionReader.cpp.

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

1064  {
1065  FunctionMap::const_iterator it1;
1066  FunctionVariableMap::const_iterator it2;
1067  std::string vName = boost::to_upper_copy(pName);
1068 
1069  // Check function exists
1070  if ((it1 = m_functions.find(vName)) != m_functions.end())
1071  {
1072  pair<std::string, int> key(pVariable,pDomain);
1073  pair<std::string, int> defkey("*",pDomain);
1074  bool varExists =
1075  (it2 = it1->second.find(key)) != it1->second.end() ||
1076  (it2 = it1->second.find(defkey)) != it1->second.end();
1077  return varExists;
1078  }
1079  return false;
1080  }
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 909 of file SessionReader.cpp.

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

Definition at line 863 of file SessionReader.cpp.

865  {
866 
867  GloSysSolnInfoList::const_iterator iter =
868  GetGloSysSolnList().find(pVariable);
869  if(iter == GetGloSysSolnList().end())
870  {
871  return false;
872  }
873 
874  std::string vProperty = boost::to_upper_copy(pProperty);
875 
876  GloSysInfoMap::const_iterator iter1 = iter->second.find(vProperty);
877  if(iter1 == iter->second.end())
878  {
879  return false;
880  }
881 
882  return true;
883  }
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 636 of file SessionReader.cpp.

637  {
638  std::string vName = boost::to_upper_copy(pName);
639  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
640  return (paramIter != m_parameters.end());
641  }
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 759 of file SessionReader.cpp.

Referenced by GetSolverInfoAsEnum().

760  {
761  std::string vName = boost::to_upper_copy(pName);
762  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
763  return (infoIter != m_solverInfo.end());
764  }
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 1303 of file SessionReader.cpp.

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

628  {
629  m_comm->Finalise();
630  }
CommSharedPtr m_comm
Communication object.
BndRegionOrdering Nektar::LibUtilities::SessionReader::GetBndRegionOrdering ( ) const

Definition at line 1373 of file SessionReader.cpp.

1374  {
1375  return m_bndRegOrder;
1376  }
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 615 of file SessionReader.cpp.

616  {
617  return m_comm;
618  }
CommSharedPtr m_comm
Communication object.
CompositeOrdering Nektar::LibUtilities::SessionReader::GetCompositeOrdering ( ) const

Definition at line 1368 of file SessionReader.cpp.

1369  {
1370  return m_compOrder;
1371  }
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 506 of file SessionReader.cpp.

References ASSERTL1.

507  {
508  ASSERTL1(m_xmlDoc, "XML Document not defined.");
509  return *m_xmlDoc;
510  }
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:191
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 535 of file SessionReader.cpp.

References ASSERTL0.

536  {
537  std::string vPath = boost::to_upper_copy(pPath);
538  std::vector<std::string> st;
539  boost::split(st, vPath, boost::is_any_of("\\/ "));
540  ASSERTL0(st.size() > 0, "No path given in XML element request.");
541 
542  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
543  ASSERTL0(vReturn, std::string("Cannot find element '")
544  + st[0] + std::string("'."));
545  for (int i = 1; i < st.size(); ++i)
546  {
547  vReturn = vReturn->FirstChildElement(st[i].c_str());
548  ASSERTL0(vReturn, std::string("Cannot find element '")
549  + st[i] + std::string("'."));
550  }
551  return vReturn;
552  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1294 of file SessionReader.cpp.

1295  {
1296  return m_exprEvaluator;
1297  }
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 580 of file SessionReader.cpp.

581  {
582  return m_filenames;
583  }
std::vector< std::string > m_filenames
Filenames.
const FilterMap & Nektar::LibUtilities::SessionReader::GetFilters ( ) const

Definition at line 1339 of file SessionReader.cpp.

1340  {
1341  return m_filters;
1342  }
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 1086 of file SessionReader.cpp.

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

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

References ASSERTL0.

1133  {
1134  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1135  return GetFunction(pName, m_variables[pVar],pDomain);
1136  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1197 of file SessionReader.cpp.

References ASSERTL0.

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

References ASSERTL0.

1243  {
1244  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1245  return GetFunctionFilename(pName, m_variables[pVar],pDomain);
1246  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1252 of file SessionReader.cpp.

References ASSERTL0.

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

References ASSERTL0.

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

References ASSERTL0.

1188  {
1189  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1190  return GetFunctionType(pName, m_variables[pVar],pDomain);
1191  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 889 of file SessionReader.cpp.

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

890  {
891  GloSysSolnInfoList::const_iterator iter;
892 
893  ASSERTL0( (iter = GetGloSysSolnList().find(pVariable)) !=
894  GetGloSysSolnList().end(),
895  "Failed to find variable in GlobalSysSolnInfoList");
896 
897  std::string vProperty = boost::to_upper_copy(pProperty);
898  GloSysInfoMap::const_iterator iter1;
899 
900  ASSERTL0( (iter1 = iter->second.find(vProperty)) != iter->second.end(),
901  "Failed to find property: " + vProperty + " in GlobalSysSolnInfoList");
902 
903  return iter1->second;
904  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 652 of file SessionReader.cpp.

References ASSERTL0.

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

Returns the session name of the loaded XML document.

Definition at line 589 of file SessionReader.cpp.

590  {
591  return m_sessionName;
592  }
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 599 of file SessionReader.cpp.

References Nektar::LibUtilities::PortablePath().

600  {
601  std::string dirname = m_sessionName + "_xml";
602  fs::path pdirname(dirname);
603 
604  std::string vFilename = "P" + boost::lexical_cast<std::string>(m_comm->GetRowComm()->GetRank());
605  fs::path pFilename(vFilename);
606 
607  fs::path fullpath = pdirname / pFilename;
608 
609  return PortablePath(fullpath);
610  }
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 708 of file SessionReader.h.

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

Returns the value of the specified solver info property.

Definition at line 770 of file SessionReader.cpp.

References ASSERTL1.

Referenced by GetSolverInfoAsEnum().

772  {
773  std::string vProperty = boost::to_upper_copy(pProperty);
774  SolverInfoMap::const_iterator iter = m_solverInfo.find(vProperty);
775 
776  ASSERTL1(iter != m_solverInfo.end(),
777  "Unable to find requested property: " + pProperty);
778 
779  return iter->second;
780  }
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:191
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 554 of file SessionReader.h.

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

556  {
557  std::string vName = boost::to_upper_copy(pName);
559  "Solver info '" + pName + "' not defined.");
560 
561  std::string vValue = GetSolverInfo(vName);
563  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
564  GetSolverInfoEnums().end(),
565  "Enum for SolverInfo property '" + pName + "' not found.");
566 
568  ASSERTL0((y = x->second.find(vValue)) != x->second.end(),
569  "Value of SolverInfo property '" + pName +
570  "' is invalid.");
571 
572  return T(y->second);
573  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1326 of file SessionReader.cpp.

References ASSERTL0.

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

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

584  {
585  std::string vName = boost::to_upper_copy(pName);
586 
588  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
589  GetSolverInfoEnums().end(),
590  "Enum for property '" + pName + "' not found.");
591 
593  ASSERTL0((y = x->second.find(pValue)) != x->second.end(),
594  "Value of property '" + pValue + "' is invalid.");
595  return T(y->second);
596  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1012 of file SessionReader.cpp.

References ASSERTL0.

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

Returns the names of all variables.

Definition at line 1035 of file SessionReader.cpp.

1036  {
1037  return m_variables;
1038  }
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  // Verify SOLVERINFO values
284 
285  // In verbose mode, print out parameters and solver info sections
286  if (m_verbose && m_comm)
287  {
288  if (m_comm->GetRank() == 0 && m_parameters.size() > 0)
289  {
290  cout << "Parameters:" << endl;
292  for (x = m_parameters.begin(); x != m_parameters.end(); ++x)
293  {
294  cout << "\t" << x->first << " = " << x->second << endl;
295  }
296  cout << endl;
297  }
298 
299  if (m_comm->GetRank() == 0 && m_solverInfo.size() > 0)
300  {
301  cout << "Solver Info:" << endl;
303  for (x = m_solverInfo.begin(); x != m_solverInfo.end(); ++x)
304  {
305  cout << "\t" << x->first << " = " << x->second << endl;
306  }
307  cout << endl;
308  }
309  }
310  }
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.
void VerifySolverInfo()
Check values of solver info options are valid.
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 1381 of file SessionReader.cpp.

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

1384  {
1385  if (pFilename.size() > 3 &&
1386  pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1387  {
1388  ifstream file(pFilename.c_str(),
1389  ios_base::in | ios_base::binary);
1390  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1391  stringstream ss;
1392  io::filtering_streambuf<io::input> in;
1393  in.push(io::gzip_decompressor());
1394  in.push(file);
1395  try
1396  {
1397  io::copy(in, ss);
1398  ss >> (*pDoc);
1399  }
1400  catch (io::gzip_error& e)
1401  {
1402  ASSERTL0(false,
1403  "Error: File '" + pFilename + "' is corrupt.");
1404  }
1405  }
1406  else if (pFilename.size() > 4 &&
1407  pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1408  {
1409  fs::path pdirname(pFilename);
1410  boost::format pad("P%1$07d.xml");
1411  pad % m_comm->GetRank();
1412  fs::path pRankFilename(pad.str());
1413  fs::path fullpath = pdirname / pRankFilename;
1414 
1415  ifstream file(PortablePath(fullpath).c_str());
1416  ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1417  file >> (*pDoc);
1418  }
1419  else
1420  {
1421  ifstream file(pFilename.c_str());
1422  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1423  file >> (*pDoc);
1424  }
1425  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 920 of file SessionReader.cpp.

924  {
925  std::string vName = boost::to_upper_copy(pName);
926  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
927  if(iter != m_geometricInfo.end())
928  {
929  pVar = iter->second;
930  }
931  else
932  {
933  pVar = pDefault;
934  }
935  }
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 941 of file SessionReader.cpp.

945  {
946  std::string vName = boost::to_upper_copy(pName);
947  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
948  if(iter != m_geometricInfo.end())
949  {
950  if (iter->second == "TRUE")
951  {
952  pVar = true;
953  }
954  else
955  {
956  pVar = false;
957  }
958  }
959  else
960  {
961  pVar = pDefault;
962  }
963  }
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 969 of file SessionReader.cpp.

973  {
974  std::string vName = boost::to_upper_copy(pName);
975  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
976  if(iter != m_geometricInfo.end())
977  {
978  pVar = std::atoi(iter->second.c_str());
979  }
980  else
981  {
982  pVar = pDefault;
983  }
984  }
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 668 of file SessionReader.cpp.

References ASSERTL0.

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

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

Load a double precision parameter.

Definition at line 701 of file SessionReader.cpp.

References ASSERTL0.

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

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

804  {
805  std::string vName = boost::to_upper_copy(pName);
806  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
807  if(infoIter != m_solverInfo.end())
808  {
809  pVar = infoIter->second;
810  }
811  else
812  {
813  pVar = pDefault;
814  }
815  }
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 990 of file SessionReader.cpp.

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

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

846  {
847  if (DefinesSolverInfo(pName))
848  {
849  std::string vName = boost::to_upper_copy(pName);
850  SolverInfoMap::const_iterator iter = m_solverInfo.find(vName);
851  if(iter != m_solverInfo.end())
852  {
853  return true;
854  }
855  }
856  return false;
857  }
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 543 of file SessionReader.h.

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

References ASSERTL0.

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

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

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

References ASSERTL0.

1484  {
1485  // Check we actually have a document loaded.
1486  ASSERTL0(m_xmlDoc, "No XML document loaded.");
1487 
1488  // Look for all data in CONDITIONS block.
1489  TiXmlHandle docHandle(m_xmlDoc);
1490  TiXmlElement* e;
1491  e = docHandle.FirstChildElement("NEKTAR").
1492  FirstChildElement("CONDITIONS").Element();
1493 
1494  // Read the various sections of the CONDITIONS block
1495  ReadParameters (e);
1496  ReadSolverInfo (e);
1498  ReadExpressions (e);
1499  ReadVariables (e);
1500  ReadFunctions (e);
1501 
1502  e = docHandle.FirstChildElement("NEKTAR").
1503  FirstChildElement("FILTERS").Element();
1504 
1505  ReadFilters(e);
1506  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 2536 of file SessionReader.cpp.

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

Parse the session name.

Definition at line 468 of file SessionReader.cpp.

References ASSERTL0.

470  {
471  ASSERTL0(!filenames.empty(),
472  "At least one filename expected.");
473 
474  std::string retval = "";
475 
476  // First input file defines the session name
477  std::string fname = filenames[0];
478 
479  // If loading a pre-partitioned mesh, remove _xml extension
480  if (fname.size() > 4 &&
481  fname.substr(fname.size() - 4, 4) == "_xml")
482  {
483  retval = fname.substr(0, fname.find_last_of("_"));
484  }
485  // otherwise remove the .xml extension
486  else if (fname.size() > 4 &&
487  fname.substr(fname.size() - 4, 4) == ".xml")
488  {
489  retval = fname.substr(0, fname.find_last_of("."));
490  }
491  // If compressed .xml.gz, remove both extensions
492  else if (fname.size() > 7 &&
493  fname.substr(fname.size() - 7, 7) == ".xml.gz")
494  {
495  retval = fname.substr(0, fname.find_last_of("."));
496  retval = retval.substr(0, retval.find_last_of("."));
497  }
498 
499  return retval;
500  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1816 of file SessionReader.cpp.

References ASSERTL0.

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

Partitions the mesh when running in parallel.

Definition at line 1536 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.

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

References ASSERTL0, and Nektar::iterator.

2158  {
2159  m_expressions.clear();
2160 
2161  if (!conditions)
2162  {
2163  return;
2164  }
2165 
2166  TiXmlElement *expressionsElement =
2167  conditions->FirstChildElement("EXPRESSIONS");
2168 
2169  if (expressionsElement)
2170  {
2171  TiXmlElement *expr = expressionsElement->FirstChildElement("E");
2172 
2173  while (expr)
2174  {
2175  stringstream tagcontent;
2176  tagcontent << *expr;
2177  ASSERTL0(expr->Attribute("NAME"),
2178  "Missing NAME attribute in expression "
2179  "definition: \n\t'" + tagcontent.str() + "'");
2180  std::string nameString = expr->Attribute("NAME");
2181  ASSERTL0(!nameString.empty(),
2182  "Expressions must have a non-empty name: \n\t'"
2183  + tagcontent.str() + "'");
2184 
2185  ASSERTL0(expr->Attribute("VALUE"),
2186  "Missing VALUE attribute in expression "
2187  "definition: \n\t'" + tagcontent.str() + "'");
2188  std::string valString = expr->Attribute("VALUE");
2189  ASSERTL0(!valString.empty(),
2190  "Expressions must have a non-empty value: \n\t'"
2191  + tagcontent.str() + "'");
2192 
2193  ExpressionMap::iterator exprIter
2194  = m_expressions.find(nameString);
2195  ASSERTL0(exprIter == m_expressions.end(),
2196  std::string("Expression '") + nameString
2197  + std::string("' already specified."));
2198 
2199  m_expressions[nameString] = valString;
2200  expr = expr->NextSiblingElement("E");
2201  }
2202  }
2203  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 2495 of file SessionReader.cpp.

References ASSERTL0.

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

Reads the FUNCTIONS section of the XML document.

Definition at line 2292 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.

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

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

2024  {
2025  GetGloSysSolnList().clear();
2026 
2027  if (!conditions)
2028  {
2029  return;
2030  }
2031 
2032  TiXmlElement *GlobalSys =
2033  conditions->FirstChildElement("GLOBALSYSSOLNINFO");
2034 
2035  if(!GlobalSys)
2036  {
2037  return;
2038  }
2039 
2040  TiXmlElement *VarInfo = GlobalSys->FirstChildElement("V");
2041 
2042  while (VarInfo)
2043  {
2044  std::stringstream tagcontent;
2045  tagcontent << *VarInfo;
2046  ASSERTL0(VarInfo->Attribute("VAR"),
2047  "Missing VAR attribute in GobalSysSolnInfo XML "
2048  "element: \n\t'" + tagcontent.str() + "'");
2049 
2050  std::string VarList = VarInfo->Attribute("VAR");
2051  ASSERTL0(!VarList.empty(),
2052  "VAR attribute must be non-empty in XML element:\n\t'"
2053  + tagcontent.str() + "'");
2054 
2055  // generate a list of variables.
2056  std::vector<std::string> varStrings;
2058  VarList.c_str(),varStrings);
2059 
2060  ASSERTL0(valid,"Unable to process list of variable in XML "
2061  "element \n\t'" + tagcontent.str() + "'");
2062 
2063  if(varStrings.size())
2064  {
2065  TiXmlElement *SysSolnInfo = VarInfo->FirstChildElement("I");
2066 
2067  while (SysSolnInfo)
2068  {
2069  tagcontent.clear();
2070  tagcontent << *SysSolnInfo;
2071  // read the property name
2072  ASSERTL0(SysSolnInfo->Attribute("PROPERTY"),
2073  "Missing PROPERTY attribute in "
2074  "GlobalSysSolnInfo for variable(s) '"
2075  + VarList + "' in XML element: \n\t'"
2076  + tagcontent.str() + "'");
2077 
2078  std::string SysSolnProperty =
2079  SysSolnInfo->Attribute("PROPERTY");
2080 
2081  ASSERTL0(!SysSolnProperty.empty(),
2082  "GlobalSysSolnIno properties must have a "
2083  "non-empty name for variable(s) : '"
2084  + VarList + "' in XML element: \n\t'"
2085  + tagcontent.str() + "'");
2086 
2087  // make sure that solver property is capitalised
2088  std::string SysSolnPropertyUpper =
2089  boost::to_upper_copy(SysSolnProperty);
2090 
2091  // read the value
2092  ASSERTL0(SysSolnInfo->Attribute("VALUE"),
2093  "Missing VALUE attribute in GlobalSysSolnInfo "
2094  "for variable(s) '" + VarList
2095  + "' in XML element: \n\t"
2096  + tagcontent.str() + "'");
2097 
2098  std::string SysSolnValue =
2099  SysSolnInfo->Attribute("VALUE");
2100  ASSERTL0(!SysSolnValue.empty(),
2101  "GlobalSysSolnInfo properties must have a "
2102  "non-empty value for variable(s) '"
2103  + VarList + "' in XML element: \n\t'"
2104  + tagcontent.str() + "'");
2105 
2106  // Store values under variable map.
2107  for(int i = 0; i < varStrings.size(); ++i)
2108  {
2110  if ((x = GetGloSysSolnList().find(varStrings[i])) ==
2111  GetGloSysSolnList().end())
2112  {
2113  (GetGloSysSolnList()[varStrings[i]])[
2114  SysSolnPropertyUpper] = SysSolnValue;
2115  }
2116  else
2117  {
2118  x->second[SysSolnPropertyUpper] = SysSolnValue;
2119  }
2120  }
2121 
2122  SysSolnInfo = SysSolnInfo->NextSiblingElement("I");
2123  }
2124  VarInfo = VarInfo->NextSiblingElement("V");
2125  }
2126  }
2127 
2128  if (m_verbose && GetGloSysSolnList().size() > 0 && m_comm)
2129  {
2130  if(m_comm->GetRank() == 0)
2131  {
2132  cout << "GlobalSysSoln Info:" << endl;
2133 
2135  for (x = GetGloSysSolnList().begin();
2136  x != GetGloSysSolnList().end();
2137  ++x)
2138  {
2139  cout << "\t Variable: " << x->first << endl;
2140 
2142  for (y = x->second.begin(); y != x->second.end(); ++y)
2143  {
2144  cout << "\t\t " << y->first << " = " << y->second
2145  << endl;
2146  }
2147  }
2148  cout << endl;
2149  }
2150  }
2151  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1863 of file SessionReader.cpp.

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

1864  {
1865  m_parameters.clear();
1866 
1867  if (!conditions)
1868  {
1869  return;
1870  }
1871 
1872  TiXmlElement *parametersElement = conditions->FirstChildElement(
1873  "PARAMETERS");
1874 
1875  // See if we have parameters defined. They are optional so we go on
1876  // if not.
1877  if (parametersElement)
1878  {
1879  TiXmlElement *parameter =
1880  parametersElement->FirstChildElement("P");
1881 
1882  ParameterMap caseSensitiveParameters;
1883 
1884  // Multiple nodes will only occur if there is a comment in
1885  // between definitions.
1886  while (parameter)
1887  {
1888  stringstream tagcontent;
1889  tagcontent << *parameter;
1890  TiXmlNode *node = parameter->FirstChild();
1891 
1892  while (node && node->Type() != TiXmlNode::TINYXML_TEXT)
1893  {
1894  node = node->NextSibling();
1895  }
1896 
1897  if (node)
1898  {
1899  // Format is "paramName = value"
1900  std::string line = node->ToText()->Value(), lhs, rhs;
1901 
1902  try {
1903  ParseEquals(line, lhs, rhs);
1904  }
1905  catch (...)
1906  {
1907  ASSERTL0(false, "Syntax error in parameter "
1908  "expression '" + line
1909  + "' in XML element: \n\t'"
1910  + tagcontent.str() + "'");
1911  }
1912 
1913  // We want the list of parameters to have their RHS
1914  // evaluated, so we use the expression evaluator to do
1915  // the dirty work.
1916  if (!lhs.empty() && !rhs.empty())
1917  {
1918  NekDouble value=0.0;
1919  try
1920  {
1921  LibUtilities::Equation expession(
1922  GetSharedThisPtr(), rhs);
1923  value = expession.Evaluate();
1924  }
1925  catch (const std::runtime_error &)
1926  {
1927  ASSERTL0(false,
1928  "Error evaluating parameter expression"
1929  " '" + rhs + "' in XML element: \n\t'"
1930  + tagcontent.str() + "'");
1931  }
1933  caseSensitiveParameters[lhs] = value;
1934  boost::to_upper(lhs);
1935  m_parameters[lhs] = value;
1936  }
1937  }
1938 
1939  parameter = parameter->NextSiblingElement();
1940  }
1941  }
1942  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 1948 of file SessionReader.cpp.

References ASSERTL0.

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

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

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

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

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

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

Referenced by main().

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

References GetSolverInfoDefaults().

661  {
662  std::string vName = boost::to_upper_copy(pName);
663  GetSolverInfoDefaults()[vName] = pValue;
664  return pValue;
665  }
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 625 of file SessionReader.h.

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

627  {
628  std::string vEnum = boost::to_upper_copy(pEnum);
630  if ((x = GetSolverInfoEnums().find(vEnum)) ==
631  GetSolverInfoEnums().end())
632  {
633  GetSolverInfoEnums()[vEnum] = EnumMap();
634  x = GetSolverInfoEnums().find(vEnum);
635  }
636  x->second[pString] = pEnumValue;
637  return pString;
638  }
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 737 of file SessionReader.cpp.

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

Set a double precision parameter.

Definition at line 747 of file SessionReader.cpp.

749  {
750  std::string vName = boost::to_upper_copy(pName);
751  m_parameters[vName] = pVar;
752  }
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 785 of file SessionReader.cpp.

References ASSERTL1, and Nektar::iterator.

787  {
788  std::string vProperty = boost::to_upper_copy(pProperty);
789  SolverInfoMap::iterator iter = m_solverInfo.find(vProperty);
790 
791  ASSERTL1(iter != m_solverInfo.end(),
792  "Unable to find requested property: " + pProperty);
793 
794  iter->second = pValue;
795  }
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:191
void Nektar::LibUtilities::SessionReader::SetTag ( const std::string &  pName,
const std::string &  pValue 
)

Sets a specified tag.

Definition at line 1314 of file SessionReader.cpp.

1317  {
1318  std::string vName = boost::to_upper_copy(pName);
1319  m_tags[vName] = pValue;
1320  }
void Nektar::LibUtilities::SessionReader::SetUpXmlDoc ( void  )

Definition at line 2648 of file SessionReader.cpp.

2649  {
2651  }
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 1024 of file SessionReader.cpp.

References ASSERTL0.

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

Substitutes expressions defined in the XML document.

Definition at line 1358 of file SessionReader.cpp.

References Nektar::iterator.

1359  {
1360  ExpressionMap::iterator exprIter;
1361  for (exprIter = m_expressions.begin();
1362  exprIter != m_expressions.end(); ++exprIter)
1363  {
1364  boost::replace_all(pExpr, exprIter->first, exprIter->second);
1365  }
1366  }
ExpressionMap m_expressions
Expressions.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void Nektar::LibUtilities::SessionReader::VerifySolverInfo ( )
private

Check values of solver info options are valid.

Definition at line 2626 of file SessionReader.cpp.

References ASSERTL0.

2627  {
2628  SolverInfoMap::const_iterator x;
2629  for (x = m_solverInfo.begin(); x != m_solverInfo.end(); ++x)
2630  {
2631  std::string solverProperty = x->first;
2632  std::string solverValue = x->second;
2633 
2634  EnumMapList::const_iterator propIt =
2635  GetSolverInfoEnums().find(solverProperty);
2636  if (propIt != GetSolverInfoEnums().end())
2637  {
2638  EnumMap::const_iterator valIt =
2639  propIt->second.find(solverValue);
2640  ASSERTL0(valIt != propIt->second.end(),
2641  "Value '" + solverValue + "' is not valid for "
2642  "property '" + solverProperty + "'");
2643  }
2644  }
2645  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
SolverInfoMap m_solverInfo
Solver information properties.
static EnumMapList & GetSolverInfoEnums()
String to enumeration map for Solver Info parameters.

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.