added: support any number of scalar solutions in anasol
necessary for e.g. boussinesq where you have a velocity, a pressure and a temperature field.
This commit is contained in:
@@ -18,8 +18,17 @@
|
|||||||
#include "tinyxml.h"
|
#include "tinyxml.h"
|
||||||
|
|
||||||
|
|
||||||
|
AnaSol::AnaSol(RealFunc* s1, VecFunc* s2, VecFunc* v1,
|
||||||
|
TensorFunc* v2, STensorFunc* v3) :
|
||||||
|
vecSol(v1), vecSecSol(v2), stressSol(v3)
|
||||||
|
{
|
||||||
|
if (s1) scalSol.push_back(s1);
|
||||||
|
if (s2) scalSecSol.push_back(s2);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
AnaSol::AnaSol (std::istream& is, const int nlines, bool scalarSol)
|
AnaSol::AnaSol (std::istream& is, const int nlines, bool scalarSol)
|
||||||
: scalSol(0), scalSecSol(0), vecSol(0), vecSecSol(0), stressSol(0)
|
: vecSol(0), vecSecSol(0), stressSol(0)
|
||||||
{
|
{
|
||||||
size_t pos = 0;
|
size_t pos = 0;
|
||||||
std::string variables;
|
std::string variables;
|
||||||
@@ -39,7 +48,7 @@ AnaSol::AnaSol (std::istream& is, const int nlines, bool scalarSol)
|
|||||||
std::string primary = function.substr(pos+8);
|
std::string primary = function.substr(pos+8);
|
||||||
IFEM::cout <<"\tPrimary="<< primary << std::endl;
|
IFEM::cout <<"\tPrimary="<< primary << std::endl;
|
||||||
if (scalarSol)
|
if (scalarSol)
|
||||||
scalSol = new EvalFunction((variables+primary).c_str());
|
scalSol.push_back(new EvalFunction((variables+primary).c_str()));
|
||||||
else
|
else
|
||||||
vecSol = new VecFuncExpr(primary,variables);
|
vecSol = new VecFuncExpr(primary,variables);
|
||||||
}
|
}
|
||||||
@@ -49,7 +58,7 @@ AnaSol::AnaSol (std::istream& is, const int nlines, bool scalarSol)
|
|||||||
std::string secondary = function.substr(pos+10);
|
std::string secondary = function.substr(pos+10);
|
||||||
IFEM::cout <<"\tSecondary="<< secondary << std::endl;
|
IFEM::cout <<"\tSecondary="<< secondary << std::endl;
|
||||||
if (scalarSol)
|
if (scalarSol)
|
||||||
scalSecSol = new VecFuncExpr(secondary,variables);
|
scalSecSol.push_back(new VecFuncExpr(secondary,variables));
|
||||||
else
|
else
|
||||||
vecSecSol = new TensorFuncExpr(secondary,variables);
|
vecSecSol = new TensorFuncExpr(secondary,variables);
|
||||||
}
|
}
|
||||||
@@ -65,7 +74,7 @@ AnaSol::AnaSol (std::istream& is, const int nlines, bool scalarSol)
|
|||||||
|
|
||||||
|
|
||||||
AnaSol::AnaSol (const TiXmlElement* elem, bool scalarSol)
|
AnaSol::AnaSol (const TiXmlElement* elem, bool scalarSol)
|
||||||
: scalSol(0), scalSecSol(0), vecSol(0), vecSecSol(0), stressSol(0)
|
: vecSol(0), vecSecSol(0), stressSol(0)
|
||||||
{
|
{
|
||||||
std::string variables;
|
std::string variables;
|
||||||
|
|
||||||
@@ -83,17 +92,18 @@ AnaSol::AnaSol (const TiXmlElement* elem, bool scalarSol)
|
|||||||
std::string primary = prim->FirstChild()->Value();
|
std::string primary = prim->FirstChild()->Value();
|
||||||
IFEM::cout <<"\tPrimary="<< primary << std::endl;
|
IFEM::cout <<"\tPrimary="<< primary << std::endl;
|
||||||
if (scalarSol)
|
if (scalarSol)
|
||||||
scalSol = new EvalFunction((variables+primary).c_str());
|
scalSol.push_back(new EvalFunction((variables+primary).c_str()));
|
||||||
else
|
else
|
||||||
vecSol = new VecFuncExpr(primary,variables);
|
vecSol = new VecFuncExpr(primary,variables);
|
||||||
}
|
}
|
||||||
|
|
||||||
prim = elem->FirstChildElement("scalarprimary");
|
prim = elem->FirstChildElement("scalarprimary");
|
||||||
if (prim && prim->FirstChild())
|
while (prim && prim->FirstChild())
|
||||||
{
|
{
|
||||||
std::string primary = prim->FirstChild()->Value();
|
std::string primary = prim->FirstChild()->Value();
|
||||||
IFEM::cout <<"\tScalar Primary="<< primary << std::endl;
|
IFEM::cout <<"\tScalar Primary="<< primary << std::endl;
|
||||||
scalSol = new EvalFunction((variables+primary).c_str());
|
scalSol.push_back(new EvalFunction((variables+primary).c_str()));
|
||||||
|
prim = prim->NextSiblingElement("scalarprimary");
|
||||||
}
|
}
|
||||||
|
|
||||||
const TiXmlElement* sec = elem->FirstChildElement("secondary");
|
const TiXmlElement* sec = elem->FirstChildElement("secondary");
|
||||||
@@ -102,17 +112,18 @@ AnaSol::AnaSol (const TiXmlElement* elem, bool scalarSol)
|
|||||||
std::string secondary = sec->FirstChild()->Value();
|
std::string secondary = sec->FirstChild()->Value();
|
||||||
IFEM::cout <<"\tSecondary="<< secondary << std::endl;
|
IFEM::cout <<"\tSecondary="<< secondary << std::endl;
|
||||||
if (scalarSol)
|
if (scalarSol)
|
||||||
scalSecSol = new VecFuncExpr(secondary,variables);
|
scalSecSol.push_back(new VecFuncExpr(secondary,variables));
|
||||||
else
|
else
|
||||||
vecSecSol = new TensorFuncExpr(secondary,variables);
|
vecSecSol = new TensorFuncExpr(secondary,variables);
|
||||||
}
|
}
|
||||||
|
|
||||||
sec = elem->FirstChildElement("scalarsecondary");
|
sec = elem->FirstChildElement("scalarsecondary");
|
||||||
if (sec && sec->FirstChild())
|
while (sec && sec->FirstChild())
|
||||||
{
|
{
|
||||||
std::string secondary = sec->FirstChild()->Value();
|
std::string secondary = sec->FirstChild()->Value();
|
||||||
IFEM::cout <<"\tScalar Secondary="<< secondary << std::endl;
|
IFEM::cout <<"\tScalar Secondary="<< secondary << std::endl;
|
||||||
scalSecSol = new VecFuncExpr(secondary,variables);
|
scalSecSol.push_back(new VecFuncExpr(secondary,variables));
|
||||||
|
sec = sec->NextSiblingElement("scalarsecondary");
|
||||||
}
|
}
|
||||||
|
|
||||||
const TiXmlElement* stress = elem->FirstChildElement("stress");
|
const TiXmlElement* stress = elem->FirstChildElement("stress");
|
||||||
|
|||||||
@@ -16,6 +16,7 @@
|
|||||||
|
|
||||||
#include "Function.h"
|
#include "Function.h"
|
||||||
#include <iostream>
|
#include <iostream>
|
||||||
|
#include <vector>
|
||||||
|
|
||||||
class TiXmlElement;
|
class TiXmlElement;
|
||||||
|
|
||||||
@@ -27,8 +28,8 @@ class TiXmlElement;
|
|||||||
class AnaSol
|
class AnaSol
|
||||||
{
|
{
|
||||||
protected:
|
protected:
|
||||||
RealFunc* scalSol; //!< Primary scalar solution field
|
std::vector<RealFunc*> scalSol; //!< Primary scalar solution fields
|
||||||
VecFunc* scalSecSol; //!< Secondary scalar solution field (gradient field)
|
std::vector<VecFunc*> scalSecSol; //!< Secondary scalar solution fields (gradient fields)
|
||||||
|
|
||||||
VecFunc* vecSol; //!< Primary vector solution field
|
VecFunc* vecSol; //!< Primary vector solution field
|
||||||
TensorFunc* vecSecSol; //!< Secondary solution field (vector gradient field)
|
TensorFunc* vecSecSol; //!< Secondary solution field (vector gradient field)
|
||||||
@@ -45,25 +46,24 @@ public:
|
|||||||
//! \details It is assumed that all the arguments are pointers to dynamically
|
//! \details It is assumed that all the arguments are pointers to dynamically
|
||||||
//! allocated objects, as the class destructor will attempt to delete them.
|
//! allocated objects, as the class destructor will attempt to delete them.
|
||||||
AnaSol(RealFunc* s1 = 0, VecFunc* s2 = 0,
|
AnaSol(RealFunc* s1 = 0, VecFunc* s2 = 0,
|
||||||
VecFunc* v1 = 0, TensorFunc* v2 = 0, STensorFunc* v3 = 0)
|
VecFunc* v1 = 0, TensorFunc* v2 = 0, STensorFunc* v3 = 0);
|
||||||
: scalSol(s1), scalSecSol(s2), vecSol(v1), vecSecSol(v2), stressSol(v3) {}
|
|
||||||
|
|
||||||
//! \brief Constructor initializing the primary and secondary solution fields.
|
//! \brief Constructor initializing the primary and secondary solution fields.
|
||||||
//! \param[in] s Primary scalar solution field
|
//! \param[in] s Primary scalar solution field
|
||||||
//! \param[in] sigma Symmetric stress tensor field
|
//! \param[in] sigma Symmetric stress tensor field
|
||||||
AnaSol(RealFunc* s, STensorFunc* sigma)
|
AnaSol(RealFunc* s, STensorFunc* sigma)
|
||||||
: scalSol(s), scalSecSol(0), vecSol(0), vecSecSol(0), stressSol(sigma) {}
|
: vecSol(0), vecSecSol(0), stressSol(sigma) { if (s) scalSol.push_back(s); }
|
||||||
|
|
||||||
//! \brief Constructor initializing the primary and secondary solution fields.
|
//! \brief Constructor initializing the primary and secondary solution fields.
|
||||||
//! \param[in] v Primary vector solution field
|
//! \param[in] v Primary vector solution field
|
||||||
//! \param[in] sigma Symmetric stress tensor field
|
//! \param[in] sigma Symmetric stress tensor field
|
||||||
AnaSol(VecFunc* v, STensorFunc* sigma)
|
AnaSol(VecFunc* v, STensorFunc* sigma)
|
||||||
: scalSol(0), scalSecSol(0), vecSol(v), vecSecSol(0), stressSol(sigma) {}
|
: vecSol(v), vecSecSol(0), stressSol(sigma) {}
|
||||||
|
|
||||||
//! \brief Constructor initializing the symmetric stress tensor field only.
|
//! \brief Constructor initializing the symmetric stress tensor field only.
|
||||||
//! \param[in] sigma Symmetric stress tensor field
|
//! \param[in] sigma Symmetric stress tensor field
|
||||||
AnaSol(STensorFunc* sigma)
|
AnaSol(STensorFunc* sigma)
|
||||||
: scalSol(0), scalSecSol(0), vecSol(0), vecSecSol(0), stressSol(sigma) {}
|
: vecSol(0), vecSecSol(0), stressSol(sigma) {}
|
||||||
|
|
||||||
//! \brief Constructor initializing expression functions by parsing a stream.
|
//! \brief Constructor initializing expression functions by parsing a stream.
|
||||||
//! \param is The file stream to read function definition from
|
//! \param is The file stream to read function definition from
|
||||||
@@ -78,11 +78,13 @@ public:
|
|||||||
//! \brief The destructor frees the analytical solution fields.
|
//! \brief The destructor frees the analytical solution fields.
|
||||||
virtual ~AnaSol()
|
virtual ~AnaSol()
|
||||||
{
|
{
|
||||||
if (scalSol) delete scalSol;
|
for (auto& it : scalSol)
|
||||||
if (scalSecSol) delete scalSecSol;
|
delete it;
|
||||||
if (vecSol) delete vecSol;
|
for (auto& it : scalSecSol)
|
||||||
if (vecSecSol) delete vecSecSol;
|
delete it;
|
||||||
if (stressSol) delete stressSol;
|
delete vecSol;
|
||||||
|
delete vecSecSol;
|
||||||
|
delete stressSol;
|
||||||
}
|
}
|
||||||
|
|
||||||
//! \brief Checks whether a scalar solution is defined.
|
//! \brief Checks whether a scalar solution is defined.
|
||||||
@@ -90,9 +92,9 @@ public:
|
|||||||
{
|
{
|
||||||
if (stressSol && !vecSecSol && !vecSol)
|
if (stressSol && !vecSecSol && !vecSol)
|
||||||
return 3;
|
return 3;
|
||||||
else if (scalSecSol)
|
else if (!scalSecSol.empty())
|
||||||
return 2;
|
return 2;
|
||||||
else if (scalSol)
|
else if (!scalSol.empty())
|
||||||
return 1;
|
return 1;
|
||||||
else
|
else
|
||||||
return 0;
|
return 0;
|
||||||
@@ -112,9 +114,11 @@ public:
|
|||||||
}
|
}
|
||||||
|
|
||||||
//! \brief Returns the scalar solution, if any.
|
//! \brief Returns the scalar solution, if any.
|
||||||
RealFunc*getScalarSol() const { return scalSol; }
|
RealFunc* getScalarSol(size_t idx = 0) const
|
||||||
|
{ return scalSol.size() <= idx ? nullptr : scalSol[idx]; }
|
||||||
//! \brief Returns the secondary scalar solution, if any.
|
//! \brief Returns the secondary scalar solution, if any.
|
||||||
VecFunc* getScalarSecSol() const { return scalSecSol; }
|
VecFunc* getScalarSecSol(size_t idx = 0) const
|
||||||
|
{ return scalSecSol.size() <= idx ? nullptr : scalSecSol[idx]; }
|
||||||
|
|
||||||
//! \brief Returns the vector solution, if any.
|
//! \brief Returns the vector solution, if any.
|
||||||
VecFunc* getVectorSol() const { return vecSol; }
|
VecFunc* getVectorSol() const { return vecSol; }
|
||||||
|
|||||||
Reference in New Issue
Block a user