G+Smo
24.08.0
Geometry + Simulation Modules
|
Container class for a set of geometry patches and their topology, that is, the interface connections and outer boundary faces.
T | coefficient type |
Public Types | |
typedef memory::shared_ptr < gsMultiPatch > | Ptr |
Shared pointer for gsMultiPatch. | |
typedef memory::unique_ptr < gsMultiPatch > | uPtr |
Unique pointer for gsMultiPatch. | |
Public Member Functions | |
gsMatrix< index_t > | active (const gsMatrix< T > &u) const |
Returns the indices of active (nonzero) functions at points u, as a list of indices. More... | |
virtual void | active_into (const gsMatrix< T > &u, gsMatrix< index_t > &result) const |
Indices of active (non-zero) function(s) for each point. More... | |
void | addAutoBoundaries () |
Make all patch sides which are not yet declared as interface or boundary to a boundary. | |
void | addBoundary (index_t p, boxSide s, std::string l="") |
Set side s of box p to a boundary. | |
void | addBoundary (const patchSide &ps) |
Set patch side ps to a boundary. | |
void | addBox (index_t i=1) |
Add i new boxes. | |
template<class T > | |
gsProperty< T > | addBoxProperty (const std::string &name, T t=T()) |
Adds a box property. More... | |
void | addInterface (index_t p1, boxSide s1, index_t p2, boxSide s2, std::string l="") |
Add an interface between side s1 of box p1 and side s2 of box p2. | |
void | addInterface (const boundaryInterface &bi) |
Add an interface described by bi. | |
GISMO_DEPRECATED void | addInterface (gsGeometry< T > *g1, boxSide s1, gsGeometry< T > *g2, boxSide s2) |
Add an interface joint between side s1 of geometry g1 side s2 of geometry g2. More... | |
index_t | addPatch (typename gsGeometry< T >::uPtr g) |
Add a patch from a gsGeometry<T>::uPtr. | |
index_t | addPatch (const gsGeometry< T > &g) |
Add a patch by copying argument. | |
void | addPatchBoundary (gsGeometry< T > *g, boxSide s) |
Add side s of patch g to the outer boundary of the domain. | |
std::vector< std::vector < patchComponent > > | allComponents (bool combineCorners=false) const |
Returns all components representing the topology. More... | |
gsMultiPatch< T > | approximateLinearly (index_t nsamples) const |
Computes linear approximation of the patches using nsamples per direction. | |
std::vector< gsBasis< T > * > | basesCopy (bool numeratorOnly=false) const |
Makes a deep copy of all bases and puts them in a vector. More... | |
gsBasis< T > & | basis (const size_t i) const |
Return the basis of the i-th patch. | |
const gsBasis< T > & | basis (const index_t k) const |
Helper which casts and returns the k-th piece of this function set as a gsBasis. | |
const_biterator | bBegin () const |
biterator | bBegin () |
const_iterator | begin () const |
iterator | begin () |
const_biterator | bEnd () const |
biterator | bEnd () |
const boundaryInterface & | bInterface (int i) const |
Access i-th boundary interface. | |
const bContainer & | boundaries () const |
Return the vector of boundaries. | |
bContainer | boundaries (const std::string l) const |
Return the vector of boundaries with label l. | |
void | boundingBox (gsMatrix< T > &result) const |
Returns a bounding box for the multipatch domain. The output result is a matrix with two columns, corresponding to two points: the lower and upper corner of the bounding box. | |
void | checkConsistency () const |
Check that boundaries and interfaces are consistent. | |
void | clear () |
Clear (delete) all patches. | |
void | clearAll () |
Clear all boxes, boundary and interface data. | |
void | clearTopology () |
Clear all boundary and interface data. | |
gsBoxTopology * | clone () const |
Clone function. Used to make a copy of the object. | |
uPtr | clone () |
Clone methode. Produceds a deep copy inside a uPtr. | |
void | closeGaps (T tol=1e-4) |
Attempt to close gaps between the interfaces. Assumes that the topology is computed, ie. computeTopology() has been called. | |
short_t | coDim () const |
Co-dimension of the geometry (must match for all patches). | |
gsMatrix< T > | coefs () const |
Returns the coefficient matrix of the multi-patch geometry. | |
index_t | coefsSize () const |
Return the number of coefficients (control points) | |
virtual void | compute (const gsMatrix< T > &in, gsFuncData< T > &out) const |
Computes function data. More... | |
bool | computeTopology (T tol=1e-4, bool cornersOnly=false, bool tjunctions=false) |
Attempt to compute interfaces and boundaries automatically. More... | |
void | constructBoundaryRep () |
Construct the boundary representation. | |
void | constructBoundaryRep (const std::string l) |
Construct the boundary representation of sides with label l. | |
void | constructInterfaceRep (const std::string l) |
Construct the interface representation of sides with label l. | |
void | degreeElevate (short_t const elevationSteps=1, short_t const dir=-1) |
Elevate the degree of all patches by elevationSteps, preserves smoothness. | |
void | degreeIncrease (short_t const elevationSteps=1, short_t const dir=-1) |
Increase the degree of all patches by elevationSteps, preserves multiplicity. | |
void | degreeReduce (int elevationSteps=1) |
Reduce the degree of all patches by elevationSteps. | |
gsMatrix< T > | deriv (const gsMatrix< T > &u) const |
Evaluate the derivatives,. More... | |
gsMatrix< T > | deriv2 (const gsMatrix< T > &u) const |
Evaluates the second derivatives of active (i.e., non-zero) functions at points u. More... | |
virtual void | deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Second derivatives. More... | |
virtual void | deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
First derivatives. More... | |
std::string | detail () const |
Prints the object as a string with extended details. | |
short_t | dim () const |
Dimension of the boxes. | |
short_t | domainDim () const |
Dimension of the (source) domain. More... | |
bool | empty () const |
Returns true if gsMultiPatch is empty. | |
const_iterator | end () const |
iterator | end () |
gsMatrix< T > | eval (const gsMatrix< T > &u) const |
Evaluate the function,. More... | |
virtual void | eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluates the function(s). More... | |
std::vector< gsMatrix< T > > | evalAllDers (const gsMatrix< T > &u, int n, bool sameElement=false) const |
Evaluate all derivatives upto order n,. More... | |
virtual void | evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result, bool sameElement=false) const |
Evaluate the nonzero functions and their derivatives up to order n at points u into result. More... | |
InterfacePtr | findInterface (const index_t b1, const index_t b2) const |
size_t | findPatchIndex (gsGeometry< T > *g) const |
Search for the given geometry and return its patch index. | |
void | fixOrientation () |
Provides positive orientation for all patches. | |
const gsFunction< T > & | function (const index_t k) const |
Helper which casts and returns the k-th piece of this function set as a gsFunction. | |
short_t | geoDim () const |
Dimension of the geometry (must match for all patches). | |
BoundaryPtr | getBoundary (const patchSide &ps) |
Returns a pointer to the boundary stored in this class, given ps. | |
template<class T > | |
gsProperty< T > | getBoxProperty (const std::string &name) const |
Adds a box property. More... | |
bool | getCornerList (const patchCorner &start, std::vector< patchCorner > &cornerList) const |
void | getEVs (std::vector< std::vector< patchCorner > > &cornerLists, bool boundaries=false) const |
bool | getInterface (const patchSide &ps, boundaryInterface &result) const |
gsAffineFunction< T > | getMapForInterface (const boundaryInterface &bi, T scaling=0) const |
construct the affine map that places bi.first() next to bi.second() and identifies the two matching sides. The optional scaling specifies the scaling of first() in the direction orthogonal to the interface. By default it preserves the size of bi.first() in that direction More... | |
gsDofMapper | getMapper (T tol) const |
Used to get a mapper with unique vertices. | |
int | getMaxValence () const |
returns the maximal valence of a vertex of this topology. | |
bool | getNeighbour (const patchSide &ps, patchSide &result, int &ii) const |
bool | getNeighbour (const patchSide &ps, patchSide &result) const |
void | getOVs (std::vector< std::vector< patchCorner > > &cornerLists) const |
gsMultiPatch () | |
Default empty constructor. | |
gsMultiPatch (const gsMultiPatch &other) | |
Copy constructor (makes deep copy) | |
gsMultiPatch (PatchContainer &patches) | |
Create from a vector of patches. | |
gsMultiPatch (const gsGeometry< T > &geo) | |
Create a single-patch instance. | |
gsMultiPatch (PatchContainer &patches, const std::vector< patchSide > &boundary, const std::vector< boundaryInterface > &interfaces) | |
Create from patches and boundary/interface information. | |
std::vector< T > | HausdorffDistance (const gsMultiPatch< T > &other, const index_t nsamples=1000, const T accuracy=1e-6, const bool directed=false) |
Construct the interface representation. | |
const_iiterator | iBegin () const |
iiterator | iBegin () |
const_iiterator | iEnd () const |
iiterator | iEnd () |
const ifContainer & | interfaces () const |
Return the vector of interfaces. | |
ifContainer | interfaces (const std::string l) const |
Return the vector of interfaces with label l. | |
bool | isBoundary (const patchSide &ps) const |
Is the given patch side ps set to a boundary? | |
bool | isBoundary (index_t p, boxSide s) |
Returns true if side s on patch p is a boundary. | |
bool | isClosed () |
Returns true if the multipatch object is a closed manifold (ie. it has no boundaries) | |
bool | isInterface (const patchSide &ps) const |
Is the given patch side ps set to an interface? | |
void | locatePoints (const gsMatrix< T > &points, gsVector< index_t > &pids, gsMatrix< T > &preim, const T accuracy=1e-6) const |
For each point in points, locates the parametric coordinates of the point. More... | |
void | locatePoints (const gsMatrix< T > &points, index_t pid1, gsVector< index_t > &pid2, gsMatrix< T > &preim) const |
For each point in points located on patch pid1, locates the parametric coordinates of the point. More... | |
size_t | nBoundary () const |
Number of boundaries. | |
index_t | nBoxes () const |
Number of boxes. | |
size_t | nInterfaces () const |
Number of interfaces. | |
size_t | nPatches () const |
Number of patches. | |
index_t | nPieces () const |
Number of pieces in the domain of definition. | |
index_t | numBoxProperties () const |
Returns the number of assigned box properties. More... | |
gsMultiPatch & | operator= (gsMultiPatch other) |
Assignment operator (uses copy-and-swap idiom) | |
const gsGeometry< T > & | operator[] (size_t i) const |
Return the i-th patch. | |
gsMatrix< T > | parameterRange (int i=0) const |
Returns the range of parameter. | |
short_t | parDim () const |
Dimension of the parameter domain (must match for all patches). | |
gsGeometry< T > & | patch (size_t i) const |
Return the i-th patch. | |
PatchContainer const & | patches () const |
Returns a vector of patches // to do : replace by copies. | |
void | permute (const std::vector< short_t > &perm) |
Permutes the patches according to perm. | |
const gsGeometry< T > & | piece (const index_t i) const |
Returns the piece(s) of the function(s) at subdomain k. | |
gsMatrix< T > | pointOn (const patchCorner &pc) |
Get coordinates of the patchCorner pc in the physical domain. | |
gsMatrix< T > | pointOn (const patchSide &ps) |
Get coordinates of a central point of the the patchSide ps in the physical domain. More... | |
std::ostream & | print (std::ostream &os) const |
Prints the object as a string. | |
bool | repairInterface (const boundaryInterface &bi) |
Checks if the interface is fully matching, and if not, repairs it. More... | |
void | repairInterfaces () |
Checks if all patch-interfaces are fully matching, and if not, repairs them, i.e., makes them fully matching. More... | |
void | setDim (short_t i) |
Set the dimension of the boxes. | |
index_t | size () const |
size More... | |
void | swap (gsMultiPatch &other) |
Swap with another gsMultiPatch. | |
void | swap (gsBoxTopology &other) |
Swap with another gsBoxTopology. | |
short_t | targetDim () const |
Dimension of the target space. More... | |
gsSurfMesh | toMesh () const |
Creates a surface mesh out of this multipatch. | |
void | uniformCoarsen (int numKnots=1) |
Coarsen uniformly all patches by removing numKnots in each knot-span. | |
void | uniformRefine (int numKnots=1, int mul=1) |
Refine uniformly all patches by inserting numKnots in each knot-span with multipliplicity mul. | |
gsMultiPatch< T > | uniformSplit (index_t dir=-1) const |
Splits each patch uniformly in each direction (if dir = -1) into two new patches, giving a total number of 2^d new patches per patch. If dir is a parametric direction, then it only splits in that one direction. This method allocated new space for each new geometry, the original one stays unchanged. | |
~gsMultiPatch () | |
Destructor. | |
Protected Attributes | |
bContainer | m_boundary |
List of boundaries of the boxes. | |
gsProperty_container | m_boxProp |
List of properties for each box. | |
short_t | m_dim |
Dimension of the boxes held. | |
ifContainer | m_interfaces |
List of intefaces between boxes. | |
index_t | nboxes |
Number of boxes held. | |
Returns the indices of active (nonzero) functions at points u, as a list of indices.
Indices of active (non-zero) function(s) for each point.
The columns are sorted in increasing order, if on a point there are less active then the number of rows in the result matrix (some other point has more actives) then the rest of the column is filled with 0s.
u | |
result |
Reimplemented in gsHTensorBasis< d, T >, gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBSplineBasis< 1, T >, gsTensorBSplineBasis< d, T >, gsTensorBSplineBasis< domainDim+1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTHBSplineBasis< d, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsFunction< T >, gsConstantBasis< T >, gsLegendreBasis< T >, gsMonomialBasis< T >, and gsMvLegendreBasis< T >.
|
inlineinherited |
Adds a box property.
[in] | name | The name of the property |
[in] | t | The value of the property |
T | The type of the property |
void addInterface | ( | gsGeometry< T > * | g1, |
boxSide | s1, | ||
gsGeometry< T > * | g2, | ||
boxSide | s2 | ||
) |
Add an interface joint between side s1 of geometry g1 side s2 of geometry g2.
|
inherited |
Returns all components representing the topology.
Each entry of the outer vector represents one component (patch-interior, face, edge, corner, etc.). Since the components refering to one interface can be addressed as belonging to different patches, each component itself is represented by an inner vector which contains all patchComponent objects that refer to the particular component.
combineCorners | If this is set, all corners are treated as one component |
std::vector< gsBasis< T > * > basesCopy | ( | bool | numeratorOnly = false | ) | const |
Makes a deep copy of all bases and puts them in a vector.
numeratorOnly | If true, and the bases are derived from gsRationalBasis, then only the source bases (numerators) are returned |
|
inlineinherited |
Get a const-iterator to the beginning of the boundaries
|
inlineinherited |
Get an iterator to the beginning of the boundaries
|
inline |
Get a const-iterator to the patches
|
inline |
Get an iterator to the beginning of the patches
|
inlineinherited |
Get a const-iterator to the end of the boundaries
|
inlineinherited |
Get an iterator to the end of the boundaries
|
virtualinherited |
Computes function data.
This function evaluates the functions and their derivatives at the points in and writes them in the corresponding fields of out. Which field to write (and what to compute) is controlled by the out.flags (see also gsFuncData).
The input points in are expected to be compatible with the implementation/representation of the function, i.e. they should be points inside the domain of definitition of the function
[in] | in | |
[out] | out |
Reimplemented in gsGeometry< T >, and gsConstantFunction< T >.
bool computeTopology | ( | T | tol = 1e-4 , |
bool | cornersOnly = false , |
||
bool | tjunctions = false |
||
) |
Attempt to compute interfaces and boundaries automatically.
tol | The tolerance to test for matching points |
cornersOnly | When set to true an interface is accepted if the patch corners match, even if the parameterization does not agree |
Evaluate the derivatives,.
Evaluates the second derivatives of active (i.e., non-zero) functions at points u.
See documentation for deriv2_into() (the one without input parameter coefs) for details.
[in] | u | Evaluation points in columns. |
Second derivatives.
For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:
\[ \left[ \begin{array}{ccccc} \partial_{1}\partial_{1}f_1(p_1) & \partial_{1}\partial_{1}f_1(p_2) & \ldots & \partial_{1}\partial_{1}f_1(p_N)\\ \partial_{2}\partial_{2}f_1(p_1) & \partial_{2}\partial_{2}f_1(p_2) & \ldots & \partial_{2}\partial_{2}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}\partial_{k}f_1(p_1) & \partial_{k}\partial_{k}f_1(p_2) & \ldots & \partial_{k}\partial_{k}f_1(p_N)\\ \partial_{1}\partial_{2}f_1(p_1) & \partial_{1}\partial_{2}f_1(p_2) & \ldots & \partial_{1}\partial_{2}f_1(p_N)\\ \partial_{1}\partial_{3}f_1(p_1) & \partial_{1}\partial_{3}f_1(p_2) & \ldots & \partial_{1}\partial_{3}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{1}\partial_{k}f_1(p_1) & \partial_{1}\partial_{k}f_1(p_2) & \ldots & \partial_{1}\partial_{k}f_1(p_N)\\ \partial_{2}\partial_{3}f_1(p_1) & \partial_{2}\partial_{3}f_1(p_2) & \ldots & \partial_{2}\partial_{3}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{2}\partial_{k}f_1(p_1) & \partial_{2}\partial_{k}f_1(p_2) & \ldots & \partial_{2}\partial_{k}f_1(p_N)\\ \partial_{3}\partial_{4}f_1(p_1) & \partial_{3}\partial_{4}f_1(p_2) & \ldots & \partial_{3}\partial_{4}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1(p_1) & \partial_{k-1}\partial_{k}f_1(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1(p_N)\\ \partial_{1}\partial_{1}f_2(p_1) & \partial_{1}\partial_{1}f_2(p_2) & \ldots & \partial_{1}\partial_{1}f_2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_S(p_1) & \partial_{k-1}\partial_{k}f_S(p_2) & \ldots & \partial_{k-1}\partial_{k}f_S(p_N)\\ \end{array} \right] \]
For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^{m}\) the format is:
\[ \left[ \begin{array}{ccccc} \partial_{1}\partial_{1}f_1^{(1)}(p_1) & \partial_{1}\partial_{1}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{1}f_1^{(1)}(p_N)\\ \partial_{2}\partial_{2}f_1^{(1)}(p_1) & \partial_{2}\partial_{2}f_1^{(1)}(p_2) & \ldots & \partial_{2}\partial_{2}f_1^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}\partial_{k}f_1^{(1)}(p_1) & \partial_{k}\partial_{k}f_1^{(1)}(p_2) & \ldots & \partial_{k}\partial_{k}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{2}f_1^{(1)}(p_1) & \partial_{1}\partial_{2}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{2}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{3}f_1^{(1)}(p_1) & \partial_{1}\partial_{3}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{3}f_1^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1{(1)}(p_1) & \partial_{k-1}\partial_{k}f_1^{(1)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{1}f_1^{(2)}(p_1) & \partial_{1}\partial_{1}f_1^{(2)}(p_2) & \ldots & \partial_{1}\partial_{1}f_1^{(2)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1^{(m)}(p_1) & \partial_{k-1}\partial_{k}f_1^{(m)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1^{(m)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_S^{(m)}(p_1) & \partial_{k-1}\partial_{k}f_S^{(m)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_S^{(m)}(p_N)\\ \end{array} \right] \]
where \( f^{(i)}_j\) is the \(i\)-th component of function \( f_j\) of the set.
u | |
result |
Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsGeometry< T >, gsFunctionExpr< T >, gsFunction< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsPreCICEFunction< T >, gsLegendreBasis< T >, gsAffineFunction< T >, gsHBSplineBasis< d, T >, gsMappedSingleSpline< d, T >, gsMonomialBasis< T >, gsMvLegendreBasis< T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.
First derivatives.
For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:
\[ \left[ \begin{array}{ccccc} \partial_{1}f_1(p_1) & \partial_{1}f_1(p_2) & \ldots & \partial_{1}f_1(p_N)\\ \partial_{2}f_1(p_1) & \partial_{2}f_1(p_2) & \ldots & \partial_{2}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1(p_1) & \partial_{k}f_1(p_2) & \ldots & \partial_{k}f_1(p_N)\\ \partial_{1}f_2(p_1) & \partial_{1}f_2(p_2) & \ldots & \partial_{1}f_2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_S(p_1) & \partial_{k}f_S(p_2) & \ldots & \partial_{k}f_S(p_N)\\ \end{array} \right] \]
For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^{m}\) the format is:
\[ \left[ \begin{array}{ccccc} \partial_{1}f_1^1(p_1) & \partial_{1}f_1^1(p_2) & \ldots & \partial_{1}f_1^1(p_N)\\ \partial_{2}f_1^1(p_1) & \partial_{1}f_2^1(p_2) & \ldots & \partial_{1}f_2^1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1^1(p_1) & \partial_{k}f_1^1(p_2) & \ldots & \partial_{k}f_1^1(p_N)\\ \partial_{1}f_1^2(p_1) & \partial_{1}f_1^2(p_2) & \ldots & \partial_{1}f_1^2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1^2(p_1) & \partial_{k}f_1^2(p_2) & \ldots & \partial_{k}f_1^2(p_N)\\ \partial_{1}f_2^1(p_1) & \partial_{1}f_2^1(p_2) & \ldots & \partial_{1}f_2^1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_S^{(m)}(p_1) & \partial_{k}f_S^{(m)}(p_2) & \ldots & \partial_{k}f_S^{(m)}(p_N) \end{array} \right] \]
where \(f^{(i)}_j\) is the \(i\)-th component of function \(f_j\) of the set.
u | |
result |
Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTensorBSplineBasis< 1, T >, gsGeometry< T >, gsFunctionExpr< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsFunction< T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsHBSplineBasis< d, T >, gsLegendreBasis< T >, gsMappedSingleSpline< d, T >, gsMvLegendreBasis< T >, gsMonomialBasis< T >, gsConstantBasis< T >, gsGeometrySlice< T >, gsSquaredDistance< T >, gsBasisFun< T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.
|
inlinevirtual |
Dimension of the (source) domain.
Implements gsFunctionSet< T >.
|
inline |
Get a const iterator to the end of the patches
|
inline |
Get an iterator to the end of the patches
Evaluate the function,.
Evaluates the function(s).
For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:
\[ \left[ \begin{array}{ccccc} f_1(p_1) & f_1(p_2) & \ldots & f_1(p_N)\\ f_2(p_1) & f_2(p_2) & \ldots & f_2(p_N)\\ \vdots & \vdots & & \vdots\\ f_S(p_1) & f_S(p_2) & \ldots & f_S(p_N) \end{array} \right] \]
For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^m\) the format is:
\[ \left[ \begin{array}{ccccc} f_1^1(p_1) & f_1^{(1)}(p_2) & \ldots & f_1^{(1)}(p_N)\\ f_1^2(p_1) & f_1^{(2)}(p_2) & \ldots & f_1^{(2)}(p_N)\\ \vdots & \vdots & & \vdots\\ f_1^{(m)}(p_1) & f_1^{(m)}(p_2) & \ldots & f_1^{(m)}(p_N)\\ f_2^{(1)}(p_1) & f_2^{(1)}(p_2) & \ldots & f_2^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ f_S^{(m)}(p_1) & f_S^{(m)}(p_2) & \ldots & f_S^{(m)}(p_N) \end{array} \right] \]
where \(f^{(i)}_j\) is the \(i\)-th component of function \(f_j\) of the set.
u | |
result |
Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTHBSplineBasis< d, T >, gsPatchIdField< T >, gsMaterialMatrixEvalSingle< T, out >, gsParamField< T >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsShellStressFunction< T >, gsMaterialMatrixIntegrateSingle< T, out >, gsNormalField< T >, gsFsiLoad< T >, gsFunctionExpr< T >, gsGeometry< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsPiecewiseFunction< T >, gsConstantFunction< T >, gsJacDetField< T >, gsDetFunction< T >, gsRemapInterface< T >, gsFunction< T >, gsGradientField< T >, gsCPPInterface< T >, gsHBSplineBasis< d, T >, gsAffineFunction< T >, gsPreCICEFunction< T >, gsLegendreBasis< T >, gsMvLegendreBasis< T >, gsMappedSingleSpline< d, T >, gsConstantBasis< T >, gsMonomialBasis< T >, gsCauchyStressFunction< T >, gsBasisFun< T >, gsGeometrySlice< T >, gsGeometryTransform< T >, gsFuncCoordinate< T >, gsAbsError< T >, gsElementErrorPlotter< T >, and gsSquaredDistance< T >.
|
inherited |
Evaluate all derivatives upto order n,.
|
virtualinherited |
Evaluate the nonzero functions and their derivatives up to order n at points u into result.
The derivatives (the 0-th derivative is the function value) are stored in a result. result is a std::vector, where result[i] is a gsMatrix which contains the i-th derivatives.
The entries in result[0], result[1], and result[2] are ordered as in eval_into(), deriv_into(), and deriv2_into(), respectively. For i > 2, the derivatives are stored in lexicographical order, e.g. for order i = 3 and dimension 2 the derivatives are stored as follows: \( \partial_{xxx}, \, \partial_{xxy}, \, \partial_{xyy}, \, \partial_{yyy}.\, \)
[in] | u | Evaluation points, each column corresponds to one evaluation point. |
[in] | n | All derivatives up to order n are computed and stored in result. |
[in,out] | result | See above for format. |
Reimplemented in gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsGeometry< T >, gsMappedSingleBasis< d, T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsMappedSingleSpline< d, T >, and gsSquaredDistance< T >.
|
inherited |
Returns a pointer to the interface between boxes b1 and b2, if one exists, otherwise it returns a null pointer
|
inlineinherited |
Adds a box property.
[in] | name | The name of the property |
T | The type of the property |
|
inherited |
takes a patchCorner start and gives back all other patchCorners, that represent the same point in the vector cornerList
|
inherited |
gives back all the extraordinary vertices (3 faces or more than 4) of the topology each EV is represented by a vector of patchCorners, which represent the same vertex all the vectors are put in the vector cornerLists. It will only find vertices on the inside. CAREFUL: works only for 2D
|
inlineinherited |
set result to the associated interface of ps, returns false if it is a boundary patchSide
gsAffineFunction< T > getMapForInterface | ( | const boundaryInterface & | bi, |
T | scaling = 0 |
||
) | const |
construct the affine map that places bi.first() next to bi.second() and identifies the two matching sides. The optional scaling specifies the scaling of first() in the direction orthogonal to the interface. By default it preserves the size of bi.first() in that direction
bi | |
scaling |
|
inherited |
gives back all the ordinary vertices (4 faces) of the topology each OV is represented by a vector of patchCorners, which represent the same vertex all the vectors are put in the vector cornerLists It will only find vertices on the inside. CAREFUL: works only for 2D
|
inlineinherited |
Get a const-iterator to the interfaces
|
inlineinherited |
Get an iterator to the beginning of the interfaces
|
inlineinherited |
Get a const iterator to the end of the interfaces
|
inlineinherited |
Get an iterator to the end of the interfaces
void locatePoints | ( | const gsMatrix< T > & | points, |
gsVector< index_t > & | pids, | ||
gsMatrix< T > & | preim, | ||
const T | accuracy = 1e-6 |
||
) | const |
For each point in points, locates the parametric coordinates of the point.
points | |
pids | vector containing for each point the patch id where it belongs (or -1 if not found) |
preim | in each column, the parametric coordinates of the corresponding point in the patch |
void locatePoints | ( | const gsMatrix< T > & | points, |
index_t | pid1, | ||
gsVector< index_t > & | pid2, | ||
gsMatrix< T > & | preim | ||
) | const |
For each point in points located on patch pid1, locates the parametric coordinates of the point.
pid2 | vector containing for each point the patch id where it belongs (or -1 if not found) |
preim | in each column, the parametric coordinates of the corresponding point in the patch |
|
inlineinherited |
Returns the number of assigned box properties.
Get coordinates of a central point of the the patchSide ps in the physical domain.
The central point in the physical domain is the the midpoint on the parameter domain, mapped to the physical domain
bool repairInterface | ( | const boundaryInterface & | bi | ) |
Checks if the interface is fully matching, and if not, repairs it.
|
inline |
Checks if all patch-interfaces are fully matching, and if not, repairs them, i.e., makes them fully matching.
|
inlinevirtual |
size
Reimplemented from gsFunctionSet< T >.
|
inlinevirtual |
Dimension of the target space.
Reimplemented from gsFunctionSet< T >.