Skip to content
Snippets Groups Projects
Commit 8bccc9f0 authored by Elias Pipping's avatar Elias Pipping Committed by Elias Pipping
Browse files

Make frictional boundary term assembly into a func.

parent 3fa3365e
No related branches found
No related tags found
No related merge requests found
......@@ -137,6 +137,26 @@ void assemble_neumann(GridType const &grid, FEBasis const &feBasis,
true); // resize the output vector and zero all of its entries
}
// Assembles constant 1-function on frictional boundary in nodalIntegrals
template <class GridType, class LocalVectorType, class FEBasis>
void assemble_frictional(
GridType const &grid, FEBasis const &feBasis,
Dune::BitSetVector<1> const &frictionalNodes,
std::vector<Dune::FieldVector<double, 1>> &nodalIntegrals) {
BoundaryPatch<typename GridType::LeafGridView> frictionalBoundary(
grid.leafView(), frictionalNodes);
ConstantFunction<LocalVectorType, Dune::FieldVector<double, 1>>
constantOneFunction(1);
NeumannBoundaryAssembler<GridType, Dune::FieldVector<double, 1>>
frictionalBoundaryAssembler(constantOneFunction);
BoundaryFunctionalAssembler<FEBasis> boundaryFunctionalAssembler(
feBasis, frictionalBoundary);
boundaryFunctionalAssembler.assemble(
frictionalBoundaryAssembler, nodalIntegrals,
true); // resize the output vector and zero all of its entries
}
int main() {
try {
typedef Dune::FieldVector<double, dim> SmallVector;
......@@ -190,20 +210,9 @@ int main() {
neumannNodes, f);
{
// constant 1-function on frictional boundary
BoundaryPatch<GridType::LeafGridView> frictionalBoundary(grid.leafView(),
frictionalNodes);
std::vector<Dune::FieldVector<double, 1>> nodalIntegrals;
ConstantFunction<SmallVector, Dune::FieldVector<double, 1>>
constantOneFunction(1);
NeumannBoundaryAssembler<GridType, Dune::FieldVector<double, 1>>
frictionalBoundaryAssembler(constantOneFunction);
BoundaryFunctionalAssembler<P1Basis> boundaryFunctionalAssembler(
p1Basis, frictionalBoundary);
boundaryFunctionalAssembler.assemble(
frictionalBoundaryAssembler, nodalIntegrals,
true); // resize the output vector and zero all of its entries
assemble_frictional<GridType, SmallVector, P1Basis>(
grid, p1Basis, frictionalNodes, nodalIntegrals);
// TODO: populate on S_F
std::vector<double> normalStress;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment