G+Smo
24.08.0
Geometry + Simulation Modules
|
Provides a mapping between the corresponding sides of two patches sharing an interface.
Have two patches \( \Omega_1 = G_1( \widehat \Omega_1 )\) and \( \Omega_2 = G_2( \widehat \Omega_2 )\) with geometry functions \( G_1 \) and \( G_2 \) and parameter domains \( \widehat \Omega_1 \) and \( \widehat \Omega_2 \).
This class represents an interface
\[ \Gamma := G_1( \widehat S_1 ) \cap G_2( \widehat S_2 ), \]
where \( \widehat S_1 \) and \( \widehat S_2 \) are sides of the corresponding parameter domains (i.e., whole edges or whole faces). The indices of the patches, as well as the corresponding sides (as boxSide) are prescribed by the incoming boundaryInterface.
The pre-images of the interface \( \Gamma \) are \( \widehat \Gamma_1 := G_1^{-1}( \Gamma ) \) and \( \widehat \Gamma_2 := G_2^{-1}( \Gamma ) \).
We say the interface is matching if
\[ \Gamma = G_1( \widehat S_1 ) = G_2( \widehat S_2 ). \]
In this case, the pre-images of the interface \( \Gamma \) are the whole sides, i.e., \( \widehat \Gamma_1 = \widehat S_1 \) and \( \widehat \Gamma_2 = \widehat S_2 \). Otherwise, the pre-images are only subsets of the sides.
We say the interface is affine if \( G_1^{-1} \circ G_2 \) or, euqivalently \( G_2^{-1} \circ G_1 \) is an affine-linear function. In this case, this class uses that affine-linear function.
If the mapping is found out not to be affine-linear, this class constructs an approximation of the interface map. This is available only if \( \Omega_1 \) and \( \Omega_2 \) are two dimensional, i.e., the interface itself is one dimensional. We refer to this as the 2D case.
Public Types | |
enum | { neverAffine, alwaysAffine } |
typedef memory::shared_ptr < gsRemapInterface > | Ptr |
Shared pointer for gsRemapInterface. | |
typedef memory::unique_ptr < gsRemapInterface > | uPtr |
Unique pointer for gsRemapInterface. | |
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... | |
void | active_into (const gsMatrix< T > &u, gsMatrix< index_t > &result) const |
Indices of active (non-zero) function(s) for each point. More... | |
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 std::vector< std::vector < T > > & | breakPoints () const |
Returns the break points used in makeDomainIterator. | |
uPtr | clone () |
Clone methode. Produceds a deep copy inside a uPtr. | |
virtual void | compute (const gsMatrix< T > &in, gsFuncData< T > &out) const |
Computes function data. More... | |
virtual void | computeMap (gsMapData< T > &InOut) const |
Computes map function data. More... | |
gsFuncCoordinate< T > | coord (const index_t c) const |
Returns the scalar function giving the i-th coordinate of this function. | |
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 T | distanceL2 (gsFunction< T > const &) const |
Computes the L2-distance between this function and the field and a function func. | |
virtual short_t | domainDim () const |
Returns parameter dimension of the domains. | |
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 interface map. 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... | |
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. | |
gsRemapInterface (const gsMultiPatch< T > &mp, const gsMultiBasis< T > &mb, const boundaryInterface &bi, const gsOptionList &opt=defaultOptions()) | |
Constructor. More... | |
const gsFunction< T >::Ptr & | interfaceMap () const |
Returns the interface map iff affine, otherwise the fitting curve. More... | |
virtual void | invertPoints (const gsMatrix< T > &points, gsMatrix< T > &result, const T accuracy=1e-6, const bool useInitialPoint=false) const |
bool | isAffine () const |
Returns true iff the interface is affine. | |
bool | isMatching () const |
Returns true iff the interface is matching. | |
gsDomainIterator< T >::uPtr | makeDomainIterator () const |
Returns a domain iterator. More... | |
int | newtonRaphson (const gsVector< T > &value, gsVector< T > &arg, bool withSupport=true, const T accuracy=1e-6, int max_loop=100, T damping_factor=1) const |
virtual index_t | nPieces () const |
Number of pieces in the domain of definition. | |
virtual gsMatrix< T > | parameterCenter () const |
Returns a "central" point inside inside the parameter domain. | |
gsMatrix< T > | parameterCenter (const boxCorner &bc) const |
Get coordinates of the boxCorner bc in the parameter domain. | |
gsMatrix< T > | parameterCenter (const boxSide &bs) const |
Get coordinates of the midpoint of the boxSide bs in the parameter domain. | |
virtual const gsFunction & | piece (const index_t k) const |
Returns the piece(s) of the function(s) at subdomain k. | |
virtual std::ostream & | print (std::ostream &os) const |
Prints the state of the object. | |
void | recoverPoints (gsMatrix< T > &xyz, gsMatrix< T > &uv, index_t k, const T accuracy=1e-6) const |
index_t | size () const |
size More... | |
virtual short_t | targetDim () const |
Dimension of the target space. More... | |
Evaluation functions | |
These functions allow one to evaluate the function as well as its derivatives at one or more points in the parameter space. See also Evaluation members. | |
virtual void | eval_component_into (const gsMatrix< T > &u, const index_t comp, gsMatrix< T > &result) const |
Evaluate the function for component comp in the target dimension at points u into result. | |
virtual void | deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluate derivatives of the function \(f:\mathbb{R}^n\rightarrow\mathbb{R}^m\) at points u into result. More... | |
virtual void | jacobian_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Computes for each point u a block of result containing the Jacobian matrix. | |
void | div_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Computes for each point u a block of result containing the divergence matrix. | |
gsMatrix< T > | jacobian (const gsMatrix< T > &u) const |
virtual void | deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const |
Evaluate second derivatives of the function at points u into result. More... | |
virtual void | hessian_into (const gsMatrix< T > &u, gsMatrix< T > &result, index_t coord=0) const |
virtual gsMatrix< T > | hessian (const gsMatrix< T > &u, index_t coord=0) const |
virtual gsMatrix< T > | laplacian (const gsMatrix< T > &u) const |
Evaluate the Laplacian at points u. More... | |
Static Public Member Functions | |
static gsOptionList | defaultOptions () |
Private Member Functions | |
void | constructBreaks () |
Constructs the breakpoints m_breakpoints. | |
void | constructFittingCurve (index_t numSamplePoints, index_t intervalsOfFittingCurve, index_t degreeOfFittingCurve) |
Constructs the fitting curve m_intfMap in the non-affine case. | |
void | constructInterfaceBox () |
Computes the box which represents the intersection of sides of incoming patches. | |
T | estimateReparamError (index_t steps) const |
Estimates error between geometry1 and geometry2 after repearametrization on physical domain. | |
Private Attributes | |
const gsBasis< T > * | m_b1 |
Basis on first patch. | |
const gsBasis< T > * | m_b2 |
Basis on second patch. | |
boundaryInterface | m_bi |
Corresponding boundary interface. | |
std::vector< std::vector< T > > | m_breakpoints |
Union of breakpoints of both bases. | |
T | m_equalityTolerance |
Tolerance for considering points to be equal. | |
const gsGeometry< T > * | m_g1 |
Geometry of first patch. | |
const gsGeometry< T > * | m_g2 |
Geometry of second patch. | |
gsFunction< T >::Ptr | m_intfMap |
Iff affine, interface map, otherwise the fitting curve. | |
bool | m_isAffine |
True iff the interface is affine. | |
bool | m_isMatching |
True iff the interface is matching. | |
T | m_newtonTolerance |
Tolerance for Newton solver. | |
gsMatrix< T > | m_parameterBounds1 |
The bounds of the box that represents \( \widehat \Gamma_1 \). More... | |
gsMatrix< T > | m_parameterBounds2 |
The bounds of the box that represents \( \widehat \Gamma_2 \). More... | |
Related Functions | |
(Note that these are not member functions.) | |
template<class T > | |
std::ostream & | operator<< (std::ostream &os, const gsRemapInterface< T > &remapIf) |
Prints the state of the object. More... | |
anonymous enum |
gsRemapInterface | ( | const gsMultiPatch< T > & | mp, |
const gsMultiBasis< T > & | mb, | ||
const boundaryInterface & | bi, | ||
const gsOptionList & | opt = defaultOptions() |
||
) |
Constructor.
mp | The multi-patch object |
mb | The multi-basis object |
bi | The boundary interface (specifying the patches and their patchSide s) |
opt | Options; see below: |
Option | Meaning |
---|---|
CheckAffine | The number of interior points (per direction) in point grid used to |
check if the mapping is affine. If set to i, the grid consits of i interior | |
points and the two boundary points per direction, so \( (i+2)^d \) points. | |
For alwaysAffine (=0) or neverAffine (-1), no checks are performed; default: 1 | |
NumSamplePoints | Number of sample points for computing fitting curve (if not affine); default: 11 |
IntervalsOfFittingCurve | Number of knot space in the fitting curve (if not affine); default: 5 |
DegreeOfFittingCurve | Spline degree of fitting curve (if not affine); default: 3 |
EqualityTolerance | Points are considered equal iff difference is smaller; default 1e-5 |
NewtonTolerance | Tolerance for Newton solvers (should be significantly smaller than EqualityTolerance); |
default: 1e-10 |
Returns the indices of active (nonzero) functions at points u, as a list of indices.
|
inlinevirtualinherited |
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 from gsFunctionSet< T >.
|
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 >.
|
virtualinherited |
Computes map function data.
This function evaluates the functions and their derivatives at the points InOut.points and writes them in the corresponding fields of InOut. Which field to write (and what to compute) is controlled by the InOut.flags (see also gsMapData). This is intended for parametrizations only and it works on functions sets of cardinality 1 only.
[in,out] | InOut |
|
inlinestatic |
default options
See constructor for their meaning
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. |
Evaluate second derivatives of the function at points u into result.
Let n be the dimension of the source space ( n = domainDim() ).
Let m be the dimension of the image/target space ( m = targetDim() ).
Let N denote the number of evaluation points.
[in] | u | gsMatrix of size n x N, where each column of u represents one evaluation point. |
[out] | result | gsMatrix of size (S*m) x N, where S=n*(n+1)/2. Each column in result corresponds to one point (i.e., one column in u) and contains the following values (for n=3, m=3): \( (\partial_{xx} f^{(1)}, \partial_{yy} f^{(1)}, \partial_{zz} f^{(1)}, \partial_{xy} f^{(1)}, \partial_{xz} f^{(1)}, \partial_{yz} f^{(1)}, \partial_{xx} f^{(2)},\ldots,\partial_{yz} f^{(3)} )^T\) |
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsGeometry< T >, gsFunctionExpr< T >, gsConstantFunction< T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsMappedSingleSpline< d, T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.
Evaluate derivatives of the function \(f:\mathbb{R}^n\rightarrow\mathbb{R}^m\) at points u into result.
Let n be the dimension of the source space ( n = domainDim() ).
Let m be the dimension of the image/target space ( m = targetDim() ).
Let N denote the number of evaluation points.
Let \( f:\mathbb R^2 \rightarrow \mathbb R^3 \), i.e., \( f(x,y) = ( f^{(1)}(x,y), f^{(2)}(x,y), f^{(3)}(x,y) )^T\),
and let \( u = ( u_1, \ldots, u_N) = ( (x_1,y_1)^T, \ldots, (x_N, y_N)^T )\).
Then, result is of the form
\[ \left[ \begin{array}{cccc} \partial_x f^{(1)}(u_1) & \partial_x f^{(1)}(u_2) & \ldots & \partial_x f^{(1)}(u_N) \\ \partial_y f^{(1)}(u_1) & \partial_y f^{(1)}(u_2) & \ldots & \partial_y f^{(1)}(u_N) \\ \partial_x f^{(2)}(u_1) & \partial_x f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_y f^{(2)}(u_1) & \partial_y f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_x f^{(3)}(u_1) & \partial_x f^{(3)}(u_2) & \ldots & \partial_x f^{(3)}(u_N)\\ \partial_y f^{(3)}(u_1) & \partial_y f^{(3)}(u_2) & \ldots & \partial_y f^{(3)}(u_N) \end{array} \right] \]
[in] | u | gsMatrix of size n x N, where each column of u represents one evaluation point. |
[out] | result | gsMatrix of size (n * m) x N. Each row of result corresponds to one component in the target space and contains the gradients for each evaluation point, as row vectors, one after the other (see above for details on the format). |
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsGeometry< T >, gsFunctionExpr< T >, gsConstantFunction< T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsMappedSingleSpline< d, T >, gsGeometrySlice< T >, gsSquaredDistance< T >, gsBasisFun< T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.
Evaluate the function,.
Evaluates the interface map.
Interfaces map \( \widehat \Gamma_1 \rightarrow \widehat \Gamma_2 \) that represents \( G_2^{-1} \circ G_1 \)
Implements gsFunction< 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 >.
|
inlinevirtualinherited |
Evaluates the Hessian (matrix of second partial derivatives) of coordinate coord at points u.
|
inline |
Returns the interface map iff affine, otherwise the fitting curve.
Interfaces map \( \widehat \Gamma_1 \rightarrow \widehat \Gamma_2 \) that represents \( G_2^{-1} \circ G_1 \)
|
virtualinherited |
Takes the physical points and computes the corresponding parameter values. If the point cannot be inverted (eg. is not part of the geometry) the corresponding parameter values will be undefined
Evaluate the Laplacian at points u.
By default uses central finite differences with h=0.00001
Reimplemented in gsFunctionExpr< T >.
gsDomainIterator< T >::uPtr makeDomainIterator | ( | ) | const |
Returns a domain iterator.
The domain iterator lives on \( \widehat \Gamma_1 \). Its break points are the union of the brakpoints of the basis on \( \widehat \Omega_1 \) and the breakpoints of the basis on \( \widehat \Omega_2 \), mapped to \( \widehat \Omega_1 \).
|
inherited |
Newton-Raphson method to find a solution of the equation f(arg) = value with starting vector arg. If the point cannot be inverted the corresponding parameter values will be undefined
|
inherited |
Recovers a point on the (geometry) together with its parameters uv, assuming that the k-th coordinate of the point xyz is not known (and has a random value as input argument).
|
inlinevirtualinherited |
size
Reimplemented from gsFunctionSet< T >.
Reimplemented in gsPiecewiseFunction< T >.
|
inlinevirtualinherited |
Dimension of the target space.
Reimplemented in gsPatchIdField< T >, gsGeometry< T >, gsParamField< T >, gsNormalField< T >, gsMultiBasis< T >, gsMultiBasis< real_t >, gsMultiPatch< T >, gsMultiPatch< real_t >, gsFsiLoad< T >, gsMaterialMatrixEvalSingle< T, out >, gsMaterialMatrixIntegrateSingle< T, out >, gsFunctionExpr< T >, gsJacDetField< T >, gsDetFunction< T >, gsConstantFunction< T >, gsShellStressFunction< T >, gsGradientField< T >, gsPiecewiseFunction< T >, gsSquaredDistance< T >, gsAffineFunction< T >, gsPreCICEFunction< T >, gsMappedSingleSpline< d, T >, gsBasisFun< T >, gsAbsError< T >, gsGeometrySlice< T >, gsGeometryTransform< T >, gsFuncCoordinate< T >, and gsCauchyStressFunction< T >.
|
related |
Prints the state of the object.
|
private |
The bounds of the box that represents \( \widehat \Gamma_1 \).
If \( \widehat \Gamma_1 = (a,b)\times (c,d) \times (e,f) \), the matrix has the structure
\[ \begin{pmatrix} a & c & e \\ b & d & f \end{pmatrix} \]
|
private |
The bounds of the box that represents \( \widehat \Gamma_2 \).