G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsGeoUtils.h
Go to the documentation of this file.
1
15#pragma once
16
17#include <gsCore/gsMultiPatch.h>
18
19namespace gismo
20{
21
22class gsParaviewCollection;
23
24//-----------------------------------//
25//--------- Mesh Analysis -----------//
26//-----------------------------------//
27
29template <class T>
30void plotGeometry(gsMultiPatch<T> const & domain, std::string fileName, index_t numSamples);
31
33template <class T>
34void plotGeometry(const gsMultiPatch<T> & domain, std::string const & fileName,
35 gsParaviewCollection & collection, index_t step);
36
39template <class T>
40void plotDeformation(const gsMultiPatch<T> & initDomain, const std::vector<gsMultiPatch<T> > & displacements,
41 std::string fileName, index_t numSamplingPoints = 10000);
42
44template <class T>
45void plotDeformation(const gsMultiPatch<T> & initDomain, const gsMultiPatch<T> & displacement,
46 std::string const & fileName, gsParaviewCollection & collection, index_t step);
47
51template <class T>
52index_t checkGeometry(gsMultiPatch<T> const & domain);
53
57template <class T>
58index_t checkDisplacement(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & displacement);
59
61template <class T>
62T normL2(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & solution);
63
66template <class T>
67T geometryJacRatio(gsMultiPatch<T> const & domain);
68
71template <class T>
72T displacementJacRatio(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & displacement);
73
76template <class T>
77void genSamplingPoints(const gsVector<T> & lower, const gsVector<T> & upper,
78 const gsQuadRule<T> & quRule, gsMatrix<T> & points);
79
81template <class T>
82T patchLength(const gsGeometry<T> & geo, short_t dir = 0);
83
85template <class T>
86T curveLength(const gsGeometry<T> & geo);
87
89template <class T>
90gsVector<unsigned> distributePoints(const gsGeometry<T> & geo, unsigned numPoints);
91
92//-----------------------------------//
93//----------- Modelling--------------//
94//-----------------------------------//
95
100template<class T>
101typename gsGeometry<T>::uPtr simplifyCurve(gsGeometry<T> const & curve,
102 index_t additionalPoints = 0, index_t degree = 0,
103 index_t numSamples = 1000);
104
106template<class T>
107T curveDistance(gsGeometry<T> const & curveA,
108 gsGeometry<T> const & curveB,
109 index_t numSamples = 1000);
110
113template <class T>
114typename gsGeometry<T>::uPtr fittingDirichlet(gsMatrix<T> const & params,
115 gsMatrix<T> const & points,
116 gsBasis<T> const & basis);
117
122template<class T>
123typename gsGeometry<T>::uPtr genPatchInterpolation(gsGeometry<T> const & A, gsGeometry<T> const & B,
124 index_t deg, index_t num, bool xiDir = false);
125
130template<class T>
131typename gsGeometry<T>::uPtr genPatchScaling(gsGeometry<T> const & boundary,
132 index_t deg, index_t num,
133 T scaling, gsVector<T> const & center);
134
137template<class T>
139 gsMatrix<T> const & A, gsMatrix<T> const & B,
140 index_t iA = 0, index_t iB = 0);
141
144template<class T>
146 T radius = 1., T x0 = 0., T y0 = 0.,
147 T angle0 = 0., T arcAngle = 2*EIGEN_PI);
148
151template<class T>
152typename gsGeometry<T>::uPtr genCircle(gsBasis<T> & basis,
153 T radius = 1., T x0 = 0., T y0 = 0.,
154 T angle0 = 0., T arcAngle = 2*EIGEN_PI);
155
159template<class T>
160typename gsGeometry<T>::uPtr genQuad(index_t xiDeg, index_t xiNum, index_t etaDeg, index_t etaNum,
161 gsMatrix<T> const & A, gsMatrix<T> const & B,
162 gsMatrix<T> const & C, gsMatrix<T> const & D,
163 index_t iA = 0, index_t iB = 0, index_t iC = 0, index_t iD = 0);
164
167template<class T>
168typename gsGeometry<T>::uPtr genSphere(index_t xiDeg, index_t xiNum, index_t etaDeg, index_t etaNum,
169 T xi0 = 0., T xi1 = 2*EIGEN_PI,
170 T eta0 = -EIGEN_PI/2, T eta1 = EIGEN_PI/2);
171
174template<class T>
175typename gsGeometry<T>::uPtr genSphere(gsKnotVector<T> & xiKnots, gsKnotVector<T> & etaKnots,
176 T xi0 = 0., T xi1 = 2*EIGEN_PI,
177 T eta0 = -EIGEN_PI/2, T eta1 = EIGEN_PI/2);
178
180template<class T>
181typename gsGeometry<T>::uPtr genCylinder(gsGeometry<T> const & base,
182 index_t deg, index_t num, T height);
183
185template<class T>
186typename gsGeometry<T>::uPtr genScrew(gsGeometry<T> const & base,
187 index_t deg, index_t num,
188 T height, T pitch, T x0 = 0., T y0 = 0.);
189
191template<class T>
192typename gsGeometry<T>::uPtr genSpring(gsGeometry<T> const & crossSection,
193 T springRadius = 6.0, T springPitch = 2.60258,
194 index_t numQuarterSegments = 12, bool nurbs = false);
195
197template<class T>
198void genMuscleMP(gsGeometry<T> const & muscleSurface, gsMultiPatch<T> & result);
199
200//----------------------------------------//
201//----------- Auxiliary functions --------//
202//----------------------------------------//
203
205template<class T>
206inline T combine(T a, T b, T x) { return a*(1-x)+b*x; }
207
211template<class T>
212gsMatrix<T> combine(gsMatrix<T> const & A, gsMatrix<T> const & B, T x,
213 index_t iA = 0, index_t iB = 0, bool cols = false);
214
219template <class T>
220T distance(gsMatrix<T> const & A, gsMatrix<T> const & B,
221 index_t i = 0, index_t j = 0, bool cols = false);
222
223} // namespace ends
224
memory::unique_ptr< gsGeometry > uPtr
Unique pointer for gsGeometry.
Definition gsGeometry.h:100
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
Provides declaration of the MultiPatch class.
The G+Smo namespace, containing all definitions for the library.
T geometryJacRatio(gsMultiPatch< T > const &domain)
Returns min(Jacobian determinant) divided by max(Jacobian determinant); samples the Jacobian elementw...
Definition gsGeoUtils.hpp:336
gsGeometry< T >::uPtr genSphere(index_t xiDeg, index_t xiNum, index_t etaDeg, index_t etaNum, T xi0=0., T xi1=2 *EIGEN_PI, T eta0=-EIGEN_PI/2, T eta1=EIGEN_PI/2)
generates a tensor product B-spline spherical patch with radius 1 and center at 0 given the degrees a...
Definition gsGeoUtils.hpp:855
void genSamplingPoints(const gsVector< T > &lower, const gsVector< T > &upper, const gsQuadRule< T > &quRule, gsMatrix< T > &points)
Generates a matrix of sampling points for a given parametric element; includes quadrature points for ...
Definition gsGeoUtils.hpp:422
T normL2(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &solution)
@ Compute norm of the isogeometric solution
Definition gsGeoUtils.hpp:295
gsVector< unsigned > distributePoints(const gsGeometry< T > &geo, unsigned numPoints)
distributes sampling points according to the length of the patch in each parametric direction
Definition gsGeoUtils.hpp:536
gsGeometry< T >::uPtr genLine(index_t deg, index_t num, gsMatrix< T > const &A, gsMatrix< T > const &B, index_t iA=0, index_t iB=0)
Definition gsGeoUtils.hpp:788
gsGeometry< T >::uPtr fittingDirichlet(gsMatrix< T > const &params, gsMatrix< T > const &points, gsBasis< T > const &basis)
fits a given parametrized point cloud with a curve using a given basis; the resulting curve interpola...
Definition gsGeoUtils.hpp:623
gsGeometry< T >::uPtr simplifyCurve(gsGeometry< T > const &curve, index_t additionalPoints=0, index_t degree=0, index_t numSamples=1000)
generates a simplified curve by fitting with the coarsest basis of the same degree; then reparametriz...
Definition gsGeoUtils.hpp:562
T displacementJacRatio(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &displacement)
Returns min(Jacobian determinant) divided by max(Jacobian determinant) for geo+disp samples the Jacob...
Definition gsGeoUtils.hpp:377
T curveLength(const gsGeometry< T > &geo)
compute curve length
Definition gsGeoUtils.hpp:510
index_t checkDisplacement(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &displacement)
Checks whether the deformed configuration is bijective, i.e. det(Jac(geo+disp)) > 0; returns -1 if ye...
Definition gsGeoUtils.hpp:245
gsGeometry< T >::uPtr genCylinder(gsGeometry< T > const &base, index_t deg, index_t num, T height)
generates a 3D tensor product B-spline cylindrical patch
Definition gsGeoUtils.hpp:897
gsGeometry< T >::uPtr genSpring(gsGeometry< T > const &crossSection, T springRadius=6.0, T springPitch=2.60258, index_t numQuarterSegments=12, bool nurbs=false)
generates a 3D NURBS spring using provided geometry as a cross-section
Definition gsGeoUtils.hpp:985
T combine(T a, T b, T x)
compute a convex combintation (1-x)a+xb
Definition gsGeoUtils.h:206
T curveDistance(gsGeometry< T > const &curveA, gsGeometry< T > const &curveB, index_t numSamples=1000)
returns a distance in L2 sense between two curves parametrized from 0 to 1
Definition gsGeoUtils.hpp:603
gsGeometry< T >::uPtr genCircle(index_t deg, index_t num, T radius=1., T x0=0., T y0=0., T angle0=0., T arcAngle=2 *EIGEN_PI)
Definition gsGeoUtils.hpp:813
gsGeometry< T >::uPtr genPatchScaling(gsGeometry< T > const &boundary, index_t deg, index_t num, T scaling, gsVector< T > const &center)
generates a tensor product B-spline bdry south | front patch by scaling a given geometry object \ / |...
Definition gsGeoUtils.hpp:776
T patchLength(const gsGeometry< T > &geo, short_t dir=0)
compute length of a patch in a given parametric direction as a mean of all boundary edges correspondi...
Definition gsGeoUtils.hpp:439
gsGeometry< T >::uPtr genScrew(gsGeometry< T > const &base, index_t deg, index_t num, T height, T pitch, T x0=0., T y0=0.)
generates a 3D tensor product B-spline screw-like patch
Definition gsGeoUtils.hpp:915
T distance(gsMatrix< T > const &A, gsMatrix< T > const &B, index_t i=0, index_t j=0, bool cols=false)
compute a distance between the point number in the set and the point number <j> in the set ; by def...
void plotDeformation(const gsMultiPatch< T > &initDomain, const std::vector< gsMultiPatch< T > > &displacements, std::string fileName, index_t numSamplingPoints=10000)
Definition gsGeoUtils.hpp:95
gsGeometry< T >::uPtr genQuad(index_t xiDeg, index_t xiNum, index_t etaDeg, index_t etaNum, gsMatrix< T > const &A, gsMatrix< T > const &B, gsMatrix< T > const &C, gsMatrix< T > const &D, index_t iA=0, index_t iB=0, index_t iC=0, index_t iD=0)
generates a quad patch given by its four C—D corners with the following orientation; | | the points a...
Definition gsGeoUtils.hpp:843
gsGeometry< T >::uPtr genPatchInterpolation(gsGeometry< T > const &A, gsGeometry< T > const &B, index_t deg, index_t num, bool xiDir=false)
Definition gsGeoUtils.hpp:687
index_t checkGeometry(gsMultiPatch< T > const &domain)
Checks whether configuration is bijective, i.e. det(Jac(geo)) > 0; returns -1 if yes or the number of...
Definition gsGeoUtils.hpp:200
void genMuscleMP(gsGeometry< T > const &muscleSurface, gsMultiPatch< T > &result)
This is more of a script than a function. I use it to generate a multi-parametrization for the biceps...
Definition gsGeoUtils.hpp:1038
void plotGeometry(gsMultiPatch< T > const &domain, std::string fileName, index_t numSamples)
Plots the mesh and the jacobian (if <numSamples> > 0) to Paraview.
Definition gsGeoUtils.hpp:43