Medusa  1.1
Coordinate Free Mehless Method implementation
mm::UnknownShape< vec_t > Class Template Reference

#include <UnknownShape_fwd.hpp>

Detailed Description

template<typename vec_t>
class mm::UnknownShape< vec_t >

This class represents an unknown domain shape.

It can be used for domains when shape is not known or important, like in simple examples:

UnknownShape<Vec2d> sh;
DomainDiscretization<Vec2d> d(sh);
std::vector<Vec2d> pos = {{-2.3, 3.4}, {-2.1, 3.4}, {0, 1.4}, {-1.4, 0.001}};
std::vector<int> types = {2, -3, 1, -1};
std::vector<Vec2d> normals = {{-1.0, 0}, {0.0, 1.0}};
d.addInternalNode(pos[0], types[0]);
d.addBoundaryNode(pos[1], types[1], normals[0]);
d.addInternalNode(pos[2], types[2]);
d.addBoundaryNode(pos[3], types[3], normals[1]);

This is the shape used as a placeholder when loading domain discretizations.

// Load domain from file `test_domain_load.h5`. It is stored in the `domain` group.
// Shape of `d` is UnknownShape.
HDF file("test/testdata/test_domain_load.h5", HDF::READONLY);
auto d = DomainDiscretization<Vec2d>::load(file, "domain");

Definition at line 27 of file UnknownShape_fwd.hpp.

+ Inheritance diagram for mm::UnknownShape< vec_t >:
+ Collaboration diagram for mm::UnknownShape< vec_t >:

Public Member Functions

bool contains (const vec_t &) const override
 Return true if point is not more than margin() outside the domain. More...
 
bool hasContains () const override
 Return true if shape has contains() method implemented. More...
 
std::pair< vec_t, vec_t > bbox () const override
 Return the bounding box of the domain. More...
 
DomainDiscretization< vec_t > discretizeBoundaryWithDensity (const std::function< scalar_t(vec_t)> &, int) const override
 Discretizes boundary with given density and fill engine. More...
 
std::ostream & print (std::ostream &ostream) const override
 Output information about this shape to given output stream os. More...
 
UnknownShapeclone () const override
 Polymorphic clone pattern. More...
 
scalar_t margin () const
 Returns current margin. More...
 
virtual void setMargin (scalar_t margin)
 Sets domain margin to margin. More...
 
void toggleMargin ()
 Toggles the margin from positive to negative. More...
 
ShapeUnion< vec_t > add (const DomainShape &other) const
 Returns a shape representing a union of *this and other. More...
 
ShapeUnion< vec_t > operator+ (const DomainShape &other) const
 Operator form of DomainShape::add. More...
 
ShapeDifference< vec_t > subtract (const DomainShape &other) const
 Returns a shape representing a difference of *this and other. More...
 
ShapeDifference< vec_t > operator- (const DomainShape &other) const
 Operator form of DomainShape::subtract. More...
 
virtual std::pair< bool, vec_t > projectPointToBoundary (const vec_t &point, const vec_t &unit_normal) const
 Project point to boundary using bisection along the line define by unit_normal. More...
 
virtual DomainDiscretization< vec_t > discretizeBoundaryWithStep (scalar_t step, int type) const
 Returns a discretization of the boundary of this shape with approximately uniform distance step between nodes. More...
 
DomainDiscretization< vec_t > discretizeBoundaryWithStep (scalar_t step) const
 Returns a discretization of the boundary of this shape with approximately uniform distance step between nodes. More...
 
virtual DomainDiscretization< vec_t > discretizeWithStep (scalar_t step, int internal_type, int boundary_type) const
 Returns a discretization of this shape with approximately uniform distance step between nodes. More...
 
DomainDiscretization< vec_t > discretizeWithStep (scalar_t step) const
 discretizeWithStep but with default types as assigned by the shape. More...
 
virtual DomainDiscretization< vec_t > discretizeWithDensity (const std::function< scalar_t(vec_t)> &dr, int internal_type, int boundary_type) const
 Returns a discretization of the domain with spatially variable step. More...
 
template<typename func_t , typename fill_t >
DomainDiscretization< vec_t > discretizeWithDensity (const func_t &dr, const fill_t &fill, int internal_type, int boundary_type) const
 Overload for fill engine. More...
 
DomainDiscretization< vec_t > discretizeWithDensity (const std::function< scalar_t(vec_t)> &dr) const
 Overload with default types. More...
 
template<typename func_t , typename fill_t >
DomainDiscretization< vec_t > discretizeWithDensity (const func_t &dr, const fill_t &fill) const
 Overload for fill engine with default types. More...
 
DomainDiscretization< vec_t > discretizeBoundaryWithDensity (const std::function< scalar_t(vec_t)> &dr) const
 Overload with default type. More...
 
TranslatedShape< vec_t > translate (const vec_t &a)
 Translate the shape by given vector a. More...
 
RotatedShape< vec_t > rotate (const Eigen::Matrix< scalar_t, dim, dim > &Q)
 Transform the shape by given orthogonal matrix Q. More...
 
RotatedShape< vec_t > rotate (scalar_t angle)
 2D version of rotate accepting an angle. More...
 

Public Types

enum  { dim = vec_t::dim }
 Store dimension of the domain. More...
 
typedef vec_t vector_t
 Vector data type used in computations. More...
 
typedef vec_t::Scalar scalar_t
 Scalar data type used in computation. More...
 

Protected Attributes

scalar_t margin_
 Tolerance for the geometric operation of the domain. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<typename vec_t >
anonymous enum
inherited

Store dimension of the domain.

Enumerator
dim 

Dimensionality of the domain.

Definition at line 57 of file DomainShape_fwd.hpp.

Member Function Documentation

◆ add()

template<typename vec_t >
ShapeUnion< vec_t > mm::DomainShape< vec_t >::add ( const DomainShape< vec_t > &  other) const
inherited

Returns a shape representing a union of *this and other.

Definition at line 66 of file DomainShape.hpp.

◆ bbox()

template<typename vec_t >
std::pair<vec_t, vec_t> mm::UnknownShape< vec_t >::bbox ( ) const
inlineoverridevirtual

Return the bounding box of the domain.

Bounding box is returned in format bbox() == {{mx, my, ...}, {MX, MY, ...}}, such that mx <= Mx and my <= My etc.\ and that the whole domain is contained in the cuboid [mx, my, ...] x [Mx, My, ...].

Implements mm::DomainShape< vec_t >.

Definition at line 36 of file UnknownShape_fwd.hpp.

◆ clone()

template<typename vec_t >
UnknownShape* mm::UnknownShape< vec_t >::clone ( ) const
inlineoverridevirtual

Polymorphic clone pattern.

Implements mm::DomainShape< vec_t >.

Definition at line 51 of file UnknownShape_fwd.hpp.

◆ contains()

template<typename vec_t >
bool mm::UnknownShape< vec_t >::contains ( const vec_t &  point) const
inlineoverridevirtual

Return true if point is not more than margin() outside the domain.

Implements mm::DomainShape< vec_t >.

Definition at line 32 of file UnknownShape_fwd.hpp.

◆ discretizeBoundaryWithDensity() [1/2]

template<typename vec_t >
DomainDiscretization<vec_t> mm::UnknownShape< vec_t >::discretizeBoundaryWithDensity ( const std::function< scalar_t(vec_t)> &  dr,
int  type 
) const
inlineoverridevirtual

Discretizes boundary with given density and fill engine.

If type is 0, the underlying shape provides the default.

Implements mm::DomainShape< vec_t >.

Definition at line 41 of file UnknownShape_fwd.hpp.

◆ discretizeBoundaryWithDensity() [2/2]

template<typename vec_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeBoundaryWithDensity ( const std::function< scalar_t(vec_t)> &  dr) const
inlineinherited

Overload with default type.

Definition at line 189 of file DomainShape_fwd.hpp.

◆ discretizeBoundaryWithStep() [1/2]

template<typename vec_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeBoundaryWithStep ( scalar_t  step) const
inlineinherited

Returns a discretization of the boundary of this shape with approximately uniform distance step between nodes.

Node type are decided by the underlying shape.

Definition at line 126 of file DomainShape_fwd.hpp.

◆ discretizeBoundaryWithStep() [2/2]

template<typename vec_t >
virtual DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeBoundaryWithStep ( scalar_t  step,
int  type 
) const
inlinevirtualinherited

Returns a discretization of the boundary of this shape with approximately uniform distance step between nodes.

step must be positive. Added nodes are of type type, which must be non-positive. Value 0 indicates that types are dependant on the implementation of concrete shape.

Reimplemented in mm::PolygonShape< vec_t >, mm::TranslatedShape< vec_t >, mm::ShapeUnion< vec_t >, mm::BoxShape< vec_t >, mm::ShapeDifference< vec_t >, and mm::BallShape< vec_t >.

Definition at line 118 of file DomainShape_fwd.hpp.

◆ discretizeWithDensity() [1/4]

template<typename vec_t >
template<typename func_t , typename fill_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeWithDensity ( const func_t &  dr,
const fill_t &  fill 
) const
inlineinherited

Overload for fill engine with default types.

Definition at line 177 of file DomainShape_fwd.hpp.

◆ discretizeWithDensity() [2/4]

template<typename vec_t >
template<typename func_t , typename fill_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeWithDensity ( const func_t &  dr,
const fill_t &  fill,
int  internal_type,
int  boundary_type 
) const
inlineinherited

Overload for fill engine.

Definition at line 162 of file DomainShape_fwd.hpp.

◆ discretizeWithDensity() [3/4]

template<typename vec_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeWithDensity ( const std::function< scalar_t(vec_t)> &  dr) const
inlineinherited

Overload with default types.

Definition at line 170 of file DomainShape_fwd.hpp.

◆ discretizeWithDensity() [4/4]

template<typename vec_t >
DomainDiscretization< vec_t > mm::DomainShape< vec_t >::discretizeWithDensity ( const std::function< scalar_t(vec_t)> &  dr,
int  internal_type,
int  boundary_type 
) const
virtualinherited

Returns a discretization of the domain with spatially variable step.

Parameters
drFunction giving desired internodal distance at each point.
internal_typeUser supplied type of internal nodes. Must be non-negative.
boundary_typeUser supplied type of boundary nodes. Must be non-positive. If any of the types is 0, the underlying shape provides the default.
Returns
Discretization with nodes distributed according to dr.

Reimplemented in mm::TranslatedShape< vec_t >, mm::ShapeUnion< vec_t >, and mm::ShapeDifference< vec_t >.

Definition at line 77 of file DomainShape.hpp.

◆ discretizeWithStep() [1/2]

template<typename vec_t >
DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeWithStep ( scalar_t  step) const
inlineinherited

discretizeWithStep but with default types as assigned by the shape.

Definition at line 144 of file DomainShape_fwd.hpp.

◆ discretizeWithStep() [2/2]

template<typename vec_t >
virtual DomainDiscretization<vec_t> mm::DomainShape< vec_t >::discretizeWithStep ( scalar_t  step,
int  internal_type,
int  boundary_type 
) const
inlinevirtualinherited

Returns a discretization of this shape with approximately uniform distance step between nodes.

step must be positive.

Parameters
stepDesired internodal distance.
internal_typeUser supplied type of internal nodes. Must be non-negative.
boundary_typeUser supplied type of boundary nodes. Must be non-positive. If any of the types is 0, the underlying shape provides the default.

Reimplemented in mm::TranslatedShape< vec_t >, mm::ShapeUnion< vec_t >, mm::BoxShape< vec_t >, mm::ShapeDifference< vec_t >, and mm::BallShape< vec_t >.

Definition at line 137 of file DomainShape_fwd.hpp.

◆ hasContains()

template<typename vec_t >
bool mm::UnknownShape< vec_t >::hasContains ( ) const
inlineoverridevirtual

Return true if shape has contains() method implemented.

Reimplemented from mm::DomainShape< vec_t >.

Definition at line 34 of file UnknownShape_fwd.hpp.

◆ margin()

template<typename vec_t >
scalar_t mm::DomainShape< vec_t >::margin ( ) const
inlineinherited

Returns current margin.

Definition at line 72 of file DomainShape_fwd.hpp.

◆ operator+()

template<typename vec_t >
ShapeUnion<vec_t> mm::DomainShape< vec_t >::operator+ ( const DomainShape< vec_t > &  other) const
inlineinherited

Operator form of DomainShape::add.

See also
add

Definition at line 81 of file DomainShape_fwd.hpp.

◆ operator-()

template<typename vec_t >
ShapeDifference<vec_t> mm::DomainShape< vec_t >::operator- ( const DomainShape< vec_t > &  other) const
inlineinherited

Operator form of DomainShape::subtract.

See also
subtract

Definition at line 86 of file DomainShape_fwd.hpp.

◆ print()

template<typename vec_t >
std::ostream& mm::UnknownShape< vec_t >::print ( std::ostream &  os) const
inlineoverridevirtual

Output information about this shape to given output stream os.

Implements mm::DomainShape< vec_t >.

Definition at line 47 of file UnknownShape_fwd.hpp.

◆ projectPointToBoundary()

template<typename vec_t >
std::pair< bool, vec_t > mm::DomainShape< vec_t >::projectPointToBoundary ( const vec_t &  point,
const vec_t &  unit_normal 
) const
virtualinherited

Project point to boundary using bisection along the line define by unit_normal.

Definition at line 21 of file DomainShape.hpp.

◆ rotate() [1/2]

template<typename vec_t >
RotatedShape< vec_t > mm::DomainShape< vec_t >::rotate ( const Eigen::Matrix< scalar_t, dim, dim > &  Q)
inherited

Transform the shape by given orthogonal matrix Q.

Note
It is usually faster to first discretize the domain and than transform the whole discretization using DomainDiscretization::rotate.

Definition at line 98 of file DomainShape.hpp.

◆ rotate() [2/2]

template<typename vec_t >
RotatedShape< vec_t > mm::DomainShape< vec_t >::rotate ( scalar_t  angle)
inherited

2D version of rotate accepting an angle.

Definition at line 89 of file DomainShape.hpp.

◆ setMargin()

template<typename vec_t >
virtual void mm::DomainShape< vec_t >::setMargin ( scalar_t  margin)
inlinevirtualinherited

Sets domain margin to margin.

Reimplemented in mm::PolygonShape< vec_t >.

Definition at line 75 of file DomainShape_fwd.hpp.

◆ subtract()

template<typename vec_t >
ShapeDifference< vec_t > mm::DomainShape< vec_t >::subtract ( const DomainShape< vec_t > &  other) const
inherited

Returns a shape representing a difference of *this and other.

Definition at line 71 of file DomainShape.hpp.

◆ toggleMargin()

template<typename vec_t >
void mm::DomainShape< vec_t >::toggleMargin ( )
inlineinherited

Toggles the margin from positive to negative.

Definition at line 77 of file DomainShape_fwd.hpp.

◆ translate()

template<typename vec_t >
TranslatedShape< vec_t > mm::DomainShape< vec_t >::translate ( const vec_t &  a)
inherited

Translate the shape by given vector a.

Note
It is usually faster to first discretize the domain and than translate the whole discretization using DomainDiscretization::translate.

Definition at line 84 of file DomainShape.hpp.

Member Data Documentation

◆ margin_

template<typename vec_t >
scalar_t mm::DomainShape< vec_t >::margin_
protectedinherited

Tolerance for the geometric operation of the domain.

The domain should behave as if it was margin_ thicker. Default margin is 1e-10.

Definition at line 64 of file DomainShape_fwd.hpp.


The documentation for this class was generated from the following file:
mm::HDF::READONLY
@ READONLY
Read only open, the file must exist.
Definition: HDF_fwd.hpp:87
mm::DomainDiscretization::load
static DomainDiscretization< vec_t > load(hdf5_file &file, const std::string &name)
Load discretization from a HD5 file.
Definition: DomainDiscretization.hpp:354