Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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...
 
bool GetSharedFilesystem ()
 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...
 
void TestSharedFilesystem ()
 
std::vector< std::string > ParseCommandLineArguments (int argc, char *argv[])
 Parse the program arguments and fill m_cmdLineOptions. More...
 
std::string ParseSessionName (std::vector< std::string > &filenames)
 Parse the session name. More...
 
void LoadDoc (const std::string &pFilename, TiXmlDocument *pDoc) const
 Loads an xml file into a tinyxml doc and decompresses if needed. More...
 
TiXmlDocument * MergeDoc (const std::vector< std::string > &pFilenames) const
 Creates an XML document from a list of input files. More...
 
void ParseDocument ()
 Loads and parses the specified file. More...
 
void CreateComm (int &argc, char *argv[])
 Loads the given XML document and instantiates an appropriate communication object. More...
 
void 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...
 
bool m_sharedFilesystem
 Running on a shared filesystem. 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 212 of file SessionReader.cpp.

References ASSERTL0.

217  {
218  ASSERTL0(pFilenames.size() > 0, "No filenames specified.");
219 
220  ParseCommandLineArguments(argc, argv);
221  m_xmlDoc = 0;
222  m_filenames = pFilenames;
223 
225 
226  // Create communicator
227  if (!pComm.get())
228  {
229  CreateComm(argc, argv);
230  }
231  else
232  {
233  m_comm = pComm;
234 
235  if (m_comm->GetSize() > 1)
236  {
237  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
238  "IterativeStaticCond";
239  }
240  }
241 
243 
244  // If running in parallel change the default global sys solution
245  // type.
246  if (m_comm->GetSize() > 1)
247  {
248  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
249  "IterativeStaticCond";
250  }
251  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 257 of file SessionReader.cpp.

258  {
259  delete m_xmlDoc;
260  }
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 
198 
199  // If running in parallel change the default global sys solution
200  // type.
201  if (m_comm->GetSize() > 1)
202  {
203  GetSolverInfoDefaults()["GLOBALSYSSOLN"] =
204  "IterativeStaticCond";
205  }
206  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 2622 of file SessionReader.cpp.

References ASSERTL0, and Nektar::lhs.

2623  {
2624  // Parse solver info overrides
2625  if (m_cmdLineOptions.count("solverinfo"))
2626  {
2627  std::vector<std::string> solverInfoList =
2628  m_cmdLineOptions["solverinfo"].as<
2629  std::vector<std::string> >();
2630 
2631  for (int i = 0; i < solverInfoList.size(); ++i)
2632  {
2633  std::string lhs, rhs;
2634 
2635  try
2636  {
2637  ParseEquals(solverInfoList[i], lhs, rhs);
2638  }
2639  catch (...)
2640  {
2641  ASSERTL0(false, "Parse error with command line "
2642  "option: "+solverInfoList[i]);
2643  }
2644 
2645  std::string lhsUpper = boost::to_upper_copy(lhs);
2646  m_solverInfo[lhsUpper] = rhs;
2647  }
2648  }
2649 
2650  if (m_cmdLineOptions.count("parameter"))
2651  {
2652  std::vector<std::string> parametersList =
2653  m_cmdLineOptions["parameter"].as<
2654  std::vector<std::string> >();
2655 
2656  for (int i = 0; i < parametersList.size(); ++i)
2657  {
2658  std::string lhs, rhs;
2659 
2660  try
2661  {
2662  ParseEquals(parametersList[i], lhs, rhs);
2663  }
2664  catch (...)
2665  {
2666  ASSERTL0(false, "Parse error with command line "
2667  "option: "+parametersList[i]);
2668  }
2669 
2670  std::string lhsUpper = boost::to_upper_copy(lhs);
2671 
2672  try
2673  {
2674  m_parameters[lhsUpper] =
2675  boost::lexical_cast<NekDouble>(rhs);
2676  }
2677  catch (...)
2678  {
2679  ASSERTL0(false, "Unable to convert string: "+rhs+
2680  "to double value.");
2681  }
2682  }
2683  }
2684  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1573 of file SessionReader.cpp.

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

1576  {
1577  if (argc == 0)
1578  {
1579  m_comm = GetCommFactory().CreateInstance("Serial", 0, 0);
1580  }
1581  else
1582  {
1583  string vCommModule("Serial");
1584  if (GetCommFactory().ModuleExists("ParallelMPI"))
1585  {
1586  vCommModule = "ParallelMPI";
1587  }
1588 
1589  m_comm = GetCommFactory().CreateInstance(vCommModule,argc,argv);
1590  }
1591  }
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:61
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.

References CellMLToNektar.cellml_metadata::p.

Referenced by main(), Nektar::Utilities::InputMCF::ParseFile(), Nektar::Utilities::InputNekpp::Process(), Nektar::FieldUtils::InputXml::Process(), Nektar::FieldUtils::ProcessDisplacement::Process(), Nektar::FieldUtils::ProcessInterpField::Process(), Nektar::Utilities::OutputNekpp::Process(), Nektar::FieldUtils::ProcessInterpPoints::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:51
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.

References CellMLToNektar.cellml_metadata::p.

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:51
bool Nektar::LibUtilities::SessionReader::DefinesCmdLineArgument ( const std::string &  pName) const

Checks if a specified cmdline argument has been given.

Definition at line 1395 of file SessionReader.cpp.

1397  {
1398  return (m_cmdLineOptions.find(pName) != m_cmdLineOptions.end());
1399  }
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 601 of file SessionReader.cpp.

References ASSERTL0.

602  {
603  std::string vPath = boost::to_upper_copy(pPath);
604  std::vector<std::string> st;
605  boost::split(st, vPath, boost::is_any_of("\\/ "));
606  ASSERTL0(st.size() > 0, "No path given in XML element request.");
607 
608  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
609  ASSERTL0(vReturn, std::string("Cannot find element '")
610  + st[0] + std::string("'."));
611  for (int i = 1; i < st.size(); ++i)
612  {
613  vReturn = vReturn->FirstChildElement(st[i].c_str());
614  if (!vReturn) return false;
615  }
616  return true;
617  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1091 of file SessionReader.cpp.

1092  {
1093  FunctionMap::const_iterator it1;
1094  std::string vName = boost::to_upper_copy(pName);
1095 
1096  if ((it1 = m_functions.find(vName)) != m_functions.end())
1097  {
1098  return true;
1099  }
1100  return false;
1101  }
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 1107 of file SessionReader.cpp.

1111  {
1112  FunctionMap::const_iterator it1;
1113  FunctionVariableMap::const_iterator it2;
1114  std::string vName = boost::to_upper_copy(pName);
1115 
1116  // Check function exists
1117  if ((it1 = m_functions.find(vName)) != m_functions.end())
1118  {
1119  pair<std::string, int> key(pVariable,pDomain);
1120  pair<std::string, int> defkey("*",pDomain);
1121  bool varExists =
1122  (it2 = it1->second.find(key)) != it1->second.end() ||
1123  (it2 = it1->second.find(defkey)) != it1->second.end();
1124  return varExists;
1125  }
1126  return false;
1127  }
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 956 of file SessionReader.cpp.

957  {
958  std::string vName = boost::to_upper_copy(pName);
959  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
960  return (iter != m_geometricInfo.end());
961  }
GeometricInfoMap m_geometricInfo
Geometric information properties.
bool Nektar::LibUtilities::SessionReader::DefinesGlobalSysSolnInfo ( const std::string &  variable,
const std::string &  property 
) const

Definition at line 910 of file SessionReader.cpp.

912  {
913 
914  GloSysSolnInfoList::const_iterator iter =
915  GetGloSysSolnList().find(pVariable);
916  if(iter == GetGloSysSolnList().end())
917  {
918  return false;
919  }
920 
921  std::string vProperty = boost::to_upper_copy(pProperty);
922 
923  GloSysInfoMap::const_iterator iter1 = iter->second.find(vProperty);
924  if(iter1 == iter->second.end())
925  {
926  return false;
927  }
928 
929  return true;
930  }
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 683 of file SessionReader.cpp.

684  {
685  std::string vName = boost::to_upper_copy(pName);
686  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
687  return (paramIter != m_parameters.end());
688  }
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 806 of file SessionReader.cpp.

Referenced by GetSolverInfoAsEnum().

807  {
808  std::string vName = boost::to_upper_copy(pName);
809  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
810  return (infoIter != m_solverInfo.end());
811  }
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 1350 of file SessionReader.cpp.

1351  {
1352  std::string vName = boost::to_upper_copy(pName);
1353  TagMap::const_iterator vTagIterator = m_tags.find(vName);
1354  return (vTagIterator != m_tags.end());
1355  }
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 674 of file SessionReader.cpp.

675  {
676  m_comm->Finalise();
677  }
CommSharedPtr m_comm
Communication object.
BndRegionOrdering Nektar::LibUtilities::SessionReader::GetBndRegionOrdering ( ) const

Definition at line 1420 of file SessionReader.cpp.

1421  {
1422  return m_bndRegOrder;
1423  }
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 406 of file SessionReader.h.

References m_cmdLineOptions.

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

Returns the communication object.

Definition at line 658 of file SessionReader.cpp.

659  {
660  return m_comm;
661  }
CommSharedPtr m_comm
Communication object.
CompositeOrdering Nektar::LibUtilities::SessionReader::GetCompositeOrdering ( ) const

Definition at line 1415 of file SessionReader.cpp.

1416  {
1417  return m_compOrder;
1418  }
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 549 of file SessionReader.cpp.

References ASSERTL1.

550  {
551  ASSERTL1(m_xmlDoc, "XML Document not defined.");
552  return *m_xmlDoc;
553  }
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:228
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 578 of file SessionReader.cpp.

References ASSERTL0.

579  {
580  std::string vPath = boost::to_upper_copy(pPath);
581  std::vector<std::string> st;
582  boost::split(st, vPath, boost::is_any_of("\\/ "));
583  ASSERTL0(st.size() > 0, "No path given in XML element request.");
584 
585  TiXmlElement* vReturn = m_xmlDoc->FirstChildElement(st[0].c_str());
586  ASSERTL0(vReturn, std::string("Cannot find element '")
587  + st[0] + std::string("'."));
588  for (int i = 1; i < st.size(); ++i)
589  {
590  vReturn = vReturn->FirstChildElement(st[i].c_str());
591  ASSERTL0(vReturn, std::string("Cannot find element '")
592  + st[i] + std::string("'."));
593  }
594  return vReturn;
595  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1341 of file SessionReader.cpp.

1342  {
1343  return m_exprEvaluator;
1344  }
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 623 of file SessionReader.cpp.

624  {
625  return m_filenames;
626  }
std::vector< std::string > m_filenames
Filenames.
const FilterMap & Nektar::LibUtilities::SessionReader::GetFilters ( ) const

Definition at line 1386 of file SessionReader.cpp.

1387  {
1388  return m_filters;
1389  }
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 1133 of file SessionReader.cpp.

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

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

References ASSERTL0.

1180  {
1181  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1182  return GetFunction(pName, m_variables[pVar],pDomain);
1183  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1244 of file SessionReader.cpp.

References ASSERTL0.

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

References ASSERTL0.

1290  {
1291  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1292  return GetFunctionFilename(pName, m_variables[pVar],pDomain);
1293  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1299 of file SessionReader.cpp.

References ASSERTL0.

1303  {
1304  FunctionMap::const_iterator it1;
1305  FunctionVariableMap::const_iterator it2, it3;
1306  std::string vName = boost::to_upper_copy(pName);
1307 
1308  it1 = m_functions.find(vName);
1309  ASSERTL0 (it1 != m_functions.end(),
1310  std::string("Function '") + pName
1311  + std::string("' not found."));
1312 
1313  // Check for specific and wildcard definitions
1314  pair<std::string,int> key(pVariable,pDomain);
1315  pair<std::string,int> defkey("*",pDomain);
1316  bool specific = (it2 = it1->second.find(key)) !=
1317  it1->second.end();
1318  bool wildcard = (it3 = it1->second.find(defkey)) !=
1319  it1->second.end();
1320 
1321  // Check function is defined somewhere
1322  ASSERTL0(specific || wildcard,
1323  "No such variable " + pVariable
1324  + " in domain " + boost::lexical_cast<string>(pDomain)
1325  + " defined for function " + pName
1326  + " in session file.");
1327 
1328  // If not specific, must be wildcard
1329  if (!specific)
1330  {
1331  it2 = it3;
1332  }
1333 
1334  return it2->second.m_fileVariable;
1335  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1189 of file SessionReader.cpp.

References ASSERTL0.

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

References ASSERTL0.

1235  {
1236  ASSERTL0(pVar < m_variables.size(), "Variable index out of range.");
1237  return GetFunctionType(pName, m_variables[pVar],pDomain);
1238  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 936 of file SessionReader.cpp.

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

937  {
938  GloSysSolnInfoList::const_iterator iter;
939 
940  ASSERTL0( (iter = GetGloSysSolnList().find(pVariable)) !=
941  GetGloSysSolnList().end(),
942  "Failed to find variable in GlobalSysSolnInfoList");
943 
944  std::string vProperty = boost::to_upper_copy(pProperty);
945  GloSysInfoMap::const_iterator iter1;
946 
947  ASSERTL0( (iter1 = iter->second.find(vProperty)) != iter->second.end(),
948  "Failed to find property: " + vProperty + " in GlobalSysSolnInfoList");
949 
950  return iter1->second;
951  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static GloSysSolnInfoList & GetGloSysSolnList()
GlobalSysSoln Info map.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:316
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 699 of file SessionReader.cpp.

References ASSERTL0.

701  {
702  std::string vName = boost::to_upper_copy(pName);
703  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
704 
705  ASSERTL0(paramIter != m_parameters.end(),
706  "Unable to find requested parameter: " + pName);
707 
708  return paramIter->second;
709  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
const std::string & Nektar::LibUtilities::SessionReader::GetSessionName ( ) const

Returns the session name of the loaded XML document.

Definition at line 632 of file SessionReader.cpp.

633  {
634  return m_sessionName;
635  }
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 642 of file SessionReader.cpp.

References Nektar::LibUtilities::PortablePath().

643  {
644  std::string dirname = m_sessionName + "_xml";
645  fs::path pdirname(dirname);
646 
647  std::string vFilename = "P" + boost::lexical_cast<std::string>(m_comm->GetRowComm()->GetRank());
648  fs::path pFilename(vFilename);
649 
650  fs::path fullpath = pdirname / pFilename;
651 
652  return PortablePath(fullpath);
653  }
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).
bool Nektar::LibUtilities::SessionReader::GetSharedFilesystem ( )

Returns the communication object.

Definition at line 663 of file SessionReader.cpp.

664  {
665  return m_sharedFilesystem;
666  }
bool m_sharedFilesystem
Running on a shared filesystem.
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 714 of file SessionReader.h.

715  {
716  return shared_from_this();
717  }
const std::string & Nektar::LibUtilities::SessionReader::GetSolverInfo ( const std::string &  pProperty) const

Returns the value of the specified solver info property.

Definition at line 817 of file SessionReader.cpp.

References ASSERTL1.

Referenced by GetSolverInfoAsEnum().

819  {
820  std::string vProperty = boost::to_upper_copy(pProperty);
821  SolverInfoMap::const_iterator iter = m_solverInfo.find(vProperty);
822 
823  ASSERTL1(iter != m_solverInfo.end(),
824  "Unable to find requested property: " + pProperty);
825 
826  return iter->second;
827  }
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:228
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 560 of file SessionReader.h.

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

562  {
563  std::string vName = boost::to_upper_copy(pName);
565  "Solver info '" + pName + "' not defined.");
566 
567  std::string vValue = GetSolverInfo(vName);
569  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
570  GetSolverInfoEnums().end(),
571  "Enum for SolverInfo property '" + pName + "' not found.");
572 
574  ASSERTL0((y = x->second.find(vValue)) != x->second.end(),
575  "Value of SolverInfo property '" + pName +
576  "' is invalid.");
577 
578  return T(y->second);
579  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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:316
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 1373 of file SessionReader.cpp.

References ASSERTL0.

1374  {
1375  std::string vName = boost::to_upper_copy(pName);
1376  TagMap::const_iterator vTagIterator = m_tags.find(vName);
1377  ASSERTL0(vTagIterator != m_tags.end(),
1378  "Requested tag does not exist.");
1379  return vTagIterator->second;
1380  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 587 of file SessionReader.h.

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

590  {
591  std::string vName = boost::to_upper_copy(pName);
592 
594  ASSERTL0((x = GetSolverInfoEnums().find(vName)) !=
595  GetSolverInfoEnums().end(),
596  "Enum for property '" + pName + "' not found.");
597 
599  ASSERTL0((y = x->second.find(pValue)) != x->second.end(),
600  "Value of property '" + pValue + "' is invalid.");
601  return T(y->second);
602  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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:316
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 1059 of file SessionReader.cpp.

References ASSERTL0.

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

Returns the names of all variables.

Definition at line 1082 of file SessionReader.cpp.

1083  {
1084  return m_variables;
1085  }
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 269 of file SessionReader.cpp.

References Nektar::iterator.

270  {
271  m_exprEvaluator.SetRandomSeed((m_comm->GetRank() + 1) * time(NULL));
272 
273  // Split up the communicator
274  PartitionComm();
275 
276  // Partition mesh
277  PartitionMesh();
278 
279  // Parse the XML data in #m_xmlDoc
280  ParseDocument();
281 
282  // Override SOLVERINFO and parameters with any specified on the
283  // command line.
284  CmdLineOverride();
285 
286  // Verify SOLVERINFO values
288 
289  // In verbose mode, print out parameters and solver info sections
290  if (m_verbose && m_comm)
291  {
292  if (m_comm->TreatAsRankZero() && m_parameters.size() > 0)
293  {
294  cout << "Parameters:" << endl;
296  for (x = m_parameters.begin(); x != m_parameters.end(); ++x)
297  {
298  cout << "\t" << x->first << " = " << x->second << endl;
299  }
300  cout << endl;
301  }
302 
303  if (m_comm->TreatAsRankZero() && m_solverInfo.size() > 0)
304  {
305  cout << "Solver Info:" << endl;
307  for (x = m_solverInfo.begin(); x != m_solverInfo.end(); ++x)
308  {
309  cout << "\t" << x->first << " = " << x->second << endl;
310  }
311  cout << endl;
312  }
313  }
314  }
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 1428 of file SessionReader.cpp.

References ASSERTL0, CellMLToNektar.pycml::copy(), CellMLToNektar.pycml::format, and Nektar::LibUtilities::PortablePath().

1431  {
1432  if (pFilename.size() > 3 &&
1433  pFilename.substr(pFilename.size() - 3, 3) == ".gz")
1434  {
1435  ifstream file(pFilename.c_str(),
1436  ios_base::in | ios_base::binary);
1437  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1438  stringstream ss;
1439  io::filtering_streambuf<io::input> in;
1440  in.push(io::gzip_decompressor());
1441  in.push(file);
1442  try
1443  {
1444  io::copy(in, ss);
1445  ss >> (*pDoc);
1446  }
1447  catch (io::gzip_error& e)
1448  {
1449  ASSERTL0(false,
1450  "Error: File '" + pFilename + "' is corrupt.");
1451  }
1452  }
1453  else if (pFilename.size() > 4 &&
1454  pFilename.substr(pFilename.size() - 4, 4) == "_xml")
1455  {
1456  fs::path pdirname(pFilename);
1457  boost::format pad("P%1$07d.xml");
1458  pad % m_comm->GetRank();
1459  fs::path pRankFilename(pad.str());
1460  fs::path fullpath = pdirname / pRankFilename;
1461 
1462  ifstream file(PortablePath(fullpath).c_str());
1463  ASSERTL0(file.good(), "Unable to open file: " + fullpath.string());
1464  file >> (*pDoc);
1465  }
1466  else
1467  {
1468  ifstream file(pFilename.c_str());
1469  ASSERTL0(file.good(), "Unable to open file: " + pFilename);
1470  file >> (*pDoc);
1471  }
1472  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 967 of file SessionReader.cpp.

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

992  {
993  std::string vName = boost::to_upper_copy(pName);
994  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
995  if(iter != m_geometricInfo.end())
996  {
997  if (iter->second == "TRUE")
998  {
999  pVar = true;
1000  }
1001  else
1002  {
1003  pVar = false;
1004  }
1005  }
1006  else
1007  {
1008  pVar = pDefault;
1009  }
1010  }
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 1016 of file SessionReader.cpp.

1020  {
1021  std::string vName = boost::to_upper_copy(pName);
1022  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
1023  if(iter != m_geometricInfo.end())
1024  {
1025  pVar = std::atoi(iter->second.c_str());
1026  }
1027  else
1028  {
1029  pVar = pDefault;
1030  }
1031  }
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 715 of file SessionReader.cpp.

References ASSERTL0.

717  {
718  std::string vName = boost::to_upper_copy(pName);
719  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
720  ASSERTL0(paramIter != m_parameters.end(), "Required parameter '" +
721  pName + "' not specified in session.");
722  pVar = (int)round(paramIter->second);
723  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 729 of file SessionReader.cpp.

731  {
732  std::string vName = boost::to_upper_copy(pName);
733  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
734  if(paramIter != m_parameters.end())
735  {
736  pVar = (int)round(paramIter->second);
737  }
738  else
739  {
740  pVar = pDefault;
741  }
742  }
ParameterMap m_parameters
Parameters.
void Nektar::LibUtilities::SessionReader::LoadParameter ( const std::string &  name,
NekDouble var 
) const

Load a double precision parameter.

Definition at line 748 of file SessionReader.cpp.

References ASSERTL0.

750  {
751  std::string vName = boost::to_upper_copy(pName);
752  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
753  ASSERTL0(paramIter != m_parameters.end(), "Required parameter '" +
754  pName + "' not specified in session.");
755  pVar = paramIter->second;
756  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 762 of file SessionReader.cpp.

766  {
767  std::string vName = boost::to_upper_copy(pName);
768  ParameterMap::const_iterator paramIter = m_parameters.find(vName);
769  if(paramIter != m_parameters.end())
770  {
771  pVar = paramIter->second;
772  }
773  else
774  {
775  pVar = pDefault;
776  }
777  }
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 847 of file SessionReader.cpp.

851  {
852  std::string vName = boost::to_upper_copy(pName);
853  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
854  if(infoIter != m_solverInfo.end())
855  {
856  pVar = infoIter->second;
857  }
858  else
859  {
860  pVar = pDefault;
861  }
862  }
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 1037 of file SessionReader.cpp.

1042  {
1043  std::string vName = boost::to_upper_copy(pName);
1044  GeometricInfoMap::const_iterator iter = m_geometricInfo.find(vName);
1045  if(iter != m_geometricInfo.end())
1046  {
1047  pVar = boost::iequals(iter->second, pTrueVal);
1048  }
1049  else
1050  {
1051  pVar = pDefault;
1052  }
1053  }
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 868 of file SessionReader.cpp.

873  {
874  std::string vName = boost::to_upper_copy(pName);
875  SolverInfoMap::const_iterator infoIter = m_solverInfo.find(vName);
876  if(infoIter != m_solverInfo.end())
877  {
878  pVar = boost::iequals(infoIter->second, pTrueVal);
879  }
880  else
881  {
882  pVar = pDefault;
883  }
884  }
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 890 of file SessionReader.cpp.

893  {
894  if (DefinesSolverInfo(pName))
895  {
896  std::string vName = boost::to_upper_copy(pName);
897  SolverInfoMap::const_iterator iter = m_solverInfo.find(vName);
898  if(iter != m_solverInfo.end())
899  {
900  return true;
901  }
902  }
903  return false;
904  }
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 549 of file SessionReader.h.

551  {
552  return (GetSolverInfoAsEnum<T>(name) == trueval);
553  }
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 1477 of file SessionReader.cpp.

References ASSERTL0, CellMLToNektar.cellml_metadata::p, and WARNINGL0.

1479  {
1480  ASSERTL0(pFilenames.size() > 0, "No filenames for merging.");
1481 
1482  // Read the first document
1483  TiXmlDocument *vMainDoc = new TiXmlDocument;
1484  LoadDoc(pFilenames[0], vMainDoc);
1485 
1486  TiXmlHandle vMainHandle(vMainDoc);
1487  TiXmlElement* vMainNektar =
1488  vMainHandle.FirstChildElement("NEKTAR").Element();
1489 
1490  // Read all subsequent XML documents.
1491  // For each element within the NEKTAR tag, use it to replace the
1492  // version already present in the loaded XML data.
1493  for (int i = 1; i < pFilenames.size(); ++i)
1494  {
1495  if((pFilenames[i].compare(pFilenames[i].size()-3,3,"xml") == 0)
1496  ||(pFilenames[i].compare(pFilenames[i].size()-6,6,"xml.gz") == 0))
1497  {
1498  TiXmlDocument* vTempDoc = new TiXmlDocument;
1499  LoadDoc(pFilenames[i], vTempDoc);
1500 
1501  TiXmlHandle docHandle(vTempDoc);
1502  TiXmlElement* vTempNektar;
1503  vTempNektar = docHandle.FirstChildElement("NEKTAR").Element();
1504  ASSERTL0(vTempNektar, "Unable to find NEKTAR tag in file.");
1505  TiXmlElement* p = vTempNektar->FirstChildElement();
1506 
1507  while (p)
1508  {
1509  TiXmlElement *vMainEntry =
1510  vMainNektar->FirstChildElement(p->Value());
1511 
1512  // First check if the new item is in fact blank
1513  if (!p->FirstChild() && vMainEntry)
1514  {
1515  std::string warningmsg =
1516  "File " + pFilenames[i] + " contains " +
1517  "an empty XML element " +
1518  std::string(p->Value()) +
1519  " which will be ignored.";
1520  WARNINGL0(false, warningmsg.c_str());
1521  }
1522  else
1523  {
1524  if (vMainEntry)
1525  {
1526  vMainNektar->RemoveChild(vMainEntry);
1527  }
1528  TiXmlElement *q = new TiXmlElement(*p);
1529  vMainNektar->LinkEndChild(q);
1530  }
1531  p = p->NextSiblingElement();
1532  }
1533 
1534  delete vTempDoc;
1535  }
1536  }
1537  return vMainDoc;
1538  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void LoadDoc(const std::string &pFilename, TiXmlDocument *pDoc) const
Loads an xml file into a tinyxml doc and decompresses if needed.
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:204
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 359 of file SessionReader.cpp.

References Nektar::iterator, Nektar::NekConstants::kGitBranch, Nektar::NekConstants::kGitSha1, NEKTAR_VERSION, CellMLToNektar.cellml_metadata::p, and CellMLToNektar.translators::run().

361  {
362  // List the publically visible options (listed using --help)
363  po::options_description desc("Allowed options");
364  desc.add_options()
365  ("verbose,v", "be verbose")
366  ("version,V", "print version information")
367  ("help,h", "print this help message")
368  ("solverinfo,I", po::value<vector<std::string> >(),
369  "override a SOLVERINFO property")
370  ("parameter,P", po::value<vector<std::string> >(),
371  "override a parameter")
372  ("npx", po::value<int>(),
373  "number of procs in X-dir")
374  ("npy", po::value<int>(),
375  "number of procs in Y-dir")
376  ("npz", po::value<int>(),
377  "number of procs in Z-dir")
378  ("nsz", po::value<int>(),
379  "number of slices in Z-dir")
380  ("part-only", po::value<int>(),
381  "only partition mesh into N partitions.")
382  ("part-only-overlapping", po::value<int>(),
383  "only partition mesh into N overlapping partitions.")
384  ("part-info", "Output partition information")
385  ;
386 
387  CmdLineArgMap::const_iterator cmdIt;
388  for (cmdIt = GetCmdLineArgMap().begin();
389  cmdIt != GetCmdLineArgMap().end(); ++cmdIt)
390  {
391  std::string names = cmdIt->first;
392  if (cmdIt->second.shortName != "")
393  {
394  names += "," + cmdIt->second.shortName;
395  }
396  if (cmdIt->second.isFlag)
397  {
398  desc.add_options()
399  (names.c_str(), cmdIt->second.description.c_str())
400  ;
401  }
402  else
403  {
404  desc.add_options()
405  (names.c_str(), po::value<std::string>(),
406  cmdIt->second.description.c_str())
407  ;
408  }
409  }
410 
411  // List hidden options (e.g. session file arguments are not actually
412  // specified using the input-file option by the user).
413  po::options_description hidden("Hidden options");
414  hidden.add_options()
415  ("input-file", po::value< vector<string> >(),
416  "input filename")
417  ;
418 
419  // Combine all options for the parser
420  po::options_description all("All options");
421  all.add(desc).add(hidden);
422 
423  // Session file is a positional option
424  po::positional_options_description p;
425  p.add("input-file", -1);
426 
427  // Parse the command-line options
428  po::parsed_options parsed = po::command_line_parser(argc, argv).
429  options(all).
430  positional(p).
431  allow_unregistered().
432  run();
433 
434  // Extract known options to map and update
435  po::store(parsed, m_cmdLineOptions);
436  po::notify(m_cmdLineOptions);
437 
438  // Help message
439  if (m_cmdLineOptions.count("help"))
440  {
441  cout << desc;
442  exit(0);
443  }
444 
445  // Version information
446  if (m_cmdLineOptions.count("version"))
447  {
448  cout << "Nektar++ version " << NEKTAR_VERSION;
449 
450  if (NekConstants::kGitSha1 != "GITDIR-NOTFOUND")
451  {
452  string sha1(NekConstants::kGitSha1);
453  string branch(NekConstants::kGitBranch);
454  boost::replace_all(branch, "refs/heads/", "");
455 
456  cout << " (git changeset " << sha1.substr(0, 8) << ", ";
457 
458  if (branch == "")
459  {
460  cout << "detached head";
461  }
462  else
463  {
464  cout << "head " << branch;
465  }
466 
467  cout << ")";
468  }
469 
470  cout << endl;
471  exit(0);
472  }
473 
474  // Enable verbose mode
475  if (m_cmdLineOptions.count("verbose"))
476  {
477  m_verbose = true;
478  }
479  else
480  {
481  m_verbose = false;
482  }
483 
484  // Print a warning for unknown options
485  std::vector< po::basic_option<char> >::iterator x;
486  for (x = parsed.options.begin(); x != parsed.options.end(); ++x)
487  {
488  if (x->unregistered)
489  {
490  cout << "Warning: Unknown option: " << x->string_key
491  << endl;
492  }
493  }
494 
495  // Return the vector of filename(s) given as positional options
496  if (m_cmdLineOptions.count("input-file"))
497  {
498  return m_cmdLineOptions["input-file"].as<
499  std::vector<std::string> >();
500  }
501  else
502  {
503  return std::vector<std::string>();
504  }
505  }
#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 1544 of file SessionReader.cpp.

References ASSERTL0.

1545  {
1546  // Check we actually have a document loaded.
1547  ASSERTL0(m_xmlDoc, "No XML document loaded.");
1548 
1549  // Look for all data in CONDITIONS block.
1550  TiXmlHandle docHandle(m_xmlDoc);
1551  TiXmlElement* e;
1552  e = docHandle.FirstChildElement("NEKTAR").
1553  FirstChildElement("CONDITIONS").Element();
1554 
1555  // Read the various sections of the CONDITIONS block
1556  ReadParameters (e);
1557  ReadSolverInfo (e);
1559  ReadExpressions (e);
1560  ReadVariables (e);
1561  ReadFunctions (e);
1562 
1563  e = docHandle.FirstChildElement("NEKTAR").
1564  FirstChildElement("FILTERS").Element();
1565 
1566  ReadFilters(e);
1567  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 2596 of file SessionReader.cpp.

2600  {
2601  /// Pull out lhs and rhs and eliminate any spaces.
2602  int beg = line.find_first_not_of(" ");
2603  int end = line.find_first_of("=");
2604  // Check for no parameter name
2605  if (beg == end) throw 1;
2606  // Check for no parameter value
2607  if (end != line.find_last_of("=")) throw 1;
2608  // Check for no equals sign
2609  if (end == std::string::npos) throw 1;
2610 
2611  lhs = line.substr(line.find_first_not_of(" "),
2612  end-beg);
2613  lhs = lhs .substr(0, lhs.find_last_not_of(" ")+1);
2614  rhs = line.substr(line.find_last_of("=")+1);
2615  rhs = rhs .substr(rhs.find_first_not_of(" "));
2616  rhs = rhs .substr(0, rhs.find_last_not_of(" ")+1);
2617  }
StandardMatrixTag & lhs
std::string Nektar::LibUtilities::SessionReader::ParseSessionName ( std::vector< std::string > &  filenames)
private

Parse the session name.

Definition at line 511 of file SessionReader.cpp.

References ASSERTL0.

513  {
514  ASSERTL0(!filenames.empty(),
515  "At least one filename expected.");
516 
517  std::string retval = "";
518 
519  // First input file defines the session name
520  std::string fname = filenames[0];
521 
522  // If loading a pre-partitioned mesh, remove _xml extension
523  if (fname.size() > 4 &&
524  fname.substr(fname.size() - 4, 4) == "_xml")
525  {
526  retval = fname.substr(0, fname.find_last_of("_"));
527  }
528  // otherwise remove the .xml extension
529  else if (fname.size() > 4 &&
530  fname.substr(fname.size() - 4, 4) == ".xml")
531  {
532  retval = fname.substr(0, fname.find_last_of("."));
533  }
534  // If compressed .xml.gz, remove both extensions
535  else if (fname.size() > 7 &&
536  fname.substr(fname.size() - 7, 7) == ".xml.gz")
537  {
538  retval = fname.substr(0, fname.find_last_of("."));
539  retval = retval.substr(0, retval.find_last_of("."));
540  }
541 
542  return retval;
543  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 1876 of file SessionReader.cpp.

References ASSERTL0.

1877  {
1878  if (m_comm->GetSize() > 1)
1879  {
1880  int nProcZ = 1;
1881  int nProcY = 1;
1882  int nProcX = 1;
1883  int nStripZ = 1;
1884  if (DefinesCmdLineArgument("npx")) {
1885  nProcX = GetCmdLineArgument<int>("npx");
1886  }
1887  if (DefinesCmdLineArgument("npy")) {
1888  nProcY = GetCmdLineArgument<int>("npy");
1889  }
1890  if (DefinesCmdLineArgument("npz")) {
1891  nProcZ = GetCmdLineArgument<int>("npz");
1892  }
1893  if (DefinesCmdLineArgument("nsz")) {
1894  nStripZ = GetCmdLineArgument<int>("nsz");
1895  }
1896  ASSERTL0(m_comm->GetSize() % (nProcZ*nProcY*nProcX) == 0,
1897  "Cannot exactly partition using PROC_Z value.");
1898  ASSERTL0(nProcZ % nProcY == 0,
1899  "Cannot exactly partition using PROC_Y value.");
1900  ASSERTL0(nProcY % nProcX == 0,
1901  "Cannot exactly partition using PROC_X value.");
1902 
1903  // Number of processes associated with the spectral method
1904  int nProcSm = nProcZ * nProcY * nProcX;
1905 
1906  // Number of processes associated with the spectral element
1907  // method.
1908  int nProcSem = m_comm->GetSize() / nProcSm;
1909 
1910  m_comm->SplitComm(nProcSm,nProcSem);
1911  m_comm->GetColumnComm()->SplitComm(nProcZ/nStripZ,nStripZ);
1912  m_comm->GetColumnComm()->GetColumnComm()->SplitComm(
1913  (nProcY*nProcX),nProcZ/nStripZ);
1914  m_comm->GetColumnComm()->GetColumnComm()->GetColumnComm()
1915  ->SplitComm(nProcX,nProcY);
1916  }
1917  }
bool DefinesCmdLineArgument(const std::string &pName) const
Checks if a specified cmdline argument has been given.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
CommSharedPtr m_comm
Communication object.
void Nektar::LibUtilities::SessionReader::PartitionMesh ( )
private

Partitions the mesh when running in parallel.

Definition at line 1597 of file SessionReader.cpp.

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

1598  {
1599  ASSERTL0(m_comm.get(), "Communication not initialised.");
1600 
1601  // Get row of comm, or the whole comm if not split
1602  CommSharedPtr vCommMesh = m_comm->GetRowComm();
1603  const bool isRoot = m_comm->TreatAsRankZero();
1604 
1605  // Delete any existing loaded mesh
1606  if (m_xmlDoc)
1607  {
1608  delete m_xmlDoc;
1609  }
1610 
1611  // Load file for root process only (since this is always needed)
1612  // and determine if the provided geometry has already been
1613  // partitioned. This will be the case if the user provides the
1614  // directory of mesh partitions as an input. Partitioned geometries
1615  // have the attribute
1616  // PARTITION=X
1617  // where X is the number of the partition (and should match the
1618  // process rank). The result is shared with all other processes.
1619  int isPartitioned = 0;
1620  if (isRoot)
1621  {
1623  if (DefinesElement("Nektar/Geometry"))
1624  {
1625  if (GetElement("Nektar/Geometry")->Attribute("PARTITION"))
1626  {
1627  cout << "Using pre-partitioned mesh." << endl;
1628  isPartitioned = 1;
1629  }
1630  }
1631  }
1632  GetComm()->Bcast(isPartitioned, 0);
1633 
1634  // If the mesh is already partitioned, we are done. Remaining
1635  // processes must load their partitions.
1636  if (isPartitioned)
1637  {
1638  if (!isRoot)
1639  {
1641  }
1642  return;
1643  }
1644 
1645  // Default partitioner to use is Metis. Use Scotch as default
1646  // if it is installed. Override default with command-line flags
1647  // if they are set.
1648  string vPartitionerName = "Metis";
1649  if (GetMeshPartitionFactory().ModuleExists("Scotch"))
1650  {
1651  vPartitionerName = "Scotch";
1652  }
1653  if (DefinesCmdLineArgument("use-metis"))
1654  {
1655  vPartitionerName = "Metis";
1656  }
1657  if (DefinesCmdLineArgument("use-scotch"))
1658  {
1659  vPartitionerName = "Scotch";
1660  }
1661 
1662  // Mesh has not been partitioned so do partitioning if required.
1663  // Note in the serial case nothing is done as we have already loaded
1664  // the mesh.
1665  if (DefinesCmdLineArgument("part-only")||
1666  DefinesCmdLineArgument("part-only-overlapping"))
1667  {
1668  // Perform partitioning of the mesh only. For this we insist
1669  // the code is run in serial (parallel execution is pointless).
1670  ASSERTL0(GetComm()->GetSize() == 1,
1671  "The 'part-only' option should be used in serial.");
1672 
1673  // Number of partitions is specified by the parameter.
1674  int nParts;
1676  MeshPartitionSharedPtr vPartitioner =
1678  vPartitionerName, vSession);
1679  if(DefinesCmdLineArgument("part-only"))
1680  {
1681  nParts = GetCmdLineArgument<int>("part-only");
1682  vPartitioner->PartitionMesh(nParts, true);
1683  }
1684  else
1685  {
1686  nParts = GetCmdLineArgument<int>("part-only-overlapping");
1687  vPartitioner->PartitionMesh(nParts, true, true);
1688  }
1689  vPartitioner->WriteAllPartitions(vSession);
1690  vPartitioner->GetCompositeOrdering(m_compOrder);
1691  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1692 
1693  if (isRoot && DefinesCmdLineArgument("part-info"))
1694  {
1695  vPartitioner->PrintPartInfo(std::cout);
1696  }
1697 
1698  Finalise();
1699  exit(0);
1700  }
1701  else if (vCommMesh->GetSize() > 1)
1702  {
1704  int nParts = vCommMesh->GetSize();
1705  if (m_sharedFilesystem)
1706  {
1707  CommSharedPtr vComm = GetComm();
1708  vector<unsigned int> keys, vals;
1709  int i;
1710 
1711  if (isRoot)
1712  {
1714 
1715  MeshPartitionSharedPtr vPartitioner =
1717  vPartitionerName, vSession);
1718  vPartitioner->PartitionMesh(nParts, true);
1719  vPartitioner->WriteAllPartitions(vSession);
1720  vPartitioner->GetCompositeOrdering(m_compOrder);
1721  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1722 
1723  // Communicate orderings to the other processors.
1724 
1725  // First send sizes of the orderings and boundary
1726  // regions to allocate storage on the remote end.
1727  keys.resize(2);
1728  keys[0] = m_compOrder.size();
1729  keys[1] = m_bndRegOrder.size();
1730  vComm->Bcast(keys, 0);
1731 
1732  // Construct the keys and sizes of values for composite
1733  // ordering
1735  keys.resize(m_compOrder.size());
1736  vals.resize(m_compOrder.size());
1737 
1738  for (cIt = m_compOrder.begin(), i = 0;
1739  cIt != m_compOrder.end(); ++cIt, ++i)
1740  {
1741  keys[i] = cIt->first;
1742  vals[i] = cIt->second.size();
1743  }
1744 
1745  // Send across data.
1746  vComm->Bcast(keys, 0);
1747  vComm->Bcast(vals, 0);
1748  for (cIt = m_compOrder.begin();
1749  cIt != m_compOrder.end(); ++cIt)
1750  {
1751  vComm->Bcast(cIt->second, 0);
1752  }
1753 
1754  // Construct the keys and sizes of values for composite
1755  // ordering
1757  keys.resize(m_bndRegOrder.size());
1758  vals.resize(m_bndRegOrder.size());
1759 
1760  for (bIt = m_bndRegOrder.begin(), i = 0;
1761  bIt != m_bndRegOrder.end(); ++bIt, ++i)
1762  {
1763  keys[i] = bIt->first;
1764  vals[i] = bIt->second.size();
1765  }
1766 
1767  // Send across data.
1768  vComm->Bcast(keys, 0);
1769  vComm->Bcast(vals, 0);
1770  for (bIt = m_bndRegOrder.begin();
1771  bIt != m_bndRegOrder.end(); ++bIt)
1772  {
1773  vComm->Bcast(bIt->second, 0);
1774  }
1775 
1776  if (DefinesCmdLineArgument("part-info"))
1777  {
1778  vPartitioner->PrintPartInfo(std::cout);
1779  }
1780  }
1781  else
1782  {
1783  keys.resize(2);
1784  vComm->Bcast(keys, 0);
1785 
1786  int cmpSize = keys[0];
1787  int bndSize = keys[1];
1788 
1789  keys.resize(cmpSize);
1790  vals.resize(cmpSize);
1791  vComm->Bcast(keys, 0);
1792  vComm->Bcast(vals, 0);
1793 
1794  for (int i = 0; i < keys.size(); ++i)
1795  {
1796  vector<unsigned int> tmp(vals[i]);
1797  vComm->Bcast(tmp, 0);
1798  m_compOrder[keys[i]] = tmp;
1799  }
1800 
1801  keys.resize(bndSize);
1802  vals.resize(bndSize);
1803  vComm->Bcast(keys, 0);
1804  vComm->Bcast(vals, 0);
1805 
1806  for (int i = 0; i < keys.size(); ++i)
1807  {
1808  vector<unsigned int> tmp(vals[i]);
1809  vComm->Bcast(tmp, 0);
1810  m_bndRegOrder[keys[i]] = tmp;
1811  }
1812  }
1813  }
1814  else
1815  {
1816  // Need to load mesh on non-root processes.
1817  if (!isRoot)
1818  {
1820  }
1821 
1822  // Partitioner now operates in parallel
1823  // Each process receives partitioning over interconnect
1824  // and writes its own session file to the working directory.
1825  MeshPartitionSharedPtr vPartitioner =
1827  vPartitionerName, vSession);
1828  vPartitioner->PartitionMesh(nParts, false);
1829  vPartitioner->WriteLocalPartition(vSession);
1830  vPartitioner->GetCompositeOrdering(m_compOrder);
1831  vPartitioner->GetBndRegionOrdering(m_bndRegOrder);
1832 
1833  if (DefinesCmdLineArgument("part-info") && isRoot)
1834  {
1835  vPartitioner->PrintPartInfo(std::cout);
1836  }
1837  }
1838  m_comm->Block();
1839 
1840  std::string dirname = GetSessionName() + "_xml";
1841  fs::path pdirname(dirname);
1842  boost::format pad("P%1$07d.xml");
1843  pad % m_comm->GetRowComm()->GetRank();
1844  fs::path pFilename(pad.str());
1845  fs::path fullpath = pdirname / pFilename;
1846 
1847  std::string vFilename = PortablePath(fullpath);
1848 
1849  if (m_xmlDoc)
1850  {
1851  delete m_xmlDoc;
1852  }
1853  m_xmlDoc = new TiXmlDocument(vFilename);
1854 
1855  ASSERTL0(m_xmlDoc, "Failed to create XML document object.");
1856 
1857  bool loadOkay = m_xmlDoc->LoadFile(vFilename);
1858  ASSERTL0(loadOkay, "Unable to load file: " + vFilename +
1859  ". Check XML standards compliance. Error on line: " +
1860  boost::lexical_cast<std::string>(m_xmlDoc->Row()));
1861  }
1862  else
1863  {
1865  }
1866  }
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:198
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:51
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
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.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
int GetSize(const ConstArray< OneD, NekDouble > &x)
Definition: StdExpUtil.cpp:111
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.
bool m_sharedFilesystem
Running on a shared filesystem.
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 2217 of file SessionReader.cpp.

References ASSERTL0, and Nektar::iterator.

2218  {
2219  m_expressions.clear();
2220 
2221  if (!conditions)
2222  {
2223  return;
2224  }
2225 
2226  TiXmlElement *expressionsElement =
2227  conditions->FirstChildElement("EXPRESSIONS");
2228 
2229  if (expressionsElement)
2230  {
2231  TiXmlElement *expr = expressionsElement->FirstChildElement("E");
2232 
2233  while (expr)
2234  {
2235  stringstream tagcontent;
2236  tagcontent << *expr;
2237  ASSERTL0(expr->Attribute("NAME"),
2238  "Missing NAME attribute in expression "
2239  "definition: \n\t'" + tagcontent.str() + "'");
2240  std::string nameString = expr->Attribute("NAME");
2241  ASSERTL0(!nameString.empty(),
2242  "Expressions must have a non-empty name: \n\t'"
2243  + tagcontent.str() + "'");
2244 
2245  ASSERTL0(expr->Attribute("VALUE"),
2246  "Missing VALUE attribute in expression "
2247  "definition: \n\t'" + tagcontent.str() + "'");
2248  std::string valString = expr->Attribute("VALUE");
2249  ASSERTL0(!valString.empty(),
2250  "Expressions must have a non-empty value: \n\t'"
2251  + tagcontent.str() + "'");
2252 
2253  ExpressionMap::iterator exprIter
2254  = m_expressions.find(nameString);
2255  ASSERTL0(exprIter == m_expressions.end(),
2256  std::string("Expression '") + nameString
2257  + std::string("' already specified."));
2258 
2259  m_expressions[nameString] = valString;
2260  expr = expr->NextSiblingElement("E");
2261  }
2262  }
2263  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 2555 of file SessionReader.cpp.

References ASSERTL0.

2556  {
2557  if (!filters)
2558  {
2559  return;
2560  }
2561 
2562  m_filters.clear();
2563 
2564  TiXmlElement *filter = filters->FirstChildElement("FILTER");
2565  while (filter)
2566  {
2567  ASSERTL0(filter->Attribute("TYPE"),
2568  "Missing attribute 'TYPE' for filter.");
2569  std::string typeStr = filter->Attribute("TYPE");
2570 
2571  std::map<std::string, std::string> vParams;
2572 
2573  TiXmlElement *param = filter->FirstChildElement("PARAM");
2574  while (param)
2575  {
2576  ASSERTL0(param->Attribute("NAME"),
2577  "Missing attribute 'NAME' for parameter in filter "
2578  + typeStr + "'.");
2579  std::string nameStr = param->Attribute("NAME");
2580 
2581  ASSERTL0(param->GetText(), "Empty value string for param.");
2582  std::string valueStr = param->GetText();
2583 
2584  vParams[nameStr] = valueStr;
2585 
2586  param = param->NextSiblingElement("PARAM");
2587  }
2588 
2589  m_filters.push_back(
2590  std::pair<std::string, FilterParams>(typeStr, vParams));
2591 
2592  filter = filter->NextSiblingElement("FILTER");
2593  }
2594  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
FilterMap m_filters
Filters map.
void Nektar::LibUtilities::SessionReader::ReadFunctions ( TiXmlElement *  conditions)
private

Reads the FUNCTIONS section of the XML document.

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

2353  {
2354  m_functions.clear();
2355 
2356  if (!conditions)
2357  {
2358  return;
2359  }
2360 
2361  // Scan through conditions section looking for functions.
2362  TiXmlElement *function = conditions->FirstChildElement("FUNCTION");
2363  while (function)
2364  {
2365  stringstream tagcontent;
2366  tagcontent << *function;
2367 
2368  // Every function must have a NAME attribute
2369  ASSERTL0(function->Attribute("NAME"),
2370  "Functions must have a NAME attribute defined in XML "
2371  "element: \n\t'" + tagcontent.str() + "'");
2372  std::string functionStr = function->Attribute("NAME");
2373  ASSERTL0(!functionStr.empty(),
2374  "Functions must have a non-empty name in XML "
2375  "element: \n\t'" + tagcontent.str() + "'");
2376 
2377  // Store function names in uppercase to remain case-insensitive.
2378  boost::to_upper(functionStr);
2379 
2380  // Retrieve first entry (variable, or file)
2381  TiXmlElement *variable = function->FirstChildElement();
2382 
2383  // Create new function structure with default type of none.
2384  FunctionVariableMap functionVarMap;
2385 
2386  // Process all entries in the function block
2387  while (variable)
2388  {
2389  FunctionVariableDefinition funcDef;
2390  std::string conditionType = variable->Value();
2391 
2392  // If no var is specified, assume wildcard
2393  std::string variableStr;
2394  if (!variable->Attribute("VAR"))
2395  {
2396  variableStr = "*";
2397  }
2398  else
2399  {
2400  variableStr = variable->Attribute("VAR");
2401  }
2402 
2403  // Parse list of variables
2404  std::vector<std::string> variableList;
2405  ParseUtils::GenerateOrderedStringVector(variableStr.c_str(),
2406  variableList);
2407 
2408  // If no domain string put to 0
2409  std::string domainStr;
2410  if (!variable->Attribute("DOMAIN"))
2411  {
2412  domainStr = "0";
2413  }
2414  else
2415  {
2416  domainStr = variable->Attribute("DOMAIN");
2417  }
2418 
2419  // Parse list of variables
2420  std::vector<std::string> varSplit;
2421  std::vector<unsigned int> domainList;
2422  ParseUtils::GenerateSeqVector(domainStr.c_str(), domainList);
2423 
2424  // Expressions are denoted by E
2425  if (conditionType == "E")
2426  {
2427  funcDef.m_type = eFunctionTypeExpression;
2428 
2429  // Expression must have a VALUE.
2430  ASSERTL0(variable->Attribute("VALUE"),
2431  "Attribute VALUE expected for function '"
2432  + functionStr + "'.");
2433  std::string fcnStr = variable->Attribute("VALUE");
2434 
2435  ASSERTL0(!fcnStr.empty(),
2436  (std::string("Expression for var: ")
2437  + variableStr
2438  + std::string(" must be specified.")).c_str());
2439 
2440  SubstituteExpressions(fcnStr);
2441 
2442  // set expression
2443  funcDef.m_expression = MemoryManager<Equation>
2445  }
2446 
2447  // Files are denoted by F
2448  else if (conditionType == "F")
2449  {
2450  if (variable->Attribute("TIMEDEPENDENT") &&
2451  boost::lexical_cast<bool>(variable->Attribute("TIMEDEPENDENT")))
2452  {
2453  funcDef.m_type = eFunctionTypeTransientFile;
2454  }
2455  else
2456  {
2457  funcDef.m_type = eFunctionTypeFile;
2458  }
2459 
2460  // File must have a FILE.
2461  ASSERTL0(variable->Attribute("FILE"),
2462  "Attribute FILE expected for function '"
2463  + functionStr + "'.");
2464  std::string filenameStr = variable->Attribute("FILE");
2465 
2466  ASSERTL0(!filenameStr.empty(),
2467  "A filename must be specified for the FILE "
2468  "attribute of function '" + functionStr
2469  + "'.");
2470 
2471  std::vector<std::string> fSplit;
2472  boost::split(fSplit, filenameStr, boost::is_any_of(":"));
2473 
2474  ASSERTL0(fSplit.size() == 1 || fSplit.size() == 2,
2475  "Incorrect filename specification in function "
2476  + functionStr + "'. "
2477  "Specify variables inside file as: "
2478  "filename:var1,var2");
2479 
2480  // set the filename
2481  funcDef.m_filename = fSplit[0];
2482 
2483  if (fSplit.size() == 2)
2484  {
2485  ASSERTL0(variableList[0] != "*",
2486  "Filename variable mapping not valid "
2487  "when using * as a variable inside "
2488  "function '" + functionStr + "'.");
2489 
2490  boost::split(
2491  varSplit, fSplit[1], boost::is_any_of(","));
2492  ASSERTL0(varSplit.size() == variableList.size(),
2493  "Filename variables should contain the "
2494  "same number of variables defined in "
2495  "VAR in function " + functionStr + "'.");
2496  }
2497  }
2498 
2499  // Nothing else supported so throw an error
2500  else
2501  {
2502  stringstream tagcontent;
2503  tagcontent << *variable;
2504 
2505  ASSERTL0(false,
2506  "Identifier " + conditionType + " in function "
2507  + std::string(function->Attribute("NAME"))
2508  + " is not recognised in XML element: \n\t'"
2509  + tagcontent.str() + "'");
2510  }
2511 
2512 
2513 
2514  // Add variables to function
2515  for (unsigned int i = 0; i < variableList.size(); ++i)
2516  {
2517  for(unsigned int j = 0; j < domainList.size(); ++j)
2518  {
2519  // Check it has not already been defined
2520  pair<std::string,int> key(variableList[i],domainList[j]);
2522  = functionVarMap.find(key);
2523  ASSERTL0(fcnsIter == functionVarMap.end(),
2524  "Error setting expression '" + variableList[i]
2525  + " in domain "
2526  + boost::lexical_cast<std::string>(domainList[j])
2527  + "' in function '" + functionStr + "'. "
2528  "Expression has already been defined.");
2529 
2530  if (varSplit.size() > 0)
2531  {
2532  FunctionVariableDefinition funcDef2 = funcDef;
2533  funcDef2.m_fileVariable = varSplit[i];
2534  functionVarMap[key] = funcDef2;
2535  }
2536  else
2537  {
2538  functionVarMap[key] = funcDef;
2539  }
2540  }
2541  }
2542 
2543  variable = variable->NextSiblingElement();
2544  }
2545  // Add function definition to map
2546  m_functions[functionStr] = functionVarMap;
2547  function = function->NextSiblingElement("FUNCTION");
2548  }
2549  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:143
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:79
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 2083 of file SessionReader.cpp.

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

2084  {
2085  GetGloSysSolnList().clear();
2086 
2087  if (!conditions)
2088  {
2089  return;
2090  }
2091 
2092  TiXmlElement *GlobalSys =
2093  conditions->FirstChildElement("GLOBALSYSSOLNINFO");
2094 
2095  if(!GlobalSys)
2096  {
2097  return;
2098  }
2099 
2100  TiXmlElement *VarInfo = GlobalSys->FirstChildElement("V");
2101 
2102  while (VarInfo)
2103  {
2104  std::stringstream tagcontent;
2105  tagcontent << *VarInfo;
2106  ASSERTL0(VarInfo->Attribute("VAR"),
2107  "Missing VAR attribute in GobalSysSolnInfo XML "
2108  "element: \n\t'" + tagcontent.str() + "'");
2109 
2110  std::string VarList = VarInfo->Attribute("VAR");
2111  ASSERTL0(!VarList.empty(),
2112  "VAR attribute must be non-empty in XML element:\n\t'"
2113  + tagcontent.str() + "'");
2114 
2115  // generate a list of variables.
2116  std::vector<std::string> varStrings;
2118  VarList.c_str(),varStrings);
2119 
2120  ASSERTL0(valid,"Unable to process list of variable in XML "
2121  "element \n\t'" + tagcontent.str() + "'");
2122 
2123  if(varStrings.size())
2124  {
2125  TiXmlElement *SysSolnInfo = VarInfo->FirstChildElement("I");
2126 
2127  while (SysSolnInfo)
2128  {
2129  tagcontent.clear();
2130  tagcontent << *SysSolnInfo;
2131  // read the property name
2132  ASSERTL0(SysSolnInfo->Attribute("PROPERTY"),
2133  "Missing PROPERTY attribute in "
2134  "GlobalSysSolnInfo for variable(s) '"
2135  + VarList + "' in XML element: \n\t'"
2136  + tagcontent.str() + "'");
2137 
2138  std::string SysSolnProperty =
2139  SysSolnInfo->Attribute("PROPERTY");
2140 
2141  ASSERTL0(!SysSolnProperty.empty(),
2142  "GlobalSysSolnIno properties must have a "
2143  "non-empty name for variable(s) : '"
2144  + VarList + "' in XML element: \n\t'"
2145  + tagcontent.str() + "'");
2146 
2147  // make sure that solver property is capitalised
2148  std::string SysSolnPropertyUpper =
2149  boost::to_upper_copy(SysSolnProperty);
2150 
2151  // read the value
2152  ASSERTL0(SysSolnInfo->Attribute("VALUE"),
2153  "Missing VALUE attribute in GlobalSysSolnInfo "
2154  "for variable(s) '" + VarList
2155  + "' in XML element: \n\t"
2156  + tagcontent.str() + "'");
2157 
2158  std::string SysSolnValue =
2159  SysSolnInfo->Attribute("VALUE");
2160  ASSERTL0(!SysSolnValue.empty(),
2161  "GlobalSysSolnInfo properties must have a "
2162  "non-empty value for variable(s) '"
2163  + VarList + "' in XML element: \n\t'"
2164  + tagcontent.str() + "'");
2165 
2166  // Store values under variable map.
2167  for(int i = 0; i < varStrings.size(); ++i)
2168  {
2170  if ((x = GetGloSysSolnList().find(varStrings[i])) ==
2171  GetGloSysSolnList().end())
2172  {
2173  (GetGloSysSolnList()[varStrings[i]])[
2174  SysSolnPropertyUpper] = SysSolnValue;
2175  }
2176  else
2177  {
2178  x->second[SysSolnPropertyUpper] = SysSolnValue;
2179  }
2180  }
2181 
2182  SysSolnInfo = SysSolnInfo->NextSiblingElement("I");
2183  }
2184  VarInfo = VarInfo->NextSiblingElement("V");
2185  }
2186  }
2187 
2188  if (m_verbose && GetGloSysSolnList().size() > 0 && m_comm)
2189  {
2190  if(m_comm->GetRank() == 0)
2191  {
2192  cout << "GlobalSysSoln Info:" << endl;
2193 
2195  for (x = GetGloSysSolnList().begin();
2196  x != GetGloSysSolnList().end();
2197  ++x)
2198  {
2199  cout << "\t Variable: " << x->first << endl;
2200 
2202  for (y = x->second.begin(); y != x->second.end(); ++y)
2203  {
2204  cout << "\t\t " << y->first << " = " << y->second
2205  << endl;
2206  }
2207  }
2208  cout << endl;
2209  }
2210  }
2211  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:143
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:316
void Nektar::LibUtilities::SessionReader::ReadParameters ( TiXmlElement *  conditions)
private

Reads the PARAMETERS section of the XML document.

Definition at line 1923 of file SessionReader.cpp.

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

1924  {
1925  m_parameters.clear();
1926 
1927  if (!conditions)
1928  {
1929  return;
1930  }
1931 
1932  TiXmlElement *parametersElement = conditions->FirstChildElement(
1933  "PARAMETERS");
1934 
1935  // See if we have parameters defined. They are optional so we go on
1936  // if not.
1937  if (parametersElement)
1938  {
1939  TiXmlElement *parameter =
1940  parametersElement->FirstChildElement("P");
1941 
1942  ParameterMap caseSensitiveParameters;
1943 
1944  // Multiple nodes will only occur if there is a comment in
1945  // between definitions.
1946  while (parameter)
1947  {
1948  stringstream tagcontent;
1949  tagcontent << *parameter;
1950  TiXmlNode *node = parameter->FirstChild();
1951 
1952  while (node && node->Type() != TiXmlNode::TINYXML_TEXT)
1953  {
1954  node = node->NextSibling();
1955  }
1956 
1957  if (node)
1958  {
1959  // Format is "paramName = value"
1960  std::string line = node->ToText()->Value(), lhs, rhs;
1961 
1962  try {
1963  ParseEquals(line, lhs, rhs);
1964  }
1965  catch (...)
1966  {
1967  ASSERTL0(false, "Syntax error in parameter "
1968  "expression '" + line
1969  + "' in XML element: \n\t'"
1970  + tagcontent.str() + "'");
1971  }
1972 
1973  // We want the list of parameters to have their RHS
1974  // evaluated, so we use the expression evaluator to do
1975  // the dirty work.
1976  if (!lhs.empty() && !rhs.empty())
1977  {
1978  NekDouble value=0.0;
1979  try
1980  {
1981  LibUtilities::Equation expession(
1982  GetSharedThisPtr(), rhs);
1983  value = expession.Evaluate();
1984  }
1985  catch (const std::runtime_error &)
1986  {
1987  ASSERTL0(false,
1988  "Error evaluating parameter expression"
1989  " '" + rhs + "' in XML element: \n\t'"
1990  + tagcontent.str() + "'");
1991  }
1993  caseSensitiveParameters[lhs] = value;
1994  boost::to_upper(lhs);
1995  m_parameters[lhs] = value;
1996  }
1997  }
1998 
1999  parameter = parameter->NextSiblingElement();
2000  }
2001  }
2002  }
ParameterMap m_parameters
Parameters.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 2008 of file SessionReader.cpp.

References ASSERTL0.

2009  {
2010  m_solverInfo.clear();
2012 
2013  if (!conditions)
2014  {
2015  return;
2016  }
2017 
2018  TiXmlElement *solverInfoElement =
2019  conditions->FirstChildElement("SOLVERINFO");
2020 
2021  if (solverInfoElement)
2022  {
2023  TiXmlElement *solverInfo =
2024  solverInfoElement->FirstChildElement("I");
2025 
2026  while (solverInfo)
2027  {
2028  std::stringstream tagcontent;
2029  tagcontent << *solverInfo;
2030  // read the property name
2031  ASSERTL0(solverInfo->Attribute("PROPERTY"),
2032  "Missing PROPERTY attribute in solver info "
2033  "XML element: \n\t'" + tagcontent.str() + "'");
2034  std::string solverProperty =
2035  solverInfo->Attribute("PROPERTY");
2036  ASSERTL0(!solverProperty.empty(),
2037  "PROPERTY attribute must be non-empty in XML "
2038  "element: \n\t'" + tagcontent.str() + "'");
2039 
2040  // make sure that solver property is capitalised
2041  std::string solverPropertyUpper =
2042  boost::to_upper_copy(solverProperty);
2043 
2044  // read the value
2045  ASSERTL0(solverInfo->Attribute("VALUE"),
2046  "Missing VALUE attribute in solver info "
2047  "XML element: \n\t'" + tagcontent.str() + "'");
2048  std::string solverValue = solverInfo->Attribute("VALUE");
2049  ASSERTL0(!solverValue.empty(),
2050  "VALUE attribute must be non-empty in XML "
2051  "element: \n\t'" + tagcontent.str() + "'");
2052 
2053  // Set Variable
2054  m_solverInfo[solverPropertyUpper] = solverValue;
2055  solverInfo = solverInfo->NextSiblingElement("I");
2056  }
2057  }
2058 
2059  if (m_comm && m_comm->GetRowComm()->GetSize() > 1)
2060  {
2061  ASSERTL0(
2062  m_solverInfo["GLOBALSYSSOLN"] == "IterativeFull" ||
2063  m_solverInfo["GLOBALSYSSOLN"] == "IterativeStaticCond" ||
2064  m_solverInfo["GLOBALSYSSOLN"] ==
2065  "IterativeMultiLevelStaticCond" ||
2066  m_solverInfo["GLOBALSYSSOLN"] == "XxtFull" ||
2067  m_solverInfo["GLOBALSYSSOLN"] == "XxtStaticCond" ||
2068  m_solverInfo["GLOBALSYSSOLN"] ==
2069  "XxtMultiLevelStaticCond" ||
2070  m_solverInfo["GLOBALSYSSOLN"] == "PETScFull" ||
2071  m_solverInfo["GLOBALSYSSOLN"] == "PETScStaticCond" ||
2072  m_solverInfo["GLOBALSYSSOLN"] ==
2073  "PETScMultiLevelStaticCond",
2074  "A parallel solver must be used when run in parallel.");
2075  }
2076  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 2269 of file SessionReader.cpp.

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

2270  {
2271  m_variables.clear();
2272 
2273  if (!conditions)
2274  {
2275  return;
2276  }
2277 
2278  TiXmlElement *variablesElement =
2279  conditions->FirstChildElement("VARIABLES");
2280 
2281  // See if we have parameters defined. They are optional so we go on
2282  // if not.
2283  if (variablesElement)
2284  {
2285  TiXmlElement *varElement =
2286  variablesElement->FirstChildElement("V");
2287 
2288  // Sequential counter for the composite numbers.
2289  int nextVariableNumber = -1;
2290 
2291  while (varElement)
2292  {
2293  stringstream tagcontent;
2294  tagcontent << *varElement;
2295 
2296  /// All elements are of the form: "<V ID="#"> name = value
2297  /// </V>", with ? being the element type.
2298  nextVariableNumber++;
2299 
2300  int i;
2301  int err = varElement->QueryIntAttribute("ID", &i);
2302  ASSERTL0(err == TIXML_SUCCESS,
2303  "Variables must have a unique ID number attribute "
2304  "in XML element: \n\t'" + tagcontent.str() + "'");
2305  ASSERTL0(i == nextVariableNumber,
2306  "ID numbers for variables must begin with zero and"
2307  " be sequential in XML element: \n\t'"
2308  + tagcontent.str() + "'");
2309 
2310  TiXmlNode* varChild = varElement->FirstChild();
2311  // This is primarily to skip comments that may be present.
2312  // Comments appear as nodes just like elements. We are
2313  // specifically looking for text in the body of the
2314  // definition.
2315  while(varChild && varChild->Type() != TiXmlNode::TINYXML_TEXT)
2316  {
2317  varChild = varChild->NextSibling();
2318  }
2319 
2320  ASSERTL0(varChild,
2321  "Unable to read variable definition body for "
2322  "variable with ID "
2323  + boost::lexical_cast<string>(i)
2324  + " in XML element: \n\t'"
2325  + tagcontent.str() + "'");
2326  std::string variableName = varChild->ToText()->ValueStr();
2327 
2328  std::istringstream variableStrm(variableName);
2329  variableStrm >> variableName;
2330 
2331  ASSERTL0(std::find(m_variables.begin(), m_variables.end(),
2332  variableName) == m_variables.end(),
2333  "Variable with ID "
2334  + boost::lexical_cast<string>(i)
2335  + " in XML element \n\t'" + tagcontent.str()
2336  + "'\nhas already been defined.");
2337 
2338  m_variables.push_back(variableName);
2339 
2340  varElement = varElement->NextSiblingElement("V");
2341  }
2342 
2343  ASSERTL0(nextVariableNumber > -1,
2344  "Number of variables must be greater than zero.");
2345  }
2346  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
VariableList m_variables
Variables.
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:316
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 677 of file SessionReader.h.

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

681  {
682  ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
683  CmdLineArg x;
684  x.shortName = pShortName;
685  x.description = pDescription;
686  x.isFlag = false;
687  GetCmdLineArgMap()[pName] = x;
688  return pName;
689  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 695 of file SessionReader.h.

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

699  {
700  ASSERTL0(!pName.empty(), "Empty name for cmdline argument.");
701  CmdLineArg x;
702  x.shortName = pShortName;
703  x.description = pDescription;
704  x.isFlag = true;
705  GetCmdLineArgMap()[pName] = x;
706  return pName;
707  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 664 of file SessionReader.h.

References GetSolverInfoDefaults().

667  {
668  std::string vName = boost::to_upper_copy(pName);
669  GetSolverInfoDefaults()[vName] = pValue;
670  return pValue;
671  }
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:

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

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

633  {
634  std::string vEnum = boost::to_upper_copy(pEnum);
636  if ((x = GetSolverInfoEnums().find(vEnum)) ==
637  GetSolverInfoEnums().end())
638  {
639  GetSolverInfoEnums()[vEnum] = EnumMap();
640  x = GetSolverInfoEnums().find(vEnum);
641  }
642  x->second[pString] = pEnumValue;
643  return pString;
644  }
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:316
void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
int &  var 
)

Set an integer parameter.

Definition at line 784 of file SessionReader.cpp.

785  {
786  std::string vName = boost::to_upper_copy(pName);
787  m_parameters[vName] = pVar;
788  }
ParameterMap m_parameters
Parameters.
void Nektar::LibUtilities::SessionReader::SetParameter ( const std::string &  name,
NekDouble var 
)

Set a double precision parameter.

Definition at line 794 of file SessionReader.cpp.

796  {
797  std::string vName = boost::to_upper_copy(pName);
798  m_parameters[vName] = pVar;
799  }
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 832 of file SessionReader.cpp.

References ASSERTL1, and Nektar::iterator.

834  {
835  std::string vProperty = boost::to_upper_copy(pProperty);
836  SolverInfoMap::iterator iter = m_solverInfo.find(vProperty);
837 
838  ASSERTL1(iter != m_solverInfo.end(),
839  "Unable to find requested property: " + pProperty);
840 
841  iter->second = pValue;
842  }
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:228
void Nektar::LibUtilities::SessionReader::SetTag ( const std::string &  pName,
const std::string &  pValue 
)

Sets a specified tag.

Definition at line 1361 of file SessionReader.cpp.

1364  {
1365  std::string vName = boost::to_upper_copy(pName);
1366  m_tags[vName] = pValue;
1367  }
void Nektar::LibUtilities::SessionReader::SetUpXmlDoc ( void  )

Definition at line 2708 of file SessionReader.cpp.

2709  {
2711  }
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 1071 of file SessionReader.cpp.

References ASSERTL0.

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

Substitutes expressions defined in the XML document.

Definition at line 1405 of file SessionReader.cpp.

References Nektar::iterator.

1406  {
1407  ExpressionMap::iterator exprIter;
1408  for (exprIter = m_expressions.begin();
1409  exprIter != m_expressions.end(); ++exprIter)
1410  {
1411  boost::replace_all(pExpr, exprIter->first, exprIter->second);
1412  }
1413  }
ExpressionMap m_expressions
Expressions.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void Nektar::LibUtilities::SessionReader::TestSharedFilesystem ( )
private

Definition at line 317 of file SessionReader.cpp.

References ASSERTL1, and Nektar::LibUtilities::ReduceSum.

318  {
319  m_sharedFilesystem = false;
320 
321  if (m_comm->GetSize() > 1)
322  {
323  if (m_comm->GetRank() == 0)
324  {
325  std::ofstream testfile("shared-fs-testfile");
326  testfile << "" << std::endl;
327  ASSERTL1(!testfile.fail(), "Test file creation failed");
328  testfile.close();
329  }
330  m_comm->Block();
331 
332  int exists = (bool)boost::filesystem::exists("shared-fs-testfile");
333  m_comm->AllReduce(exists, LibUtilities::ReduceSum);
334 
335  m_sharedFilesystem = (exists == m_comm->GetSize());
336 
337  if ((m_sharedFilesystem && m_comm->GetRank() == 0) ||
339  {
340  std::remove("shared-fs-testfile");
341  }
342  }
343  else
344  {
345  m_sharedFilesystem = false;
346  }
347 
348  if (m_verbose && m_comm->GetRank() == 0 && m_sharedFilesystem)
349  {
350  cout << "Shared filesystem detected" << endl;
351  }
352  }
CommSharedPtr m_comm
Communication object.
bool m_sharedFilesystem
Running on a shared filesystem.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Nektar::LibUtilities::SessionReader::VerifySolverInfo ( )
private

Check values of solver info options are valid.

Definition at line 2686 of file SessionReader.cpp.

References ASSERTL0.

2687  {
2688  SolverInfoMap::const_iterator x;
2689  for (x = m_solverInfo.begin(); x != m_solverInfo.end(); ++x)
2690  {
2691  std::string solverProperty = x->first;
2692  std::string solverValue = x->second;
2693 
2694  EnumMapList::const_iterator propIt =
2695  GetSolverInfoEnums().find(solverProperty);
2696  if (propIt != GetSolverInfoEnums().end())
2697  {
2698  EnumMap::const_iterator valIt =
2699  propIt->second.find(solverValue);
2700  ASSERTL0(valIt != propIt->second.end(),
2701  "Value '" + solverValue + "' is not valid for "
2702  "property '" + solverProperty + "'");
2703  }
2704  }
2705  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 468 of file SessionReader.h.

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

Definition at line 432 of file SessionReader.h.

Referenced by GetCmdLineArgument().

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

Communication object.

Definition at line 435 of file SessionReader.h.

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

Map of original composite ordering for parallel periodic bcs.

Definition at line 465 of file SessionReader.h.

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

Expressions.

Definition at line 449 of file SessionReader.h.

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

Analytic expression evaluator instance.

Definition at line 451 of file SessionReader.h.

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

Filenames.

Definition at line 437 of file SessionReader.h.

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

Filters map.

Definition at line 459 of file SessionReader.h.

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

Functions.

Definition at line 453 of file SessionReader.h.

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

Geometric information properties.

Definition at line 447 of file SessionReader.h.

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

Parameters.

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

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

Running on a shared filesystem.

Definition at line 463 of file SessionReader.h.

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

Solver information properties.

Definition at line 445 of file SessionReader.h.

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

Custom tags.

Definition at line 457 of file SessionReader.h.

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

Variables.

Definition at line 455 of file SessionReader.h.

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

Be verbose.

Definition at line 461 of file SessionReader.h.

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

Pointer to the loaded XML document.

Definition at line 441 of file SessionReader.h.