17#include <gsUnstructuredSplines/src/gsC1SurfGluingDataAssembler.h>
25template<
class T,
class Visitor = gsC1SurfGluingDataVisitor<T>>
26class gsC1SurfGluingData :
public gsC1SurfGD<T>
37 : gsC1SurfGD<T>(mp, mb)
54 return sol.row(0) * ( ones - points ) + sol.row(1) * points;
61 return sol.row(2) * ( ones - points ) + sol.row(3) * points;
68 return solBeta.row(0) * ( ones - points ) + solBeta.row(1) * points;
75 return solBeta.row(2) * ( ones - points ) + solBeta.row(3) * points;
82 return sol.row(4) * ( ones - points ).cwiseProduct( ones - points)
83 + sol.row(5) * ( ones - points ).cwiseProduct(points) + sol.row(6) * points;
141 mSys.
matrix().makeCompressed();
152 applyBeta(visitorBeta);
154 mSysBeta.
matrix().makeCompressed();
157 void apply(Visitor visitor)
165 const int tid = omp_get_thread_num();
166 const int nt = omp_get_num_threads();
174 const gsBasis<T> & basis = this->m_mb[0].basis(0).component(1);
177 visitor_.initialize(basis,quRule);
180 typename gsBasis<T>::domainIter domIt = basis.makeDomainIterator(boundary::none);
183 for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
185 for (; domIt->good(); domIt->next() )
189 quRule.
mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
192 visitor_.evaluate(quNodes, this->m_mp);
195 visitor_.assemble(*domIt, quWeights);
198 #pragma omp critical(localToGlobal)
199 visitor_.localToGlobal( dirichletDofs, mSys);
204 void applyBeta(Visitor visitorBeta)
211 visitor_Beta(visitorBeta);
212 const int tid = omp_get_thread_num();
213 const int nt = omp_get_num_threads();
215 &visitor_Beta = visitorBeta;
221 const gsBasis<T> & basis = this->m_mb[0].basis(0).component(1);
224 visitor_Beta.initialize(basis,quRule);
227 typename gsBasis<T>::domainIter domIt = basis.makeDomainIterator(boundary::none);
230 for ( domIt->next(tid); domIt->good(); domIt->next(nt) )
232 for (; domIt->good(); domIt->next() )
236 quRule.
mapTo( domIt->lowerCorner(), domIt->upperCorner(), quNodes, quWeights );
239 visitor_Beta.evaluateBeta(quNodes, this->m_mp, sol);
242 visitor_Beta.assembleBeta(*domIt, quWeights);
245 #pragma omp critical(localToGlobal)
246 visitor_Beta.localToGlobalBeta( dirichletDofsBeta, mSysBeta);
255 typename gsSparseSolver<T>::SimplicialLDLT solver;
257 solver.compute(mSys.
matrix());
258 sol = solver.solve(mSys.
rhs());
263 typename gsSparseSolver<T>::SimplicialLDLT solver;
265 solver.compute(mSysBeta.
matrix());
266 solBeta = solver.solve(mSysBeta.
rhs());
282 solBetaTMP(0, 0) = 0;
283 solBetaTMP(1, 0) = 0;
284 solBetaTMP(2, 0) = 0;
285 solBetaTMP(3, 0) = 0;
288 solBeta = solBetaTMP;
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 boundarySize() const
Returns the number of eliminated dofs.
Definition gsDofMapper.h:457
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
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
const gsMatrix< T > & rhs() const
Access the right hand side.
Definition gsSparseSystem.h:402
const gsDofMapper & colMapper(const index_t c) const
colMapper returns the dofMapper for column block c
Definition gsSparseSystem.h:498
void reserve(const index_t nz, const index_t numRhs)
reserve reserves the memory for the sparse matrix and the rhs.
Definition gsSparseSystem.h:327
const gsSparseMatrix< T > & matrix() const
Access the system Matrix.
Definition gsSparseSystem.h:394
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Compute the gluing data for one interface.
Visitor for the Gluing Data.
The G+Smo namespace, containing all definitions for the library.