26 gsElPoissonAssembler<T>::gsElPoissonAssembler(
const gsMultiPatch<T> & patches,
27 const gsMultiBasis<T> & basis,
28 const gsBoundaryConditions<T> & bconditions,
29 const gsFunction<T> & body_force)
31 gsPiecewiseFunction<T> rightHandSides;
32 rightHandSides.addPiece(body_force);
33 typename gsPde<T>::Ptr pde(
new gsPoissonPde<T>(patches,bconditions,rightHandSides) );
34 m_bases.push_back(basis);
35 Base::initialize(pde, m_bases, defaultOptions());
39 gsOptionList gsElPoissonAssembler<T>::defaultOptions()
41 gsOptionList opt = Base::defaultOptions();
42 opt.addReal(
"LocalStiff",
"Stiffening degree for the Jacobian-based local stiffening",0.);
47 void gsElPoissonAssembler<T>::refresh()
49 std::vector<gsDofMapper> m_dofMappers(m_bases.size());
50 m_bases[0].getMapper((dirichlet::strategy)m_options.getInt(
"DirichletStrategy"),
51 iFace::glue,m_pde_ptr->bc(),m_dofMappers[0],0,
true);
53 m_system = gsSparseSystem<T>(m_dofMappers[0]);
54 m_system.reserve(m_bases[0], m_options, m_pde_ptr->numRhs());
55 Base::computeDirichletDofs(0);
59 void gsElPoissonAssembler<T>::assemble(
bool saveEliminationMatrix)
61 m_system.matrix().setZero();
62 m_system.reserve(m_bases[0], m_options, m_pde_ptr->numRhs());
63 m_system.rhs().setZero(Base::numDofs(),m_pde_ptr->numRhs());
65 if (saveEliminationMatrix)
67 eliminationMatrix.resize(Base::numDofs(),Base::numFixedDofs());
68 eliminationMatrix.setZero();
69 eliminationMatrix.reservePerColumn(m_system.numColNz(m_bases[0],m_options));
72 gsVisitorElPoisson<T> visitor(*m_pde_ptr, saveEliminationMatrix ? &eliminationMatrix :
nullptr);
73 Base::template push<gsVisitorElPoisson<T> >(visitor);
75 m_system.matrix().makeCompressed();
77 if (saveEliminationMatrix)
79 Base::rhsWithZeroDDofs = m_system.rhs();
80 eliminationMatrix.makeCompressed();
85 void gsElPoissonAssembler<T>::constructSolution(
const gsMatrix<T> & solVector,
86 const std::vector<gsMatrix<T> > & fixedDoFs,
87 gsMultiPatch<T> & result)
const
89 Base::constructSolution(solVector,fixedDoFs,result,gsVector<index_t>::Zero(1));
Visitor class for Poisson's equation.
Provides stiffness matrix for Poisson's equations.