16#include <gsUnstructuredSplines/src/gsC1SurfGluingData.h>
17#include <gsUnstructuredSplines/src/gsC1SurfVisitorBasisEdge.h>
21 template<
class T,
class bhVisitor = gsC1SurfVisitorBasisEdge<T>>
32 gsC1SurfGluingData<T> gluingD)
33 : m_mp(mp), m_basis(basis), m_uv(uv), m_isBoundary(isBoundary), m_gD(gluingD)
40 basis_plus.elevateContinuity(1);
41 m_basis_plus = basis_plus;
46 basis_minus.degreeReduce(1);
47 m_basis_minus = basis_minus;
51 m_basis_g1 = m_basis.basis(0);
62 void assemble(index_t i, std::string typeBf);
63 inline void apply(bhVisitor & visitor, index_t i, std::string typeBf);
67 void constructSolution(
const gsMatrix<T> & solVector,
88 gsC1SurfGluingData<T> m_gD;
107 template <
class T,
class bhVisitor>
108 void gsC1SurfBasisEdge<T,bhVisitor>::setG1BasisEdge(
gsMultiPatch<T> & result)
112 index_t n_plus = m_basis_plus.size();
113 index_t n_minus = m_basis_minus.size();
118 for (index_t bfID = bfID_init; bfID < n_plus - bfID_init; bfID++)
124 assemble(bfID,
"plus");
126 typename gsSparseSolver<T>::SimplicialLDLT solver;
129 solver.compute(m_system.matrix());
130 sol = solver.solve(m_system.rhs());
132 constructSolution(sol,g1EdgeBasis);
135 for (index_t bfID = bfID_init; bfID < n_minus-bfID_init; bfID++)
143 assemble(bfID,
"minus");
145 typename gsSparseSolver<T>::SimplicialLDLT solver;
148 solver.compute(m_system.matrix());
149 sol = solver.solve(m_system.rhs());
151 constructSolution(sol,g1EdgeBasis);
154 result = g1EdgeBasis;
157 template <
class T,
class bhVisitor>
158 void gsC1SurfBasisEdge<T,bhVisitor>::constructSolution(
const gsMatrix<T> & solVector,
gsMultiPatch<T> & result, short_t unk)
const
163 const index_t dim = ( 0!=solVector.cols() ? solVector.cols() : m_ddof[0].cols() );
166 const gsDofMapper & mapper = m_system.colMapper(0);
170 sz = m_basis.basis(0).size();
172 coeffs.resize(sz, dim);
174 for (index_t i = 0; i < sz; ++i)
178 coeffs.row(i) = solVector.row(mapper.
index(i, 0));
183 coeffs.row(i) = m_ddof[0].row( mapper.
bindex(i, 0) ).head(dim);
186 result.
addPatch(m_basis_g1.basis(0).makeGeometry(
give(coeffs)));
190 template <
class T,
class bhVisitor>
191 void gsC1SurfBasisEdge<T,bhVisitor>::refresh()
198 for (index_t i = 2; i < m_basis.basis(0).component(1-m_uv).size(); i++)
200 act = m_basis.basis(0).boundaryOffset(m_uv == 0 ? 3 : 1, i);
201 map.markBoundary(0, act);
211 template <
class T,
class bhVisitor>
212 void gsC1SurfBasisEdge<T,bhVisitor>::assemble(index_t bfID, std::string typeBf)
215 const index_t nz = gsAssemblerOptions::numColNz(m_basis[0],2,1,0.333333);
216 m_system.reserve(nz, 1);
227 apply(visitor, bfID, typeBf);
229 m_system.matrix().makeCompressed();
233 template <
class T,
class bhVisitor>
234 void gsC1SurfBasisEdge<T,bhVisitor>::apply(bhVisitor & visitor, index_t bf_index, std::string typeBf)
247 const int tid = omp_get_thread_num();
248 const int nt = omp_get_num_threads();
261 visitor_.initialize(basis_g1, quRule);
266 typename gsBasis<T>::domainIter domIt = m_geo.basis(0).makeDomainIterator(boundary::none);
269 for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
271 for (; domIt->good(); domIt->next() )
275 quRule.
mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
278 visitor_.evaluate(bf_index, typeBf, basis_g1, basis_geo, basis_plus, basis_minus, patch, quNodes, m_uv, m_gD, m_isBoundary);
281 visitor_.assemble(*domIt, quWeights);
284#pragma omp critical(localToGlobal)
285 visitor_.localToGlobal(0, m_ddof, m_system);
The assembler class provides generic routines for volume and boundary integrals that are used for for...
Definition gsAssembler.h:266
gsSparseSystem< T > m_system
Global sparse linear system.
Definition gsAssembler.h:290
virtual void constructSolution(const gsMatrix< T > &solVector, gsMultiPatch< T > &result, short_t unk=0) const
Construct solution from computed solution vector for a single unknows.
Definition gsAssembler.hpp:537
void apply(ElementVisitor &visitor, size_t patchIndex=0, boxSide side=boundary::none)
Generic assembly routine for volume or boundary integrals.
Definition gsAssembler.h:668
std::vector< gsMatrix< T > > m_ddof
Definition gsAssembler.h:295
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition gsDofMapper.h:69
index_t bindex(index_t i, index_t k=0, index_t c=0) const
Returns the boundary index of local dof i of patch k.
Definition gsDofMapper.h:334
index_t boundarySize() const
Returns the number of eliminated dofs.
Definition gsDofMapper.h:457
index_t index(index_t i, index_t k=0, index_t c=0) const
Returns the global dof index associated to local dof i of patch k.
Definition gsDofMapper.h:325
bool is_free(index_t i, index_t k=0, index_t c=0) const
Returns true if local dof i of patch k is not eliminated.
Definition gsDofMapper.h:382
const gsBasis< T > & basis(const index_t k) const
Helper which casts and returns the k-th piece of this function set as a gsBasis.
Definition gsFunctionSet.hpp:33
Abstract base class representing a geometry map.
Definition gsGeometry.h:93
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
index_t addPatch(typename gsGeometry< T >::uPtr g)
Add a patch from a gsGeometry<T>::uPtr.
Definition gsMultiPatch.hpp:211
void clear()
Clear (delete) all patches.
Definition gsMultiPatch.h:388
Class representing a reference quadrature rule.
Definition gsQuadRule.h:29
virtual void mapTo(const gsVector< T > &lower, const gsVector< T > &upper, gsMatrix< T > &nodes, gsVector< T > &weights) const
Maps quadrature rule (i.e., points and weights) from the reference domain to an element.
Definition gsQuadRule.h:177
A sparse linear system indexed by sets of degrees of freedom.
Definition gsSparseSystem.h:30
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266