G+Smo
24.08.0
Geometry + Simulation Modules
|
Provides isogeometric meshing and modelling routines. More...
Go to the source code of this file.
Namespaces | |
gismo | |
The G+Smo namespace, containing all definitions for the library. | |
Functions | |
template<class T > | |
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 yes or the number of the first invalid patch; samples the Jacobian elementwise at the quadrature points and the corners. | |
template<class T > | |
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 the first invalid patch; samples the Jacobian elementwise at the quadrature points and the corners. | |
template<class T > | |
T | combine (T a, T b, T x) |
compute a convex combintation (1-x)a+xb | |
template<class T > | |
gsMatrix< T > | combine (gsMatrix< T > const &A, gsMatrix< T > const &B, T x, index_t iA=0, index_t iB=0, bool cols=false) |
compute a convex combination of two points given as ROWS of matrices numbered <iA> and <iB>; set <cols> to <true> to give points as COLUMNS | |
template<class T > | |
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 | |
template<class T > | |
T | curveLength (const gsGeometry< T > &geo) |
compute curve length | |
template<class T > | |
T | displacementJacRatio (gsMultiPatch< T > const &domain, gsMultiPatch< T > const &displacement) |
Returns min(Jacobian determinant) divided by max(Jacobian determinant) for geo+disp samples the Jacobian elementwise at the quadrature points and the corners. | |
template<class T > | |
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 default, points are given as ROWS of matrices; set <cols> to <true> to give points as COLUMNS | |
template<class T > | |
gsVector< unsigned > | distributePoints (const gsGeometry< T > &geo, unsigned numPoints) |
distributes sampling points according to the length of the patch in each parametric direction | |
template<class T > | |
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 interpolates the first and the last points | |
template<class T > | |
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) |
template<class T > | |
gsGeometry< T >::uPtr | genCircle (gsBasis< T > &basis, T radius=1., T x0=0., T y0=0., T angle0=0., T arcAngle=2 *EIGEN_PI) |
template<class T > | |
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 | |
template<class T > | |
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) |
template<class T > | |
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 model given its surface. | |
template<class T > | |
gsGeometry< T >::uPtr | genPatchInterpolation (gsGeometry< T > const &A, gsGeometry< T > const &B, index_t deg, index_t num, bool xiDir=false) |
template<class T > | |
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 \ / | | | towards a given center point; (x,y) north | back oppositely lying bdry is generated by scaling the original boundary with <scaling> coeff | |
template<class T > | |
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 are given as ROWS of matrices A—B | |
template<class T > | |
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 the element as well as the corner points. | |
template<class T > | |
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 | |
template<class T > | |
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 and number of control points in two parametric dimensions | |
template<class T > | |
gsGeometry< T >::uPtr | genSphere (gsKnotVector< T > &xiKnots, gsKnotVector< T > &etaKnots, 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 knot vectors for two parametric dimensions | |
template<class T > | |
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 | |
template<class T > | |
T | geometryJacRatio (gsMultiPatch< T > const &domain) |
Returns min(Jacobian determinant) divided by max(Jacobian determinant); samples the Jacobian elementwise at the quadrature points and the corners. | |
template<class T > | |
T | normL2 (gsMultiPatch< T > const &domain, gsMultiPatch< T > const &solution) |
@ Compute norm of the isogeometric solution | |
template<class T > | |
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 corresponding to this direction | |
template<class T > | |
void | plotDeformation (const gsMultiPatch< T > &initDomain, const std::vector< gsMultiPatch< T > > &displacements, std::string fileName, index_t numSamplingPoints=10000) |
template<class T > | |
void | plotDeformation (const gsMultiPatch< T > &initDomain, const gsMultiPatch< T > &displacement, std::string const &fileName, gsParaviewCollection &collection, index_t step) |
plot a deformed isogeometric mesh and add it to a Paraview collection | |
template<class T > | |
void | plotGeometry (gsMultiPatch< T > const &domain, std::string fileName, index_t numSamples) |
Plots the mesh and the jacobian (if <numSamples> > 0) to Paraview. | |
template<class T > | |
void | plotGeometry (const gsMultiPatch< T > &domain, std::string const &fileName, gsParaviewCollection &collection, index_t step) |
plot an isogeometric mesh and add to collection | |
template<class T > | |
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 reparametrizes it using the basis of the original curve; <additionalPoints> increases the number of degrees of freedom; <numSamples> is a number of sampling points for reparametrization | |
Provides isogeometric meshing and modelling routines.
This file is part of the G+Smo library.
This Source Code Form is subject to the terms of the Mozilla Public License, v. 2.0. If a copy of the MPL was not distributed with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
Author(s): A.Shamanskiy (2016 - ...., TU Kaiserslautern)