Class xlifepp::SuTermVector#
-
class SuTermVector : public xlifepp::Term#
-
Inheritence diagram for xlifepp::SuTermVector:
Collaboration diagram for xlifepp::SuTermVector:
handles numerical representation of single unknown linear forms (SuLinearForm).
Internal class
Public Functions
-
SuTermVector(const LcTerm<SuTermVector>&)#
-
constructor from a linear combination of SuTermVector’s
-
template<typename T>
SuTermVector(const string_t&, const Unknown*, Space*, const T&)#
-
basic constructor from a constant value
-
SuTermVector(const string_t&, const Unknown*, Space*, ValueType vt = _real, number_t n = 0, dimen_t nv = 0, bool noass = false)#
-
basic constructor from explicit data type
-
SuTermVector(const SuTermVector&)#
-
copy constructor (full copy)
-
SuTermVector(const SuTermVector&, const SuTermVector&, const SuTermVector&, const SymbolicFunction &fs, const string_t &na = "")#
-
symb_fun(SuTermVector, SuTermVector, SuTermVector)
build the SuTermVector R = fs(X1,X2,X3) where fs is a symbolic function of x_1 x_2 x_3 : R(i)=fs(X1(i),X2(i),X3(i)) i=1,n X1, X2, X3 have to be three scalar SutermVectors defined on the same space the unknown of the SuTermVector result is the unknown of X1 it can be change later using changeUnknown function
-
SuTermVector(const SuTermVector&, const SuTermVector&, const SymbolicFunction &fs, const string_t &na = "")#
-
symb_fun(SuTermVector, SuTermVector)
build the SuTermVector R = fs(X1,X2) where fs is a symbolic function of x1_ x2_ when X1 and X2 are scalar SuTermVector, performs R(i)=fs(X1(i),X2(i)) i=1,n result is a scalar SuTermVector with X1 unknown when X1 is a vector SuTermVector (dim=d) and X2 a scalar SuTermvector, performs component by component R_k(i)=fs(X1_k(i),X2(i)) i=1,n, k=1,d, result is a vector SuTermVector with X1 unknown when X1 is a scalar SuTermVector and X2 a vector SuTermvector (dim=d), performs component by component R_k(i)=fs(X1(i),X2_k(i)) i=1,n, k=1,d, result is a vector SuTermVector with X2 unknown when X1 and X2 are both vector SuTermvector, computation is not available
NB: unknown can be change later using changeUnknown function
-
SuTermVector(const SuTermVector&, const SuTermVector&, funSC2_t &f, const string_t &na = "")#
-
complex_fun(SuTermVector, SuTermVector)
-
SuTermVector(const SuTermVector&, const SuTermVector&, funSR2_t &f, const string_t &na = "")#
-
real_fun(SuTermVector, SuTermVector)
apply function to two SutermVector’s build new SutermVector Z from SuTermVectors X and Y with Z(i)=f(X(i),Y(i))
-
SuTermVector(const SuTermVector&, const SymbolicFunction&, const string_t &na = "")#
-
symb_fun(SuTermVector)
build the SuTermVector R = fs(X) where fs is a symbolic function of x when X is a scalar SuTermVector, performs R(i)=fs(X(i)) i=1,n result is a scalar SuTermVector with X unknown when X is a vector SuTermVector (dim=d), performs component by component R_k(i)=fs(X_k(i)) i=1,n, k=1,d, result is a vector SuTermVector with X unknown
-
template<typename T>
SuTermVector(const SuTermVector&, const T&)#
-
copy constructor and assign a constant value
copy constructor from SuTermVector with values set to a constant value
-
SuTermVector(const SuTermVector&, funSC1_t &f, const string_t &na = "")#
-
complex_fun(SuTermVector)
-
SuTermVector(const SuTermVector&, funSR1_t &f, const string_t &na = "")#
-
real_fun(SuTermVector)
apply function to SutermVector build new SutermVector Y from SuTermVector X with Y(i)=f(X(i)) function f has to be consistent with SuTermVector coefficients use a scalar function for scalar SuTermVector and vector function for vector SuTermVector
-
SuTermVector(const SuTermVector &sutv, const GeomDomain &dom, const string_t &na = "")#
-
SuTermVector from a SuTermVector with geometrical mapping.
SuTermVector from the restriction of a SuTermVector to a given domain dom raise an error if dom does not intersect the SutermVector domain SutermVector sutv has to be computed and its entries_p has to be non zero In case of entries_p==nullptr and scalar_entries_p!=nullptr goto vector representation before calling this.
SuTermVector from the restriction of a SuTermVector to a given domain
-
SuTermVector(const Unknown&, const std::list<const SuTermVector*>&, const string_t& = "")#
-
concatenate n scalar SuTermVector’s in a vector SutermVector
construct vector SutermVector from n scalar SuTermVector’s scalar SuTermVector’s have to be defined on the same space vu has to be a vector unknown at least n
-
SuTermVector(const Unknown&, const SuTermVector&, const SuTermVector&, const string_t& = "")#
-
concatenate 2 scalar SuTermVector’s in a vector SutermVector
construct vector SutermVector from 2 scalar SuTermVector’s scalar SuTermVector’s have to be defined on the same space vu has to be a vector unknown at least d=2
-
SuTermVector(const Unknown &u, const GeomDomain &dom, const Function &f, const Function &gf, const Function g2f, const string_t &na)#
-
interpolation value constructor for non Lagrange interpolation involving second derivative dofs
-
SuTermVector(const Unknown &u, const GeomDomain &dom, const Function &f, const Function &gf, const string_t &na)#
-
interpolation constructor for non Lagrange interpolation involving first derivative dofs
-
SuTermVector(const Unknown &u, const GeomDomain &dom, const SuTermVector &sut, const Function *fmap, bool useNearest, bool errorOnOutDom, real_t tol, const string_t &na)#
-
construct SuTermVector from a SuTermVector with geometrical mapping, the unknowns may differ u : unknown attached to new SuTermVector (say nsut) dom : domain related to new SuTermVector sut : original SuTermVector (involving unknown v and domain omega) fmap : optional mapping function from dom to omega (=0 ,identity is assumed) errorOnOutDom : control behaviour on locatig error tol : tolerance factor used in error management
for each node Pi of dom, find the nearest point Qi belonging to an element Ek of omega (projection) if dist(Qi,omega) < theTolerance (i.e Qi=Pi) : nsut[i] = sut(Pi) (v-interpolation) if dist(Qi,omega) > theTolerance compute the relative distance ei=dist(Qi,omega)/E.size() if (ei > tol) stop computation if errorOnOutDom = true else nsut[i] = 0 if (ei < tol) nsut[i] = sut(Qi) (v-interpolation), throw a warning if errorOnOutDom = true
this process involves mesh interpolation and is restricted to Lagrange interpolation spaces it is assumed that the two unknowns have same structure (scalar/vector)
-
template<typename T>
SuTermVector(const Unknown &u, const GeomDomain &dom, const T &v, const string_t &na, bool nodalFct, bool noass)#
-
constructor from constant scalar value SuTermVector from explicit real scalar function (interpolation on Lagrange Dof’s)
SuTermVector constructor from Unknown and real value, explicit vector construction from a constant value (no linear form) T may be a real/complex scalar or a real/complex vector (only in the case of a vector unknown)
-
template<typename T>
SuTermVector(const Unknown &u, const GeomDomain &dom, const Vector<T> &v, const string_t &na, bool nodalFct, bool noass)#
-
constructor from constant vector or a vector
SuTermVector constructor from Unknown and real value, explicit vector construction from a given vector (no linear form) T may be a real/complex scalar or a real/complex vector (only in the case of a vector unknown) to raise the ambiguity between u(i)=v i=1,n and u = v we test the number of components (nbc) of the unknown and the size (s) of v if nbc=s we apply u(i)=v i=1,n else we apply u = v Note: if v is too small compared to number of dofs, the omitted values are replaced by 0.
-
SuTermVector(const Unknown &u, const GeomDomain &dom, funSC_t &f, const string_t &na, bool nodalFct, bool noass)#
-
SuTermVector from explicit complex scalar function (interpolation on Lagrange Dof’s)
-
SuTermVector(const Unknown &u, const GeomDomain &dom, funVC_t &f, const string_t &na, bool nodalFct, bool noass)#
-
SuTermVector from explicit complex vector function (interpolation on Lagrange Dof’s)
-
SuTermVector(const Unknown &u, const GeomDomain &dom, funVR_t &f, const string_t &na, bool nodalFct, bool noass)#
-
SuTermVector from explicit real vector function (interpolation on Lagrange Dof’s)
-
SuTermVector(SuLinearForm *sulf = nullptr, const string_t &na = "", bool noass = false)#
-
basic constructor from SuLinearForm
-
virtual ~SuTermVector()#
-
destructor
-
VectorEntry *actual_entries() const#
-
return the actual pointer to entries (priority to scalar entry)
return the actual pointer to entries (priority to scalar entry), return 0 is no allocated pointer
-
void adjustScalarEntries(const std::vector<DofComponent>&)#
-
adjust scalar entries to cdofs numbering given
adjust scalar entries to cdofs numbering given same as extendScalarTo but different method, see general adjustScalarEntries function
-
SuTermVector &applyEssentialConditions(const Constraints&)#
-
apply constraints (essential conditions) to current SuTermVector
< apply essential condition to a SuTermVector (user tool)
-
template<typename T>
Vector<T> &asVector(Vector<T>&) const#
-
reinterpret SuTermVector as a raw Vector<T>
-
void buildSubspaces()#
-
build the required subspaces and the largest subspace
build the required subspaces and the largest subspace (possibly the whole space) this function:
travel along basic linear forms of the linear combination
identify the largest geometric domain
build the largest subspace (for SuTermVector dofs numbering)
create subspace related to each basic linear form if it does not exist
create dofs numbering between domain subspaces ands largest subspace
dofs numbering: whole space largest subspace domain subspace 1 5 2 2 3 1 3 1 4 5
-
inline std::vector<DofComponent> &cdofs()#
-
access to cdofs_r vector (r/w)
-
inline const std::vector<DofComponent> &cdofs() const#
-
access to cdofs_r vector (r)
-
void changeUnknown(const Unknown&, const Vector<number_t>&)#
-
move SuTermVector unknown to an other
move SuTermVector unknown to an other one consist in changing the unknown u of SuTermVector by the unknown v (u and v must be defined on the same Space!) with an optionnal management of components of unknown, examples u a scalar unknown and v a scalar unknown: SuTermVector if not modified u a scalar unknown and v a 3 vector unknown, ci ={2} the SuTermVector is changed to a 3 column vector: new sutv=[ [0…0] [sutv] [0…0] ] u a 3-vector unknown and v a 3-vector unknown, ci ={2 3 1} the SuTermVector is changed columns are permuted: new sutv=[ [sutv2] [sutv3] [sutv1] ] v may be equal to u in that case u a 4-vector unknown and v a 2-vector unknown, ci ={2 4} the SuTermVector is changed: new sutv=[ [sutv2] [sutv4] ] u a 4-vector unknown and v a 3-vector unknown, ci ={2 0 4} the SuTermVector is changed: new sutv=[ [sutv2] [0…0] [sutv4] ] u is replaced by v, v may be equal to u but be cautious if the result is not of the dimension of v because many operations are not supported in doubt, create a new unknown of the consistent dimension
-
virtual void clear()#
-
deallocate memory used by a SuTermVector
-
virtual void compute()#
-
compute SuTermVector
-
void compute(const LcTerm<SuTermVector>&)#
-
compute SuTermVector from a linear combination of SuTermVector’s
compute SuTermVector from a linear combination of SuTermVector’s the current SuTermVector is assumed to be void, this function fills every things if it is not the case, every thing is clear
-
template<typename T, typename K>
void computeFE(const SuLinearForm&, Vector<T>&, K&)#
-
compute SuTermVector (FE and FEext computation)
-
template<typename T, typename K>
void computeFE(const SuLinearForm &sulf, Vector<T> &v, K &vt)
-
FE computation of a scalar SuLinearForm on a == unique domain ==.
- Parameters:
-
sulf – a single unknown linear form defined on a unique domain (assumed)
v – reference to vector entries of type T
vt – to pass as template the scalar type (real or complex)
-
template<typename T, typename K>
void computeIR(const SuLinearForm &sulf, Vector<T> &res, K &vt, const std::vector<Point> &xs, const Vector<T> &U, const std::vector<Vector<real_t>> *nxs)#
-
FE computation of a scalar SuLinearForm with kernels on a == unique domain ==.
examples: v(x)=intg_gamma G(x,y)*U(y) dy v(x)=intg_gamma n_y.grad_y(G(x,y))*U(y) dy
U comes from a SuTermVector having same unknown as SuLinearForm !
Algorithm: various initialisation pre-computation of shape values on ref element, geometric data and quadrature points on each element loop on finite elements related to the integration domain | initialisation related to the element (geometric element, dof numbering, n_y) | loop on linear forms | if more than one lf, load some information related to linear form | if ny required, transmit n_y to kernel | loop on xi points (parallelized) | compute relative distance(xi,e)/diameter(e) | loop on integration methods related to the distance | compute intg_E opk(xi,y)*opu(y)dy by quadrature or specific method
Result of integral representation may be a scalar, a vector or a matrix, always stored as a vector of scalars!
Integration methods are described by the IntegrationMethods object handled by each linear form. Such object manages a list of (IntegrationMethod, FunctionPart, bound) that allows to deal with complex choices of integration methods. For instance, for point x on gamma (bound=0) : specific method for the singular part, standard quadrature method (e.g GaussLegendreRule, 5) for the regular part for point x close to gamma (bound=1) : high quadrature method (e.g GaussLegendreRule, 20) for all part of the function for point x not so far from gamma (bound=3) : quadrature method (e.g Gasus_Legendre, 10) for all part of the function for point x not far from gamma (bound=inf) : low quadrature method (e.g Gasus_Legendre, 3) for all part of the function
This function supports scalar and vector unknown, scalar and vector shape functions, and domain extension if required (non tangential operators or volumic FE with no FE trace)
- Parameters:
-
sulf – a single unknown linear form defined on a unique domain (assumed)
vt – to pass as template the scalar type (real or complex)
xs – points where to evaluate IR
nxs – pointer to vector of normals at points where to evaluate IR (null if not available)
U – vector to apply to IR; U has always a scalar representation even if it is related to a vector unknown
res – reference to vector entries of type T (result)
-
template<typename T, typename K>
void computeIROld(const SuLinearForm&, Vector<T> &res, K &vt, const std::vector<Point> &xs, const Vector<T>&, const std::vector<Vector<real_t>> *nxs = nullptr)#
-
compute scalar SuTermVector (IR computation)
-
template<typename T, typename K>
void computeSP(const SuLinearForm&, Vector<T>&, K&)#
-
compute SuTermVector (SP computation) compute integral representation using integration methods
-
template<typename T, typename K>
void computeSP(const SuLinearForm &sulf, Vector<T> &v, K &vt)
-
SP computation of a scalar SuLinearForm on a == unique domain ==.
- Parameters:
-
sulf – a single unknown linear form defined on a unique domain (assumed) and involving a spectral unknown
v – reference to vector entries of type T
vt – to pass as template the scalar type (real or complex)
-
void copy(const SuTermVector&)#
-
full copy of SuTermVector
-
number_t dofRank(const Dof &dof) const#
-
return rank in vector of dof
return rank in vector of a given dof
-
const GeomDomain *domain() const#
-
return geometric domain pointer attached to SuTermVector (may be 0)
-
inline VectorEntry *&entries()#
-
access to entries pointer (r/w)
-
inline const VectorEntry *entries() const#
-
access to entries pointer (r)
-
Value evaluate(const Point&, const Element*, bool useNearest = false, bool errorOnOutDom = true) const#
-
evaluate interpolated value at a point, may be knowing Element (general function)
-
inline Value evaluate(const Point &p, bool useNearest = false, bool errorOnOutDom = true) const#
-
evaluate interpolated value at a point
-
inline Value evaluate(const Point &p, const Element &elt) const#
-
< evaluate interpolated value at a point, knowing Element
-
void extendScalarTo(const std::vector<DofComponent>&, bool useDual = false)#
-
extend current scalar entries according to a DofComponent numbering vector produce the mapped SuTermVector from current GeomDomains to a GeomDomain
extend current vector to numbering given by DofComponent vector ndofs if useDual is true, use dual components dof (default = false) update scalar entries and NOT vector entries if scalar_entries = 0 move to scalar representation Note: differs from previous extendTo functions by the fact it works with component dof numbering instead of dof numbering
-
void extendTo(const Space&)#
-
extend current vector according to numbering of a space
extend current vector to numbering of space sp update entries and NOT scalar entries
-
void extendTo(const SuTermVector&)#
-
extend current vector according to numbering of an other vector
extend current vector to numbering of vector v (useful for algebraic computation on SuTermVector) update entries and NOT scalar entries if allocated
-
void initFromFunction(const Unknown &u, const GeomDomain &dom, const Function &f, const Function &gf, const Function &g2f, const string_t &na)#
-
initialize SuTermVector from functions (interpolation on non Lagrange Dof’s)
SuTermVector constructor for non Lagrange interpolation involving second derivative dofs (Argyris) when scalar unknown and scalar dof, f returns a scalar, gf the vector ( dxf, [dyf, [dzf]] ), g2f the vector ( dxxf, [dyyf, dxyf, [dzzf, dxz2f, dyzf]] ) when vector unknown or vector dof, f returns a vector, gf the vector ( dxf1, dxf2, dxf3, [dyf1, dyf2, dyf3, [dzf1,, dzf2, dzf3]] ), (not both!) g2f the vector ( dxxf1, [dyyf1, dxyf1, [dzzf1, dxz2f1, dyzf1]], dxxf2, [dyyf2, dxyf2, [dzzf2, dxz2f2, dyzf2]], …) create the vector d_i(f,gf,g2f) where d_i is the ith dof viewed as a distribution using Dof::operator(f,gf,g2f)
-
void initFromFunction(const Unknown &u, const GeomDomain &dom, const OperatorOnFunction &opf, const string_t &na, bool noass)#
-
initialize SuTermVector from function (interpolation on Lagrange Dof’s)
-
SuTermVector *interpolate(const Unknown&, const GeomDomain&)#
-
interpolate on u space and geomdomain locate element of the SuTermVector FESpace containing P
interpolate current sutermvector on u space and domain dom sutr(i) = sum_j sutv_j w_j(Mi) Mi nodes of Lagrange u-space, Mi in dom wj shape functions related to sutv return a new SuTermMatrix allocated here quite expansive !!!
-
Element *locateElementP(const Point &p, bool useNearest = false, bool errorOnOutDom = true) const#
-
locate element (pointer) of the SuTermVector FESpace containing P
-
SuTermVector *mapTo(const GeomDomain&, const Unknown&, bool useNearest = false, bool erDom = true, Function *fmap = nullptr, real_t tol = theLocateToleranceFactor) const#
-
produce the mapped SuTermVector from current GeomDomains to GeomDomain dom current SuTermVector has to be already computed
-
complex_t maxValAbs() const#
-
return the value (in complex) of component being the largest one in absolute value
return the value (in complex) of component being the largest one in absolute value travel all SCALAR components of SuTermVector using scalar entries if exists
check if 2 SuTermVector have same dofs
-
SuTermVector &merge(const SuTermVector&)#
-
merge a SuTermVector with preserving the current value when common dofs
-
inline virtual void name(const string_t &nam)#
-
update the name of a SuTermVector
-
inline number_t nbOfComponents() const#
-
number of components when each element is a vector (1 if scalar element)
-
SuTermVector *onDomain(const GeomDomain &dom) const#
-
restriction of SuTermVector to a domain, return 0 pointer if void
-
SuTermVector operator()(const ComponentOfUnknown&) const#
-
extract unknown component term vector as a SuTermVector
extract unknown component term vector as a SuTermVector, create a SuTermVector this function must be called from TermVector class in order to encapsulate SuTermVector in a TermVector
-
template<typename T>
T &operator()(const std::vector<real_t>&, T&) const#
-
compute interpolated value of unknown at point
interpolation functions compute interpolated value of unknown at point p, T has to be compliant with SuTermVector entries
-
template<typename T>
SuTermVector &operator*=(const T&)#
-
operation U*=real
operation U*=t
-
SuTermVector &operator+=(const SuTermVector&)#
-
operation U+=V
-
SuTermVector &operator-=(const SuTermVector&)#
-
operation U+=V
-
template<typename T>
SuTermVector &operator/=(const T&)#
-
operation U/=real
-
SuTermVector &operator=(const SuTermVector&)#
-
assign operator
-
virtual void print(std::ostream&) const#
-
print SuTermVector
-
void print(std::ostream&, bool r, bool h) const#
-
print SuTermVector with raw and header option
-
SuTermVector &round(real_t prec = 0.)#
-
round with a given precision (0 no round)
-
SuTermVector &roundToZero(real_t aszero = std::sqrt(theEpsilon))#
-
round to zero all coefficients close to 0 (|.
| < aszero)
-
inline void saveToFile(const string_t &fname, bool encode) const#
-
save SuTermVector to file
-
void saveToFile(const string_t &fname, number_t prec = fullPrec, bool encode = false) const#
-
save SuTermVector to file
-
inline VectorEntry *&scalar_entries()#
-
access to entries pointer (r/w)
-
inline const VectorEntry *scalar_entries() const#
-
access to entries pointer (r)
-
void setValue(const Function&, const GeomDomain&)#
-
set values of SuTermVector from Function, restricted to a domain: sut|dom=f
set dof values of GeomDomain from Function: u_i = f(M_i) for all dof i on dom, M_i point support of dof i does not change other values restricted to FE space dom has to be a subdomain of the domain of the currrent SuTermVector be cautious when Space is not a FE Lagrange space, because does not imply that u=f on dom in the meaning of function
-
void setValue(const SuTermVector&)#
-
set values of SuTermVector from SuTermVector
set dof values of SutermVector from a given SuTermVector: u_i = v_i for all dof i of the given SutermVector does not change other values
-
void setValue(const SuTermVector&, const GeomDomain&)#
-
set values of SuTermVector from SuTermVector, restricted to a domain: sut|dom= sutv|dom
set dof values of GeomDomain from SuTermVector: u_i = v_i for all dof i on dom does not change other values dom has to be a subdomain of the domain of the currrent SuTermVector be cautious when Space is not a FE Lagrange space, because does not imply that u=v on dom in the meaning of function
-
void setValue(const Value&, const GeomDomain&)#
-
set values of SuTermVector from Value, restricted to a domain: sut|dom=v
set dof values of SutermVector on a GeomDomainfrom from Value: u_i = v for all dof i on dom does not change other values dom has to be a subdomain of the domain of the currrent SuTermVector be cautious when Space is not a FE Lagrange space because does not imply that u=v on dom in the meaning of function
-
inline SuLinearForm *&sulfp()#
-
access to single unknown linear form (pointer) (writable)
-
inline SuLinearForm *sulfp() const#
-
access to single unknown linear form (pointer)
-
SuTermVector &toAbs()#
-
change to abs(U)
-
SuTermVector &toComplex()#
-
change to complex representation
-
SuTermVector &toConj()#
-
change to conj(U)
-
SuTermVector &toImag()#
-
change to imaginary part of U
-
SuTermVector &toReal()#
-
change to real part of U
-
void toScalar(bool keepVector = false)#
-
go to scalar representation
create scalar representation scalar_entries_p and cdofs_ numbering vectors
if unknown is a vector unknown a new VectorEntry is created using VectorEntry::toScalar() function and linked to the VectorEntry pointer scalar_entries_p. In that case, if keepVector is false the vector representation is destroyed
if unknowns is a scalar unknown, no new VectorEntry: scalar_entries_p = entries_p This function builds in any case, the cdofs_ numbering vector (vector of ComponentDof’s)
-
void toVector(bool keepEntries = false)#
-
return to vector representation from scalar representation
if vector unknown, create vector representation from scalar representation (have to exist) if scalar unknown, set entries_p = scalar_entries_p keepScalar: if false (default value), scalar entries are deleted
-
SuTermVector toVectorUnknown() const#
-
return vector unknown representation if component unknown SutermVector
when current SutermVector has a component unknown, return vector unknown representation if is not a component unknown SuTermVector, a copy of current SuTermVector is returned
-
SuTermVector(const LcTerm<SuTermVector>&)#