From 8bccc9f0943c50a4b787aee153ccf071e2f35c12 Mon Sep 17 00:00:00 2001 From: Elias Pipping <elias.pipping@fu-berlin.de> Date: Thu, 3 Nov 2011 14:31:39 +0100 Subject: [PATCH] Make frictional boundary term assembly into a func. --- src/one-body-sample.cc | 35 ++++++++++++++++++++++------------- 1 file changed, 22 insertions(+), 13 deletions(-) diff --git a/src/one-body-sample.cc b/src/one-body-sample.cc index 6597700b..c46797a9 100644 --- a/src/one-body-sample.cc +++ b/src/one-body-sample.cc @@ -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; -- GitLab