G+Smo  23.12.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsRemapInterface< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsRemapInterface< T >

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.

+ Inheritance diagram for gsRemapInterface< T >:
+ Collaboration diagram for gsRemapInterface< T >:

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_tactive (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) basis 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) const
 Evaluate all derivatives upto order n,. More...
 
virtual void evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result) const
 Evaluate the nonzero functions and their derivatives up to order n. If n is -1 then no computation is performed.
 
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 >::PtrinterfaceMap () 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 gsFunctionpiece (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.
 
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.
 
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.
 
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...
 

Member Enumeration Documentation

anonymous enum
Enumerator
neverAffine 

Instructs constructor not to use an affine mapping.

alwaysAffine 

Instructs constructor always to use an affine mapping.

Constructor & Destructor Documentation

gsRemapInterface ( const gsMultiPatch< T > &  mp,
const gsMultiBasis< T > &  mb,
const boundaryInterface bi,
const gsOptionList opt = defaultOptions() 
)

Constructor.

Parameters
mpThe multi-patch object
mbThe multi-basis object
biThe boundary interface (specifying the patches and their patchSide s)
optOptions; 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

Member Function Documentation

gsMatrix<index_t> active ( const gsMatrix< T > &  u) const
inlineinherited

Returns the indices of active (nonzero) functions at points u, as a list of indices.

See Also
active_into()
void active_into ( const gsMatrix< T > &  u,
gsMatrix< index_t > &  result 
) const
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.

Parameters
u
result

Reimplemented from gsFunctionSet< T >.

void compute ( const gsMatrix< T > &  in,
gsFuncData< T > &  out 
) const
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

Parameters
[in]in
[out]out

Reimplemented in gsGeometry< T >, and gsConstantFunction< T >.

void computeMap ( gsMapData< T > &  InOut) const
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.

Parameters
[in,out]InOut
static gsOptionList defaultOptions ( )
inlinestatic

default options

See constructor for their meaning

gsMatrix< T > deriv ( const gsMatrix< T > &  u) const
inherited

Evaluate the derivatives,.

See Also
deriv_into()
gsMatrix< T > deriv2 ( const gsMatrix< T > &  u) const
inherited

Evaluates the second derivatives of active (i.e., non-zero) basis at points u.

See documentation for deriv2_into() (the one without input parameter coefs) for details.

See Also
deriv2_into()
Parameters
[in]uEvaluation points in columns.
Returns
For every column of u, a column containing the second derivatives. See documentation for deriv2_into() (the one without input parameter coefs) for details.
void deriv2_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

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.

Parameters
[in]ugsMatrix of size n x N, where each column of u represents one evaluation point.
[out]resultgsMatrix 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\)
Warning
By default uses central finite differences with h=0.00001! One must override this function in derived classes to get proper results.

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsGeometry< T >, gsFunctionExpr< T >, gsConstantFunction< T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsMappedSingleSpline< d, T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.

void deriv_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

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] \]

Parameters
[in]ugsMatrix of size n x N, where each column of u represents one evaluation point.
[out]resultgsMatrix 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).
Warning
By default, gsFunction uses central finite differences with h=0.00001! One must override this function in derived classes to get proper results.

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 >.

gsMatrix< T > eval ( const gsMatrix< T > &  u) const
inherited

Evaluate the function,.

See Also
eval_into()
void eval_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtual

Evaluates the interface map.

Interfaces map \( \widehat \Gamma_1 \rightarrow \widehat \Gamma_2 \) that represents \( G_2^{-1} \circ G_1 \)

Implements gsFunction< T >.

std::vector< gsMatrix< T > > evalAllDers ( const gsMatrix< T > &  u,
int  n 
) const
inherited

Evaluate all derivatives upto order n,.

See Also
evalAllDers_into
virtual gsMatrix<T> hessian ( const gsMatrix< T > &  u,
index_t  coord = 0 
) const
inlinevirtualinherited

Evaluates the Hessian (matrix of second partial derivatives) of coordinate coord at points u.

const gsFunction<T>::Ptr& interfaceMap ( ) const
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 \)

void invertPoints ( const gsMatrix< T > &  points,
gsMatrix< T > &  result,
const T  accuracy = 1e-6,
const bool  useInitialPoint = false 
) const
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

gsMatrix< T > laplacian ( const gsMatrix< T > &  u) const
virtualinherited

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 \).

int newtonRaphson ( const gsVector< T > &  value,
gsVector< T > &  arg,
bool  withSupport = true,
const T  accuracy = 1e-6,
int  max_loop = 100,
damping_factor = 1 
) const
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

void recoverPoints ( gsMatrix< T > &  xyz,
gsMatrix< T > &  uv,
index_t  k,
const T  accuracy = 1e-6 
) const
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).

index_t size ( ) const
inlinevirtualinherited

size

Warning
gsFunction and gsGeometry have size() == 1. This should not be confused with the size eg. of gsGeometry::basis(), which is the number of basis functions in the basis
Returns
the size of the function set: the total number of functions

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsPiecewiseFunction< T >.

Friends And Related Function Documentation

std::ostream & operator<< ( std::ostream &  os,
const gsRemapInterface< T > &  remapIf 
)
related

Prints the state of the object.

Member Data Documentation

gsMatrix<T> m_parameterBounds1
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} \]

gsMatrix<T> m_parameterBounds2
private

The bounds of the box that represents \( \widehat \Gamma_2 \).

See m_parameterBounds1