diff --git a/src/one-body-sample.cc b/src/one-body-sample.cc
index 6597700be83dc29ad328014b542ebafbbed250f6..c46797a9ced89dd47b7bf0ed6180d5c514bc51db 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;