G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsGeoUtils.h
Go to the documentation of this file.
1 
15 #pragma once
16 
17 #include <gsCore/gsMultiPatch.h>
18 
19 namespace gismo
20 {
21 
22 class gsParaviewCollection;
23 
24 //-----------------------------------//
25 //--------- Mesh Analysis -----------//
26 //-----------------------------------//
27 
29 template <class T>
30 void plotGeometry(gsMultiPatch<T> const & domain, std::string fileName, index_t numSamples);
31 
33 template <class T>
34 void plotGeometry(const gsMultiPatch<T> & domain, std::string const & fileName,
35  gsParaviewCollection & collection, index_t step);
36 
39 template <class T>
40 void plotDeformation(const gsMultiPatch<T> & initDomain, const std::vector<gsMultiPatch<T> > & displacements,
41  std::string fileName, index_t numSamplingPoints = 10000);
42 
44 template <class T>
45 void plotDeformation(const gsMultiPatch<T> & initDomain, const gsMultiPatch<T> & displacement,
46  std::string const & fileName, gsParaviewCollection & collection, index_t step);
47 
51 template <class T>
52 index_t checkGeometry(gsMultiPatch<T> const & domain);
53 
57 template <class T>
58 index_t checkDisplacement(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & displacement);
59 
61 template <class T>
62 T normL2(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & solution);
63 
66 template <class T>
67 T geometryJacRatio(gsMultiPatch<T> const & domain);
68 
71 template <class T>
72 T displacementJacRatio(gsMultiPatch<T> const & domain, gsMultiPatch<T> const & displacement);
73 
76 template <class T>
77 void genSamplingPoints(const gsVector<T> & lower, const gsVector<T> & upper,
78  const gsQuadRule<T> & quRule, gsMatrix<T> & points);
79 
81 template <class T>
82 T patchLength(const gsGeometry<T> & geo, short_t dir = 0);
83 
85 template <class T>
86 T curveLength(const gsGeometry<T> & geo);
87 
89 template <class T>
90 gsVector<unsigned> distributePoints(const gsGeometry<T> & geo, unsigned numPoints);
91 
92 //-----------------------------------//
93 //----------- Modelling--------------//
94 //-----------------------------------//
95 
100 template<class T>
101 typename gsGeometry<T>::uPtr simplifyCurve(gsGeometry<T> const & curve,
102  index_t additionalPoints = 0, index_t degree = 0,
103  index_t numSamples = 1000);
104 
106 template<class T>
107 T curveDistance(gsGeometry<T> const & curveA,
108  gsGeometry<T> const & curveB,
109  index_t numSamples = 1000);
110 
113 template <class T>
114 typename gsGeometry<T>::uPtr fittingDirichlet(gsMatrix<T> const & params,
115  gsMatrix<T> const & points,
116  gsBasis<T> const & basis);
117 
122 template<class T>
123 typename gsGeometry<T>::uPtr genPatchInterpolation(gsGeometry<T> const & A, gsGeometry<T> const & B,
124  index_t deg, index_t num, bool xiDir = false);
125 
130 template<class T>
131 typename gsGeometry<T>::uPtr genPatchScaling(gsGeometry<T> const & boundary,
132  index_t deg, index_t num,
133  T scaling, gsVector<T> const & center);
134 
137 template<class T>
138 typename gsGeometry<T>::uPtr genLine(index_t deg, index_t num,
139  gsMatrix<T> const & A, gsMatrix<T> const & B,
140  index_t iA = 0, index_t iB = 0);
141 
144 template<class T>
145 typename gsGeometry<T>::uPtr genCircle(index_t deg, index_t num,
146  T radius = 1., T x0 = 0., T y0 = 0.,
147  T angle0 = 0., T arcAngle = 2*EIGEN_PI);
148 
151 template<class T>
152 typename 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 
159 template<class T>
160 typename 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 
167 template<class T>
168 typename 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 
174 template<class T>
175 typename 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 
180 template<class T>
181 typename gsGeometry<T>::uPtr genCylinder(gsGeometry<T> const & base,
182  index_t deg, index_t num, T height);
183 
185 template<class T>
186 typename 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 
191 template<class T>
192 typename 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 
197 template<class T>
198 void genMuscleMP(gsGeometry<T> const & muscleSurface, gsMultiPatch<T> & result);
199 
200 //----------------------------------------//
201 //----------- Auxiliary functions --------//
202 //----------------------------------------//
203 
205 template<class T>
206 inline T combine(T a, T b, T x) { return a*(1-x)+b*x; }
207 
211 template<class T>
212 gsMatrix<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 
219 template <class T>
220 T 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 
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
T geometryJacRatio(gsMultiPatch< T > const &domain)
Returns min(Jacobian determinant) divided by max(Jacobian determinant); samples the Jacobian elementw...
Definition: gsGeoUtils.hpp:336
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 &lt;j&gt; in the set ; by def...
T curveLength(const gsGeometry< T > &geo)
compute curve length
Definition: gsGeoUtils.hpp:510
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
#define short_t
Definition: gsConfig.h:35
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
#define index_t
Definition: gsConfig.h:32
T combine(T a, T b, T x)
compute a convex combintation (1-x)a+xb
Definition: gsGeoUtils.h:206
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 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
void plotGeometry(gsMultiPatch< T > const &domain, std::string fileName, index_t numSamples)
Plots the mesh and the jacobian (if &lt;numSamples&gt; &gt; 0) to Paraview.
Definition: gsGeoUtils.hpp:43
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
Provides declaration of the MultiPatch class.
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
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:913
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
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 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:895
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:1036
index_t checkDisplacement(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &displacement)
Checks whether the deformed configuration is bijective, i.e. det(Jac(geo+disp)) &gt; 0; returns -1 if ye...
Definition: gsGeoUtils.hpp:245
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
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:983
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 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
T normL2(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &solution)
@ Compute norm of the isogeometric solution
Definition: gsGeoUtils.hpp:295
void plotDeformation(const gsMultiPatch< T > &initDomain, const std::vector< gsMultiPatch< T > > &displacements, std::string fileName, index_t numSamplingPoints=10000)
Definition: gsGeoUtils.hpp:95
index_t checkGeometry(gsMultiPatch< T > const &domain)
Checks whether configuration is bijective, i.e. det(Jac(geo)) &gt; 0; returns -1 if yes or the number of...
Definition: gsGeoUtils.hpp:200