G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsMPBESBSplineBasis.h
Go to the documentation of this file.
1 
14 #pragma once
15 
18 #include <gsNurbs/gsKnotVector.h>
19 
20 namespace gismo
21 {
22 
31 template<short_t d, class T>
32 class gsMPBESBSplineBasis : public gsMPBESBasis<d,T>
33 {
34 public:
36  typedef memory::shared_ptr< gsMPBESBSplineBasis > Ptr;
37 
39  typedef memory::unique_ptr< gsMPBESBSplineBasis > uPtr;
40 
42  static const short_t Dim = d;
43 
44  typedef index_t indexType;
45  typedef gsMPBESBasis<d,T> Base;
47  typedef typename std::vector<BasisType *> BasisContainer;
48  typedef typename std::vector<gsBasis<T>* >::const_iterator ConstBasisIter;
49  typedef typename std::vector<gsBasis<T>* >::iterator BasisIter;
50  typedef typename std::vector<gsMatrix<T> *>::const_iterator ConstMatrixPtrIter;
51 
52 protected:
53  using Base::m_bases;
54  using Base::m_topol;
56  using Base::m_mapper;
57  using Base::m_vertices;
58  using Base::m_distances;
59  using Base::m_minDist;
60 public:
61  using Base::updateTopol;
62  using Base::degree;
63  using Base::global_coef_to_local_coef;
64  using Base::local_coef_to_global_coef;
65  using Base::exportToPatches;
66  using Base::nPatches;
68  using Base::repairPatches;
69  using Base::maxDegree;
70 protected:
71  using Base::_getPatch;
72  using Base::_getPatchIndex;
73  using Base::_initVertices;
75 
76 public:
79  { }
80 
81  gsMPBESBSplineBasis (const BasisContainer & bases, gsBoxTopology const & topol,
82  index_t increaseSmoothnessLevel=-1, index_t minEVDistance=-1);
83 
84  gsMPBESBSplineBasis (const BasisContainer & bases, gsBoxTopology const & topol,std::vector<gsMatrix<T> * > coefs,
85  index_t increaseSmoothnessLevel=-1, index_t minEVDistance=-1);
86 
87  gsMPBESBSplineBasis( gsMultiPatch<T> const & mp, index_t increaseSmoothnessLevel = -1,
88  index_t minEVDistance=-1);
89 
91 
92  gsMPBESBSplineBasis<d,T> & operator=( const gsMPBESBSplineBasis& other );
93 
95  { } //destructor
96 
97 private:
99  // functions for initializing and updating
101 
103  {
104 // T eps=1e-6;
105  bool consistent = true;
106 // index_t nrOfBases = m_bases.size();
107 // consistent = consistent && nrOfBases>0;
108 // GISMO_ASSERT(nrOfBases>0,"empty list of bases");
109 // consistent = consistent && nrOfBases==m_topol.size();
110 // GISMO_ASSERT(nrOfBases==m_topol.size(),"bases and topol not of same size.");
111 // unsigned i = 0;
112 // for(ConstBasisIter it=m_bases.begin();it!=m_bases.end();++it)
113 // {
114 // index_t deg_u = basis(i).degree(0);
115 // index_t deg_v = basis(i).degree(1);
116 // unsigned u_dim = basis(i).size(0), v_dim = basis(i).size(1);
117 // consistent = consistent && u_dim >= static_cast<unsigned>(deg_u*2);
118 // consistent = consistent && v_dim >= static_cast<unsigned>(deg_v*2);
119 // if(u_dim<static_cast<unsigned>(deg_u*2))
120 // GISMO_ERROR("patch is too small: u_dim too small");
121 // if(v_dim<static_cast<unsigned>(deg_v*2))
122 // GISMO_ERROR("patch is too small: v_dim too small");
123 // i++;
124 // }
125 // std::vector<boundaryInterface> interfaces = m_topol.interfaces();
126 // for(unsigned i=0;i<interfaces.size();i++)
127 // {
128 // patchSide first() = interfaces[i].first();
129 // patchSide second() = interfaces[i].second();
130 // index_t nrOfBasisFuncs1=basis(first().patch).size(!direction(first().side));
131 // index_t nrOfBasisFuncs2=basis(second().patch).size(!direction(second().side));
132 // consistent = consistent && nrOfBasisFuncs1==nrOfBasisFuncs2;
133 // if(!(nrOfBasisFuncs1==nrOfBasisFuncs2))
134 // GISMO_ERROR("number of basisfunctions on interface are different");
135 // index_t deg_1 = m_bases[first().patch]->degree(!direction(first().side));
136 // index_t deg_2 = m_bases[first().patch]->degree(!direction(second().side));
137 // consistent = consistent && deg_1==deg_2;
138 // if(!(deg_1==deg_2))
139 // GISMO_ERROR("degrees of patches on interface are different");
140 // gsKnotVector<T> kv1 = basis(first().patch).knots(!direction(first().side));
141 // gsKnotVector<T> kv2 = basis(second().patch).knots(!direction(second().side));
142 // if(kv1.size()!=kv2.size())
143 // GISMO_ERROR("different knot vector lengths");
144 // if(!interfaces[i].orient[0])
145 // kv2.reverse();
146 // for(index_t j = 0;j<kv1.size();j++)
147 // {
148 // if(std::abs(kv1[j]-kv2[j])>eps)
149 // {
150 // GISMO_ERROR("different knot elements.");
151 // }
152 // }
153 // }
154  return consistent;
155  }
156 
157 protected:
158  void _setMapping();
159 
160 public:
162  // general functions for interacting with this class
164 
166  GISMO_CLONE_FUNCTION(gsMPBESBSplineBasis)
167 
168  gsTensorBSplineBasis<d,T> & basis(size_t i)
169  { return static_cast<gsTensorBSplineBasis<d,T>&>(*m_bases[i]); }
170 
171  const gsTensorBSplineBasis<d,T> & basis(size_t i) const
172  { return static_cast<const gsTensorBSplineBasis<d,T>&>(*m_bases[i]); }
173 
174  unsigned basisFunctionsOnSide(const patchSide& ps) const;
175 
176  virtual bool isLocallyConnected(indexType i,indexType j) const;
177 
178 public:
180  // functions for evaluating and derivatives
182 
184  void numActive_into(const index_t patch,const gsMatrix<T> & u, gsVector<index_t>& result) const;
185 
186 public:
188  // functions for refinement
190 
191  void refine(const index_t patch,const std::vector<T>& knots_u, const std::vector<T>& knots_v,bool updateBasis = true);
192 
193  void refine(const index_t patch, gsMatrix<T> const & boxes, bool updateBasis = true);
194 
195  void refineElements(const index_t patch, std::vector<index_t> const & boxes, bool updateBasis = true);
196 
197  void refine_withCoefs(gsMatrix<T>& localCoef, const index_t patch,const std::vector<T>& knots_u, const std::vector<T>& knots_v,
198  bool updateBasis = true);
199 
200  void refine_withCoefs(gsMatrix<T>& coefs, const index_t patch, gsMatrix<T> const & boxes,
201  bool updateBasis = true);
202 
203  void refineElements_withCoefs(gsMatrix<T>& coefs, const index_t patch, std::vector<index_t> const & boxes,
204  bool updateBasis = true);
205 
206 private:
208  // private helper functions for the refinement
210 
211  void _boxesVectorToMatrix(const std::vector<index_t> & boxes, gsMatrix<T> & mat_boxes);
212 
213  void _boxToRefineKnots(const index_t patch,gsMatrix<T> const & boxes,std::vector<std::vector<T> > & refineKnots);
214 
215  // takes a knotvector, a standard vector of new knots which will be inserted, an indicator par which tells
216  // if the start or the end of the knotvector has to be checked for new knots, and a distance indicating the distance
217  // from start/end where new knots are looked for. Returns the number of new knots inserted in the knotspan, going
218  // distance knots from the beginning or end of the given knotvector.
219  unsigned getNrOfSpecialKnots(const gsKnotVector<T> kv,const std::vector<T>& new_knots,bool par,index_t distance);
220 
221  void repairPatches(std::vector<gsMatrix<T> *> & coefs,
222  index_t startFromPatch = -1);
223 
224  bool _knotsMatchNeighbours(index_t patch,std::vector<std::vector<T> >& knotsToInsert,
225  std::vector<index_t>& checkPatches,
226  T eps = 1e-10);
227 
228 protected:
230  // functions going back and forth between absolute val and parametric value for C^0 parts
232 
233  T findParameter(patchSide const & ps,patchCorner const & pc,unsigned nrBasisFuncs) const;
234 
235 }; // class gsMPBESBSplineBasis
236 
237 }
238 
A univariate Lagrange basis.
Definition: gsMPBESBSplineBasis.h:32
Knot vector for B-splines.
bool isSpecialVertex(const patchCorner &pc) const
gives back true, if the given patchCorner is a special vertex
Definition: gsMPBESBasis.hpp:362
Provides declaration of Basis abstract interface.
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...
#define short_t
Definition: gsConfig.h:35
memory::shared_ptr< gsMPBESBSplineBasis > Ptr
Shared pointer for gsMPBESBSplineBasis.
Definition: gsMPBESBSplineBasis.h:36
void repairPatches(std::vector< gsMatrix< T > * > &coefs, index_t startFromPatch=-1)
Definition: gsMPBESBSplineBasis.hpp:338
void _setMapping()
create a new mapping of the local basisfunctions
Definition: gsMPBESBSplineBasis.hpp:132
short_t m_incrSmoothnessDegree
smoothness degree that is tried to achive over patch interfaces
Definition: gsMPBESBasis.h:375
#define index_t
Definition: gsConfig.h:32
memory::unique_ptr< gsMPBESBSplineBasis > uPtr
Unique pointer for gsMPBESBSplineBasis.
Definition: gsMPBESBSplineBasis.h:39
A tensor product B-spline basis.
Definition: gsTensorBSplineBasis.h:36
std::vector< std::pair< patchCorner, bool > > m_vertices
vector storing all the inner vertices, with a flag if it is an Extraordinary vertex or an ordinary ve...
Definition: gsMPBESBasis.h:381
void refineElements(const index_t patch, std::vector< index_t > const &boxes, bool updateBasis=true)
Definition: gsMPBESBSplineBasis.hpp:202
T findParameter(patchSide const &ps, patchCorner const &pc, unsigned nrBasisFuncs) const
Definition: gsMPBESBSplineBasis.hpp:436
virtual void repairPatches(std::vector< gsMatrix< T > * > &coefs, index_t startFromPatch=-1)=0
void updateTopol()
updates the mapping of this basis (f.e. after a knot insertion)
Definition: gsMPBESBasis.h:101
static const short_t Dim
Dimension of the parameter domain.
Definition: gsMPBESBSplineBasis.h:42
gsTensorBSplineBasis< d, T > & basis(size_t i)
Clone function. Used to make a copy of a derived basis.
Definition: gsMPBESBSplineBasis.h:168
Purely abstract class gsMappedBasis, which gives means of combining basis functions to new...
Definition: gsMPBESBasis.h:46
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
void numActive_into(const index_t patch, const gsMatrix< T > &u, gsVector< index_t > &result) const
Returns the number of active (nonzero) basis functions at points u in result.
Definition: gsMPBESBSplineBasis.hpp:170
unsigned basisFunctionsOnSide(const patchSide &ps) const
Returns the amount of basis functions on a given side of a given patch.
Definition: gsMPBESBSplineBasis.hpp:142
gsMPBESBSplineBasis()
Default empty constructor.
Definition: gsMPBESBSplineBasis.h:78
bool _checkTopologyWithBases() const
Checks the gsMappedBasis for consistency.
Definition: gsMPBESBSplineBasis.h:102
void refineElements_withCoefs(gsMatrix< T > &coefs, const index_t patch, std::vector< index_t > const &boxes, bool updateBasis=true)
Definition: gsMPBESBSplineBasis.hpp:258
std::vector< distances > m_distances
vector of distances objects, that store C^0 distances from special vertices on the edges ...
Definition: gsMPBESBasis.h:383
void _initVertices()
initializes the m_vertices field
Definition: gsMPBESBasis.hpp:298
Defines a topological arrangement of a collection of &quot;boxes&quot; (e.g., parameter domains that map to phy...
Definition: gsBoxTopology.h:38
unsigned m_minDist
minimal C^0 distance from special (extraordinary) vertices, specified in basisfunctions ...
Definition: gsMPBESBasis.h:377
Provides declaration of Basis abstract interface.
A basis represents a family of scalar basis functions defined over a common parameter domain...
Definition: gsBasis.h:78
void _setDistanceOfAllVertices()
initializes the m_distances field
Definition: gsMPBESBasis.hpp:313