G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsGenericAssembler.hpp
1
14#pragma once
15
23
24namespace gismo
25{
26
27template <class T>
29{
31 options.update( gsVisitorDg<T>::defaultOptions(), gsOptionList::addIfUnknown );
32 options.update( gsVisitorNitsche<T>::defaultOptions(), gsOptionList::addIfUnknown );
33 return options;
34}
35
36template <class T>
37const gsSparseMatrix<T> & gsGenericAssembler<T>::assembleMass(const index_t patchIndex, const bool refresh)
38{
39 // Clean the sparse system
40 if (refresh)
41 {
43 m_system.reserve(m_bases[0], m_options, 1);
44 Base::computeDirichletDofs();
45 }
46
47 // Assemble mass integrals
48 if ( patchIndex < 0 )
49 {
50 this->template push<gsVisitorMass<T> >();
51 }
52 else
53 {
54 gsVisitorMass<T> visitor(*m_pde_ptr);
55 this->apply(visitor, patchIndex);
56 }
57
58 // Assembly is done, compress the matrix
59 this->finalize();
60
61 return m_system.matrix();
62}
63
64template <class T>
65const gsSparseMatrix<T> & gsGenericAssembler<T>::assembleStiffness(const index_t patchIndex, const bool refresh)
66{
67 // Clean the sparse system
68 if (refresh)
69 {
71 m_system.reserve(m_bases[0], m_options, 1);
72 Base::computeDirichletDofs();
73 }
74
75 // Assemble stiffness integrals
76 if ( patchIndex < 0 )
77 {
78 this->template push<gsVisitorGradGrad<T> >();
79 }
80 else
81 {
82 gsVisitorGradGrad<T> visitor(*m_pde_ptr);
83 this->apply(visitor, patchIndex);
84 }
85
86 // Assembly is done, compress the matrix
87 this->finalize();
88
89 return m_system.matrix();
90}
91
92template <class T>
93const gsMatrix<T> & gsGenericAssembler<T>::assembleMoments(const gsFunction<T> & func, const index_t patchIndex, const bool refresh)
94{
95 // Clean the sparse system
96 if (refresh)
97 {
99 m_system.rhs().setZero(m_system.cols(), 1);
100 }
101
102 // Assemble moment integrals
103 gsVisitorMoments<T> mom(func);
104 if (patchIndex<0)
105 this->push(mom);
106 else
107 this->apply(mom,patchIndex);
108
109 // Assembly is done, compress the matrix
110 this->finalize();
111
112 return m_system.rhs();
113}
114
115template <class T>
117{
118 GISMO_ENSURE( m_options.getInt("InterfaceStrategy")==iFace::dg,
119 "Assembling DG terms only makes sense in corresponding setting." );
120
121 // Clean the sparse system
122 if (refresh)
123 {
125 m_system.reserve(m_bases[0], m_options, 1);
126 Base::computeDirichletDofs();
127 }
128
129 gsVisitorDg<T> visitor(*m_pde_ptr);
130 this->apply(visitor, iFace);
131
132 // Assembly is done, compress the matrix
133 this->finalize();
134
135 return m_system.matrix();
136}
137
138template <class T>
140{
141 GISMO_ASSERT( bc.ctype() == "Neumann", "gsGenericAssembler::assembleNeumann: Got " << bc.ctype() << "bc.");
142
143 // Clean the sparse system
144 if (refresh)
145 {
147 m_system.rhs().setZero(m_system.cols(), 1);
148 }
149
150 gsVisitorNeumann<T> visitor(*m_pde_ptr,bc);
151 this->apply(visitor, bc.patch(), bc.side());
152
153 // Assembly is done, compress the matrix
154 this->finalize();
155
156 return m_system.matrix();
157}
158
159template <class T>
161{
162 GISMO_ASSERT( bc.ctype() == "Dirichlet", "gsGenericAssembler::assembleNitsche: Got " << bc.ctype() << "bc.");
163
164 // Clean the sparse system
165 if (refresh)
166 {
168 m_system.reserve(m_bases[0], m_options, 1);
169 Base::computeDirichletDofs();
170 }
171
172
173 gsVisitorNitsche<T> visitor(*m_pde_ptr,bc);
174 this->apply(visitor, bc.patch(), bc.side());
175
176 // Assembly is done, compress the matrix
177 this->finalize();
178
179 return m_system.matrix();
180}
181
182
183} // namespace gismo
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsAssembler.hpp:30
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
const gsMatrix< T > & assembleMoments(const gsFunction< T > &func, index_t patchIndex=-1, bool refresh=true)
Moments assembly routine.
Definition gsGenericAssembler.hpp:93
const gsSparseMatrix< T > & assembleDG(const boundaryInterface &iFace, bool refresh=true)
Assemble dG interface terms.
Definition gsGenericAssembler.hpp:116
const gsSparseMatrix< T > & assembleNitsche(const boundary_condition< T > &bc, bool refresh=true)
Assemble Nitsche terms for weakly imposing Dirichlet conditions.
Definition gsGenericAssembler.hpp:160
void refresh()
Refreshes the sparse system (deletes assembled matrices and vectors)
Definition gsGenericAssembler.h:79
const gsSparseMatrix< T > & assembleStiffness(const index_t patchIndex=-1, const bool refresh=true)
Stiffness assembly routine.
Definition gsGenericAssembler.hpp:65
const gsSparseMatrix< T > & assembleNeumann(const boundary_condition< T > &bc, bool refresh=true)
Assemble Neumann boundary terms.
Definition gsGenericAssembler.hpp:139
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsGenericAssembler.hpp:28
const gsSparseMatrix< T > & assembleMass(const index_t patchIndex=-1, bool refresh=true)
Mass assembly routine.
Definition gsGenericAssembler.hpp:37
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void update(const gsOptionList &other, updateType type=ignoreIfUnknown)
Updates the object using the data from other.
Definition gsOptionList.cpp:253
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
Visitor for adding the interface conditions for the interior penalty methods of the Poisson problem.
Definition gsVisitorDg.h:52
The visitor computes element grad-grad integrals.
Definition gsVisitorGradGrad.h:34
The visitor computes element mass integrals.
Definition gsVisitorMass.h:30
Visitor for the moment vector of a function.
Definition gsVisitorMoments.h:30
Implementation of a Neumann BC for elliptic assemblers.
Definition gsVisitorNeumann.h:34
Visitor for adding the terms associated to weak (Nitsche-type) imposition of the Dirichlet boundary c...
Definition gsVisitorNitsche.h:55
Provides generic assembler routines.
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Visitor for adding the interface conditions for the interior penalty methods of the Poisson problem.
Stiffness (grad-grad) Visitor.
Element visitor for moment vector.
Neumann conditions visitor for elliptic problems.
Weak (Nitsche-type) imposition of the Dirichlet boundary conditions.
The G+Smo namespace, containing all definitions for the library.
Mass visitor for assembling element mass matrix.
Struct which represents an interface between two patches.
Definition gsBoundary.h:650
Class that defines a boundary condition for a side of a patch for some unknown variable of a PDE.
Definition gsBoundaryConditions.h:107
const std::string & ctype() const
Returns the type of the boundary condition.
Definition gsBoundaryConditions.h:263
boxSide side() const
Returns the side to which this boundary condition refers to.
Definition gsBoundaryConditions.h:269
index_t patch() const
Returns the patch to which this boundary condition refers to.
Definition gsBoundaryConditions.h:266