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

Use the assembler class

parent 183c2a63
No related branches found
No related tags found
No related merge requests found
......@@ -46,13 +46,12 @@
#include <dune/istl/bvector.hh>
#include <dune/istl/scaledidmatrix.hh>
#include <dune/fufem/assemblers/functionalassembler.hh>
#include <dune/fufem/assemblers/assembler.hh>
#include <dune/fufem/assemblers/localassemblers/boundarymassassembler.hh>
#include <dune/fufem/assemblers/localassemblers/l2functionalassembler.hh>
#include <dune/fufem/assemblers/localassemblers/massassembler.hh>
#include <dune/fufem/assemblers/localassemblers/stvenantkirchhoffassembler.hh>
#include <dune/fufem/assemblers/localassemblers/vonmisesstressassembler.hh>
#include <dune/fufem/assemblers/operatorassembler.hh>
#include <dune/fufem/boundarypatch.hh>
#include <dune/fufem/dunepython.hh>
#include <dune/fufem/functions/basisgridfunction.hh>
......@@ -204,6 +203,8 @@ int main(int argc, char *argv[]) {
typedef P1NodalBasis<GridView, double> P1Basis;
P0Basis const p0Basis(leafView);
P1Basis const p1Basis(leafView);
Assembler<P0Basis, P0Basis> const p0Assembler(p0Basis, p0Basis);
Assembler<P1Basis, P1Basis> const p1Assembler(p1Basis, p1Basis);
// Set up the boundary
Dune::BitSetVector<dims> ignoreNodes(finestSize, false);
......@@ -260,8 +261,7 @@ int main(int argc, char *argv[]) {
MassAssembler<GridType, P1Basis::LocalFiniteElement,
P1Basis::LocalFiniteElement,
Dune::ScaledIdentityMatrix<double, dims>> const localMass;
OperatorAssembler<P1Basis, P1Basis>(p1Basis, p1Basis)
.assemble(localMass, massMatrix);
p1Assembler.assembleOperator(localMass, massMatrix);
massMatrix *= density;
}
{
......@@ -286,8 +286,8 @@ int main(int argc, char *argv[]) {
weightedGravitationalDirection);
L2FunctionalAssembler<GridType, SmallVector> gravityFunctionalAssembler(
gravityFunction);
FunctionalAssembler<P1Basis>(p1Basis)
.assemble(gravityFunctionalAssembler, gravityFunctional, true);
p1Assembler.assembleFunctional(gravityFunctionalAssembler,
gravityFunctional);
}
}
SingletonVectorType surfaceNormalStress(finestSize);
......@@ -300,8 +300,7 @@ int main(int argc, char *argv[]) {
StVenantKirchhoffAssembler<GridType, P1Basis::LocalFiniteElement,
P1Basis::LocalFiniteElement> const
localStiffness(E, nu);
OperatorAssembler<P1Basis, P1Basis>(p1Basis, p1Basis)
.assemble(localStiffness, stiffnessMatrix);
p1Assembler.assembleOperator(localStiffness, stiffnessMatrix);
}
SingletonMatrixType frictionalBoundaryMassMatrix;
......@@ -315,9 +314,8 @@ int main(int argc, char *argv[]) {
GridType, BoundaryPatch<GridView>, P1Basis::LocalFiniteElement,
P1Basis::LocalFiniteElement, SmallSingletonMatrix> const
frictionalBoundaryMassAssembler(frictionalBoundary);
OperatorAssembler<P1Basis, P1Basis>(p1Basis, p1Basis).assemble(
frictionalBoundaryMassAssembler, frictionalBoundaryMassMatrix);
p1Assembler.assembleOperator(frictionalBoundaryMassAssembler,
frictionalBoundaryMassMatrix);
}
SumNorm<VectorType> const velocityEnergyNorm(1.0, stiffnessMatrixNorm, 1.0,
massMatrixNorm);
......@@ -509,8 +507,7 @@ int main(int argc, char *argv[]) {
p1Basis, u);
VonMisesStressAssembler<GridType> localStressAssembler(
E, nu, gridDisplacement);
FunctionalAssembler<P0Basis>(p0Basis)
.assemble(localStressAssembler, vonMisesStress, true);
p0Assembler.assembleFunctional(localStressAssembler, vonMisesStress);
writeVtk<P1Basis, P0Basis, VectorType, SingletonVectorType, GridView>(
p1Basis, u, alpha, p0Basis, vonMisesStress, leafView,
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment