Skip to content
Snippets Groups Projects
Commit 6146be39 authored by Éric Béchet's avatar Éric Béchet
Browse files

weight watchers

parent fd1b5772
Branches
Tags
No related merge requests found
......@@ -132,6 +132,28 @@ class dofManager{
{
numberDof(v->getNum(), Dof::createTypeWithTwoInts(iComp, iField));
}
inline void getDofValue(std::vector<Dof> &keys,std::vector<dataVec> &Vals)
{
int ndofs=keys.size();
Vals.reserve(Vals.size()+ndofs);
for (int i=0;i<ndofs;++i) Vals.push_back(getDofValue(keys[i]));
}
inline dataVec getDofValue(Dof key) const
{
{
typename std::map<Dof, dataVec>::const_iterator it = fixed.find(key);
if (it != fixed.end()) return it->second;
}
{
std::map<Dof, int>::const_iterator it = unknown.find(key);
if (it != unknown.end())
return _current->getFromSolution(it->second);
}
return dataVec(0.0);
}
inline dataVec getDofValue(int ent, int type) const
{
Dof key(ent, type);
......
This diff is collapsed.
......@@ -11,6 +11,7 @@
#include "SVector3.h"
#include "dofManager.h"
#include "simpleFunction.h"
#include "functionSpace.h"
class GModel;
class PView;
......@@ -24,68 +25,56 @@ struct elasticField {
elasticField () : g(0), _enrichment(0),_tag(0){}
};
struct dirichletBC {
struct BoundaryCondition
{
enum location{UNDEF,ON_VERTEX,ON_EDGE,ON_FACE,ON_VOLUME};
location onWhat; // on vertices or elements
int _tag; // tag for the dofManager
int _comp; // component
groupOfElements *g; // support for this BC
BoundaryCondition() : g(0),_tag(0),onWhat(UNDEF) {};
};
struct dirichletBC : public BoundaryCondition
{
int _comp; // component
simpleFunction<double> _f;
dirichletBC () : g(0),_comp(0),_tag(0){}
dirichletBC () :BoundaryCondition(),_comp(0),_f(0){}
};
struct neumannBC {
int _tag; // tag for the dofManager
groupOfElements *g; // support for this BC
struct neumannBC : public BoundaryCondition
{
simpleFunction<SVector3> _f;
neumannBC () : g(0),_tag(0),_f(SVector3(0,0,0)){}
neumannBC () : BoundaryCondition(),_f(SVector3(0,0,0)){}
};
// an elastic solver ...
class elasticitySolver{
class elasticitySolver
{
protected:
GModel *pModel;
int _dim, _tag;
dofManager<double> *pAssembler;
FunctionSpace<SVector3> *LagSpace;
// young modulus and poisson coefficient per physical
std::vector<elasticField> elasticFields;
// imposed nodal forces
std::map<int, SVector3> nodalForces;
// imposed line forces
std::map<int, SVector3> lineForces;
std::vector<neumannBC> edgeNeu;
// imposed face forces
std::map<int, SVector3> faceForces;
std::vector<neumannBC> faceNeu;
// imposed volume forces
std::map<int, SVector3> volumeForces;
std::vector<neumannBC> volumeNeu;
// imposed nodal displacements
std::map<std::pair<int,int>, double> nodalDisplacements;
// imposed edge displacements
std::map<std::pair<int,int>, double> edgeDisplacements;
std::vector<dirichletBC> edgeDiri;
// imposed face displacements
std::map<std::pair<int,int>, double> faceDisplacements;
std::vector<dirichletBC> faceDiri;
// std::vector<dirichletBC> faceDiri;
// neumann BC
std::vector<neumannBC> allNeumann;
// dirichlet BC
std::vector<dirichletBC> allDirichlet;
public:
elasticitySolver(int tag) : _tag(tag) {}
void addNodalForces (int iNode, const SVector3 &f)
elasticitySolver(int tag) : _tag(tag),LagSpace(0),pAssembler(0) {}
virtual ~elasticitySolver()
{
nodalForces[iNode] = f;
if (LagSpace) delete LagSpace;
if (pAssembler) delete pAssembler;
}
void addNodalDisplacement(int iNode, int dir, double val)
{
nodalDisplacements[std::make_pair(iNode, dir)] = val;
}
// void addElasticConstants(double e, double nu, int physical)
// {
// elasticConstants[physical] = std::make_pair(e, nu);
// }
void readInputFile(const std::string &meshFileName);
void setMesh(const std::string &meshFileName);
virtual void solve();
void readInputFile(const std::string &meshFileName);
virtual PView *buildDisplacementView(const std::string &postFileName);
// PView *buildVonMisesView(const std::string &postFileName);
virtual PView *buildVonMisesView(const std::string &postFileName);
// std::pair<PView *, PView*> buildErrorEstimateView
// (const std::string &errorFileName, double, int);
// std::pair<PView *, PView*> buildErrorEstimateView
......@@ -93,15 +82,4 @@ class elasticitySolver{
};
// another elastic solver ...
class MyelasticitySolver : public elasticitySolver
{
public:
MyelasticitySolver(int tag) : elasticitySolver(tag) {}
virtual void solve();
};
#endif
......@@ -68,12 +68,14 @@ class FunctionSpaceBase
public:
virtual int getNumKeys(MElement *ele)=0; // if one needs the number of dofs
virtual int getKeys(MElement *ele, std::vector<Dof> &keys)=0;
virtual int getNumKeys(MVertex *ver)=0; // if one needs the number of dofs
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys)=0;
};
template<class T>
class FunctionSpace : public FunctionSpaceBase
{
protected:
public:
typedef typename TensorialTraits<T>::ValType ValType;
typedef typename TensorialTraits<T>::GradType GradType;
/* typedef typename TensorialTraits<T>::HessType HessType;
......@@ -81,6 +83,7 @@ class FunctionSpace : public FunctionSpaceBase
typedef typename TensorialTraits<T>::CurlType CurlType;*/
public:
virtual int f(MVertex *ver, std::vector<ValType> &vals) {} // used for neumann BC
virtual int f(MElement *ele, double u, double v, double w, std::vector<ValType> &vals)=0;
virtual int gradf(MElement *ele, double u, double v, double w,std::vector<GradType> &grads)=0;
// virtual groupOfElements* getSupport()=0;// probablement inutile
......@@ -90,25 +93,23 @@ class FunctionSpace : public FunctionSpaceBase
virtual int curlf(MElement *ele, double u, double v, double w,std::vector<CurlType> &curls);*/
virtual int getNumKeys(MElement *ele)=0; // if one needs the number of dofs
virtual int getKeys(MElement *ele, std::vector<Dof> &keys)=0;
virtual int getNumKeys(MVertex *ver)=0; // if one needs the number of dofs
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys)=0;
};
class ScalarLagrangeFunctionSpace : public FunctionSpace<double>
{
protected:
public:
typedef TensorialTraits<double>::ValType ValType;
typedef TensorialTraits<double>::GradType GradType;
typedef TensorialTraits<double>::HessType HessType;
typedef TensorialTraits<double>::DivType DivType;
typedef TensorialTraits<double>::CurlType CurlType;
protected:
std::vector<basisFunction*> basisFunctions;
int _iField; // field number (used to build dof keys)
Dof getLocalDof(MElement *ele, int i) const
{
return Dof(ele->getVertex(i)->getNum(), _iField);
}
public:
ScalarLagrangeFunctionSpace(int i=0):_iField(i) {}
virtual int getId(void) const {return _iField;};
......@@ -119,6 +120,7 @@ class ScalarLagrangeFunctionSpace : public FunctionSpace<double>
vals.resize(curpos+ndofs);
ele->getShapeFunctions(u, v, w, &(vals[curpos]));
};
virtual int f(MVertex *ver, std::vector<ValType> &vals) {vals.push_back(1.0);} // used for neumann BC
virtual int gradf(MElement *ele, double u, double v, double w,std::vector<GradType> &grads)
{
int ndofs= ele->getNumVertices();
......@@ -143,8 +145,13 @@ class ScalarLagrangeFunctionSpace : public FunctionSpace<double>
int ndofs= ele->getNumVertices();
keys.reserve(keys.size()+ndofs);
for (int i=0;i<ndofs;++i)
keys.push_back(getLocalDof(ele,i));
getKeys(ele->getVertex(i),keys);
}
virtual int getNumKeys(MVertex *ver) { return 1;}
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys)
{
keys.push_back(Dof(ver->getNum(), _iField));
};
};
......@@ -180,6 +187,20 @@ public :
}
virtual ~ScalarToAnyFunctionSpace() {delete ScalarFS;}
virtual int f(MVertex *ver, std::vector<ValType> &vals)
{
std::vector<double> valsd;
ScalarFS->f(ver,valsd);
int nbdofs=valsd.size();
int nbcomp=comp.size();
int curpos=vals.size();
vals.reserve(curpos+nbcomp*nbdofs);
for (int j=0;j<nbcomp;++j)
{
for (int i=0;i<nbdofs;++i) vals.push_back(multipliers[j]*valsd[i]);
}
} // used for neumann BC
virtual int f(MElement *ele, double u, double v, double w, std::vector<ValType> &vals)
{
std::vector<double> valsd;
......@@ -236,14 +257,36 @@ public :
}
}
}
virtual int getNumKeys(MVertex *ver) { return ScalarFS->getNumKeys(ver)*comp.size();}
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys)
{
int nk=ScalarFS->getNumKeys(ver);
std::vector<Dof> bufk;
bufk.reserve(nk);
ScalarFS->getKeys(ver,bufk);
int nbdofs=bufk.size();
int nbcomp=comp.size();
int curpos=keys.size();
keys.reserve(curpos+nbcomp*nbdofs);
for (int j=0;j<nbcomp;++j)
{
for (int i=0;i<nk;++i)
{
int i1,i2;
Dof::getTwoIntsFromType(bufk[i].getType(), i1,i2);
keys.push_back(Dof(bufk[i].getEntity(),Dof::createTypeWithTwoInts(comp[j],i1)));
}
}
};
};
class VectorLagrangeFunctionSpace : public ScalarToAnyFunctionSpace<SVector3>
{
protected:
static const SVector3 BasisVectors[3];
public:
enum Along { VECTOR_X=0, VECTOR_Y=1, VECTOR_Z=2 };
static const SVector3 BasisVectors[3];
typedef TensorialTraits<SVector3>::ValType ValType;
typedef TensorialTraits<SVector3>::GradType GradType;
typedef TensorialTraits<SVector3>::HessType HessType;
......@@ -305,6 +348,12 @@ class CompositeFunctionSpace : public FunctionSpace<T>
delete (*it);
}
virtual int f(MVertex *ver, std::vector<ValType> &vals)
{
for (iterFS it=_spaces.begin(); it!=_spaces.end();++it)
(*it)->f(ver,vals);
}
virtual int f(MElement *ele, double u, double v, double w,std::vector<ValType> &vals)
{
for (iterFS it=_spaces.begin(); it!=_spaces.end();++it)
......@@ -330,11 +379,26 @@ class CompositeFunctionSpace : public FunctionSpace<T>
for (iterFS it=_spaces.begin(); it!=_spaces.end();++it)
(*it)->getKeys(ele,keys);
}
virtual int getNumKeys(MVertex *ver)
{
int ndofs=0;
for (iterFS it=_spaces.begin(); it!=_spaces.end();++it)
ndofs+=(*it)->getNumKeys(ver);
return ndofs;
}
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys)
{
for (iterFS it=_spaces.begin(); it!=_spaces.end();++it)
(*it)->getKeys(ver,keys);
};
};
template<class T>
class xFemFunctionSpace : public FunctionSpace<T>
{
public:
typedef typename TensorialTraits<T>::ValType ValType;
typedef typename TensorialTraits<T>::GradType GradType;
typedef typename TensorialTraits<T>::HessType HessType;
......@@ -344,10 +408,14 @@ class xFemFunctionSpace : public FunctionSpace<T>
FunctionSpace<T>* _spacebase;
// Function<double>* enrichment;
public:
virtual int f(MVertex *ver, std::vector<ValType> &vals);
virtual int f(MElement *ele, double u, double v, double w,std::vector<ValType> &vals);
virtual int gradf(MElement *ele, double u, double v, double w,std::vector<GradType> &grads);
int getNumKeys(MElement *ele);
void getKeys(MElement *ele, std::vector<Dof> &keys);
virtual int getNumKeys(MElement *ele);
virtual void getKeys(MElement *ele, std::vector<Dof> &keys);
virtual int getNumKeys(MVertex *ver);
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys);
};
......@@ -364,10 +432,13 @@ class FilteredFunctionSpace : public FunctionSpace<T>
F &_filter;
public:
virtual int f(MVertex *ver, std::vector<ValType> &vals);
virtual int f(MElement *ele, double u, double v, double w,std::vector<ValType> &vals);
virtual int gradf(MElement *ele, double u, double v, double w,std::vector<GradType> &grads);
int getNumKeys(MElement *ele);
void getKeys(MElement *ele, std::vector<Dof> &keys);
virtual int getNumKeys(MElement *ele);
virtual void getKeys(MElement *ele, std::vector<Dof> &keys);
virtual int getNumKeys(MVertex *ver);
virtual int getKeys(MVertex *ver, std::vector<Dof> &keys);
};
......
......@@ -21,6 +21,7 @@
#include "MVertex.h"
template<class Iterator,class Assembler> void Assemble(BilinearTermBase &term,FunctionSpaceBase &space,Iterator itbegin,Iterator itend,QuadratureBase &integrator,Assembler &assembler) // symmetric
{
fullMatrix<typename Assembler::dataMat> localMatrix;
......@@ -64,6 +65,20 @@ template<class Iterator,class Assembler> void Assemble(LinearTermBase &term,Func
}
}
template<class Iterator,class Assembler> void Assemble(LinearTermBase &term,FunctionSpaceBase &space,Iterator itbegin,Iterator itend,Assembler &assembler)
{
fullVector<typename Assembler::dataMat> localVector;
std::vector<Dof> R;
for (Iterator it = itbegin;it!=itend; ++it)
{
MVertex *v = *it;
R.clear();
term.get(v,localVector);
space.getKeys(v,R);
assembler.assemble(R, localVector);
}
}
template<class Assembler> void Assemble(LinearTermBase &term,FunctionSpaceBase &space,MElement *e,QuadratureBase &integrator,Assembler &assembler)
{
fullVector<typename Assembler::dataMat> localVector;
......@@ -75,7 +90,6 @@ template<class Assembler> void Assemble(LinearTermBase &term,FunctionSpaceBase &
assembler.assemble(R, localVector);
}
template<class Assembler> void FixDofs(Assembler &assembler,std::vector<Dof> &dofs,std::vector<typename Assembler::dataVec> &vals)
{
int nbff=dofs.size();
......@@ -86,6 +100,12 @@ template<class Assembler> void FixDofs(Assembler &assembler,std::vector<Dof> &do
}
class FilterDof
{
public:
virtual bool operator()(Dof key)=0;
};
class FilterDofTrivial
{
public:
virtual bool operator()(Dof key) {return true;}
......@@ -100,18 +120,15 @@ class FilterDofComponent :public FilterDof
{
int type=key.getType();
int icomp,iphys;
Dof::getTwoIntsFromType(type, iphys, icomp);
Dof::getTwoIntsFromType(type, icomp, iphys);
if (icomp==comp) return true;
return false;
}
};
template<class Iterator,class Assembler> void FixNodalDofs(FunctionSpaceBase &space,Iterator itbegin,Iterator itend,Assembler &assembler,simpleFunction<typename Assembler::dataVec> &fct,FilterDof filter)
{
for (Iterator it=itbegin;it!=itend;++it)
template<class Assembler> void FixNodalDofs(FunctionSpaceBase &space,MElement *e,Assembler &assembler,simpleFunction<typename Assembler::dataVec> &fct,FilterDof &filter)
{
MElement *e=*it;
std::vector<MVertex*> tabV;
int nv=e->getNumVertices();
std::vector<Dof> R;
......@@ -135,8 +152,28 @@ template<class Iterator,class Assembler> void FixNodalDofs(FunctionSpaceBase &sp
}
}
}
template<class Assembler> void FixNodalDofs(FunctionSpaceBase &space,MVertex *v,Assembler &assembler,simpleFunction<typename Assembler::dataVec> &fct,FilterDof &filter)
{
std::vector<Dof> R;
space.getKeys(v,R);
for (std::vector<Dof>::iterator itd=R.begin();itd!=R.end();++itd)
{
Dof key=*itd;
if (filter(key))
{
assembler.fixDof(key, fct(v->x(),v->y(),v->z()));
}
}
}
template<class Iterator,class Assembler> void FixNodalDofs(FunctionSpaceBase &space,Iterator itbegin,Iterator itend,Assembler &assembler,simpleFunction<typename Assembler::dataVec> &fct,FilterDof &filter)
{
for (Iterator it=itbegin;it!=itend;++it)
{
FixNodalDofs(space,*it,assembler,fct,filter);
}
}
template<class Iterator,class Assembler> void NumberDofs(FunctionSpaceBase &space,Iterator itbegin,Iterator itend,Assembler &assembler)
{
......
......@@ -22,6 +22,7 @@
#include "groupOfElements.h"
// evaluation of a field ???
template<class T>
class Field {
......@@ -80,6 +81,7 @@ class LinearTermBase
{
public:
virtual void get(MElement *ele,int npts,IntPt *GP,fullVector<double> &v) =0;
virtual void get(MVertex *v,fullVector<double> &m) =0;
};
template<class S1> class LinearTerm : public LinearTermBase
......@@ -94,7 +96,7 @@ template<class S1> class LinearTerm : public LinearTermBase
class ScalarTermBase
{
public :
virtual void get(MElement *ele,int npts,IntPt *GP,double &v) =0;
virtual void get(MElement *ele,int npts,IntPt *GP,double &val) =0;
};
class ScalarTerm : public ScalarTermBase
......@@ -105,6 +107,60 @@ class ScalarTerm : public ScalarTermBase
template<class S1,class S2> class LaplaceTerm : public BilinearTerm<S1,S2>
{
public :
LaplaceTerm(S1& space1_,S2& space2_) : BilinearTerm<S1,S2>(space1_,space2_)
{}
virtual void get(MElement *ele,int npts,IntPt *GP,fullMatrix<double> &m)
{
Msg::Error("LaplaceTerm<S1,S2> w/ S1!=S2 not implemented");
}
virtual void get(MVertex *v,fullMatrix<double> &m)
{
Msg::Error("LaplaceTerm<S1,S2> w/ S1!=S2 not implemented");
}
}; // class
template<class S1> class LaplaceTerm<S1,S1> : public BilinearTerm<S1,S1> // symmetric
{
public :
LaplaceTerm(S1& space1_) : BilinearTerm<S1,S1>(space1_,space1_)
{}
virtual void get(MElement *ele,int npts,IntPt *GP,fullMatrix<double> &m)
{
int nbFF = BilinearTerm<S1,S1>::space1.getNumKeys(ele);
double jac[3][3];
m.resize(nbFF, nbFF);
m.setAll(0.);
for (int i = 0; i < npts; i++)
{
const double u = GP[i].pt[0]; const double v = GP[i].pt[1]; const double w = GP[i].pt[2];
const double weight = GP[i].weight; const double detJ = ele->getJacobian(u, v, w, jac);
std::vector<typename S1::GradType> Grads;
BilinearTerm<S1,S1>::space1.gradf(ele,u, v, w, Grads);
for (int j = 0; j < nbFF; j++)
{
for (int k = j; k < nbFF; k++)
{
double contrib=weight * detJ * dot(Grads[j],Grads[k]);
m(j,k)+=contrib;
if (j!=k) m(k,j)+=contrib;
}
}
}
// m.print("");
// exit(0);
}
}; // class
template<class S1,class S2> class IsotropicElasticTerm : public BilinearTerm<S1,S2> // non symmetric
{
protected :
......@@ -238,9 +294,9 @@ inline double dot(const double &a, const double &b)
template<class S1> class LoadTerm : public LinearTerm<S1>
{
simpleFunction<typename S1::ValType> Load;
simpleFunction<typename S1::ValType> &Load;
public :
LoadTerm(S1& space1_,simpleFunction<typename S1::ValType> Load_) :LinearTerm<S1>(space1_),Load(Load_) {};
LoadTerm(S1& space1_,simpleFunction<typename S1::ValType> &Load_) :LinearTerm<S1>(space1_),Load(Load_) {};
virtual void get(MElement *ele,int npts,IntPt *GP,fullVector<double> &m)
{
double nbFF=LinearTerm<S1>::space1.getNumKeys(ele);
......@@ -253,12 +309,29 @@ template<class S1> class LoadTerm : public LinearTerm<S1>
const double weight = GP[i].weight;const double detJ = ele->getJacobian(u, v, w, jac);
std::vector<typename S1::ValType> Vals;
LinearTerm<S1>::space1.f(ele,u, v, w, Vals);
SPoint3 p;
ele->pnt(u, v, w, p);
typename S1::ValType load=Load(p.x(),p.y(),p.z());
for (int j = 0; j < nbFF ; ++j)
{
m(j)+=dot(Vals[j],Load(u,v,w))*weight*detJ;
m(j)+=dot(Vals[j],load)*weight*detJ;
}
}
}
virtual void get(MVertex *ver,fullVector<double> &m)
{
double nbFF=LinearTerm<S1>::space1.getNumKeys(ver);
double jac[3][3];
m.resize(nbFF);
std::vector<typename S1::ValType> Vals;
LinearTerm<S1>::space1.f(ver, Vals);
typename S1::ValType load=Load(ver->x(),ver->y(),ver->z());
for (int j = 0; j < nbFF ; ++j)
{
m(j)=dot(Vals[j],load);
}
}
};
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment