22 class gsParaviewCollection;
30 void plotGeometry(gsMultiPatch<T>
const & domain, std::string fileName,
index_t numSamples);
34 void plotGeometry(
const gsMultiPatch<T> & domain, std::string
const & fileName,
35 gsParaviewCollection & collection,
index_t step);
40 void plotDeformation(
const gsMultiPatch<T> & initDomain,
const std::vector<gsMultiPatch<T> > & displacements,
41 std::string fileName,
index_t numSamplingPoints = 10000);
45 void plotDeformation(
const gsMultiPatch<T> & initDomain,
const gsMultiPatch<T> & displacement,
46 std::string
const & fileName, gsParaviewCollection & collection,
index_t step);
62 T
normL2(gsMultiPatch<T>
const & domain, gsMultiPatch<T>
const & solution);
78 const gsQuadRule<T> & quRule, gsMatrix<T> & points);
90 gsVector<unsigned>
distributePoints(
const gsGeometry<T> & geo,
unsigned numPoints);
101 typename gsGeometry<T>::uPtr
simplifyCurve(gsGeometry<T>
const & curve,
108 gsGeometry<T>
const & curveB,
115 gsMatrix<T>
const & points,
116 gsBasis<T>
const & basis);
131 typename gsGeometry<T>::uPtr
genPatchScaling(gsGeometry<T>
const & boundary,
133 T scaling, gsVector<T>
const & center);
139 gsMatrix<T>
const & A, gsMatrix<T>
const & B,
146 T radius = 1., T x0 = 0., T y0 = 0.,
147 T angle0 = 0., T arcAngle = 2*EIGEN_PI);
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);
161 gsMatrix<T>
const & A, gsMatrix<T>
const & B,
162 gsMatrix<T>
const & C, gsMatrix<T>
const & D,
169 T xi0 = 0., T xi1 = 2*EIGEN_PI,
170 T eta0 = -EIGEN_PI/2, T eta1 = EIGEN_PI/2);
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);
181 typename gsGeometry<T>::uPtr
genCylinder(gsGeometry<T>
const & base,
186 typename gsGeometry<T>::uPtr
genScrew(gsGeometry<T>
const & base,
188 T height, T pitch, T x0 = 0., T y0 = 0.);
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);
198 void genMuscleMP(gsGeometry<T>
const & muscleSurface, gsMultiPatch<T> & result);
206 inline T
combine(T a, T b, T x) {
return a*(1-x)+b*x; }
212 gsMatrix<T>
combine(gsMatrix<T>
const & A, gsMatrix<T>
const & B, T x,
220 T
distance(gsMatrix<T>
const & A, gsMatrix<T>
const & B,
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 <j> 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 <numSamples> > 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 ¢er)
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)) > 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 ¶ms, 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)) > 0; returns -1 if yes or the number of...
Definition: gsGeoUtils.hpp:200