65 typedef memory::shared_ptr< gsFunction >
Ptr;
68 typedef memory::unique_ptr< gsFunction >
uPtr;
78 GISMO_ENSURE(0==k,
"Single function of type "<<
typeid(*this).name() <<
" is defined on single subdomain, received: "<<k<<
". Is piece(k) implemented?" );
86 { result.setConstant(1,u.cols(),0); }
207 hessian_into(u,res,
coord);
227 const T accuracy = 1e-6,
228 const bool useInitialPoint =
false)
const;
230 virtual void invertPointGrid(gsGridIterator<T,0> & git,
232 const bool useInitialPoint =
false)
const;
240 bool withSupport =
true,
241 const T accuracy = 1e-6,
243 T damping_factor = 1)
const;
248 T damping_factor = 1)
const;
254 const T accuracy = 1e-6)
const;
256 void recoverPointGrid(gsGridIterator<T,0> & git,
258 index_t k,
const T accuracy = 1e-6)
const;
265 return ( S.col(0) + S.col(1) ) * (T)(0.5);
276 virtual std::ostream &
print(std::ostream &os)
const
278 os <<
"gsFunction.\n";
return os;
299 template<
int mode,
int _Dim=-1>
300 int newtonRaphson_impl(
303 const T accuracy = 1e-6,
int max_loop = 100,
304 T damping_factor = 1, T scale = 1.0)
const;
316{
return b.
print(os); }
318#ifdef GISMO_WITH_PYBIND11
323 void pybind11_init_gsFunction(pybind11::module &m);
330#ifndef GISMO_BUILD_LIB
331#include GISMO_HPP_HEADER(gsFunction.hpp)
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
Represents a certain component of a vector-valued function.
Definition gsFuncCoordinate.h:35
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition gsFunctionSet.h:219
virtual short_t domainDim() const =0
Dimension of the (source) domain.
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
virtual short_t targetDim() const
Dimension of the target space.
Definition gsFunctionSet.h:595
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
virtual void deriv2_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate second derivatives of the function at points u into result.
Definition gsFunction.hpp:130
virtual gsMatrix< T > laplacian(const gsMatrix< T > &u) const
Evaluate the Laplacian at points u.
Definition gsFunction.hpp:177
virtual T distanceL2(gsFunction< T > const &) const
Computes the L2-distance between this function and the field and a function func.
Definition gsFunction.h:220
gsFuncCoordinate< T > coord(const index_t c) const
Returns the scalar function giving the i-th coordinate of this function.
Definition gsFunction.hpp:38
void recoverPoints(gsMatrix< T > &xyz, gsMatrix< T > &uv, index_t k, const T accuracy=1e-6) const
Definition gsFunction.hpp:504
memory::unique_ptr< gsFunction > uPtr
Unique pointer for gsFunction.
Definition gsFunction.h:68
virtual gsMatrix< T > hessian(const gsMatrix< T > &u, index_t coord=0) const
Definition gsFunction.h:204
index_t size() const
size
Definition gsFunction.h:295
memory::shared_ptr< gsFunction > Ptr
Shared pointer for gsFunction.
Definition gsFunction.h:65
void active_into(const gsMatrix< T > &u, gsMatrix< index_t > &result) const
Indices of active (non-zero) function(s) for each point.
Definition gsFunction.h:85
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.
Definition gsFunction.hpp:53
void div_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Computes for each point u a block of result containing the divergence matrix.
Definition gsFunction.hpp:65
virtual void computeMap(gsMapData< T > &InOut) const
Computes map function data.
Definition gsFunction.hpp:817
virtual void invertPoints(const gsMatrix< T > &points, gsMatrix< T > &result, const T accuracy=1e-6, const bool useInitialPoint=false) const
Definition gsFunction.hpp:199
virtual std::ostream & print(std::ostream &os) const
Prints the object as a string.
Definition gsFunction.h:276
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.
Definition gsFunction.hpp:546
virtual const gsFunction & piece(const index_t k) const
Returns the piece(s) of the function(s) at subdomain k.
Definition gsFunction.h:76
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
Definition gsFunction.hpp:443
virtual void deriv_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate derivatives of the function at points u into result.
Definition gsFunction.hpp:93
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const =0
Evaluate the function at points u into result.
virtual gsMatrix< T > parameterCenter() const
Returns a "central" point inside inside the parameter domain.
Definition gsFunction.h:261
the gsMapData is a cache of pre-computed function (map) values.
Definition gsFuncData.h:349
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
This is the interface of all objects that computes functions on points like gsBasis,...
The G+Smo namespace, containing all definitions for the library.
Struct which represents a certain corner of a hyper-cube.
Definition gsBoundary.h:292