Skip to content
Snippets Groups Projects
Commit b8a6fe93 authored by Carsten Gräser's avatar Carsten Gräser
Browse files

Switch function interface in `NewPFEAssembler`

* Switch to dune-functions interface for grid functions
* Adjust `newpfeassemblertest.cc`

After this patch you can no longer pass a `Dune::VirtualFunction` or `VirtualGridFunction`,
but have to use either the `std::function` or the `Dune::Functions::GridViewFunction` interface.
parent d229e777
Branches
Tags
1 merge request!146Stop using Dune::VirtualFunction interface
......@@ -14,6 +14,7 @@
#include <dune/fufem/assemblers/localoperatorassembler.hh>
#include <dune/fufem/staticmatrixtools.hh>
#include <dune/fufem/symmetrictensor.hh>
#include <dune/fufem/functions/gridfunctionhelper.hh>
//! \todo Please doc me !
......@@ -23,7 +24,12 @@ class NewPFEAssembler : public LocalOperatorAssembler < GridType, TrialLocalFE,
private:
static const int dim = GridType::dimension;
const FunctionType& phaseField_;
using GlobalCoordinate = typename GridType::template Codim<0>::Geometry::GlobalCoordinate;
using LocalCoordinate = typename GridType::template Codim<0>::Geometry::LocalCoordinate;
using PhaseFieldRangeType = std::decay_t<std::invoke_result_t<FunctionType,GlobalCoordinate>>;
using LocalFunction = Dune::Fufem::Impl::LocalFunctionInterfaceForGrid<GridType, PhaseFieldRangeType>;
mutable LocalFunction localPhaseField_;
const AlloyType& alloy_;
const int quadOrder_;
......@@ -36,8 +42,9 @@ class NewPFEAssembler : public LocalOperatorAssembler < GridType, TrialLocalFE,
// typedef typename LocalOperatorAssembler < GridType ,T >::BoolMatrix BoolMatrix;
// typedef typename Dune::Matrix<T> LocalMatrix;
NewPFEAssembler(const FunctionType& phaseField, const AlloyType& alloy, const int quadOrder=2):
phaseField_(phaseField),
template<class FT>
NewPFEAssembler(const FT& phaseField, const AlloyType& alloy, const int quadOrder=2):
localPhaseField_(Dune::Fufem::Impl::localFunctionForGrid<GridType, PhaseFieldRangeType>(phaseField)),
alloy_(alloy),
quadOrder_(quadOrder)
{}
......@@ -67,6 +74,8 @@ class NewPFEAssembler : public LocalOperatorAssembler < GridType, TrialLocalFE,
std::vector<JacobianType> referenceGradients(tFE.localBasis().size());
std::vector<FVdim> gradients(tFE.localBasis().size());
localPhaseField_.bind(element);
// loop over quadrature points
for (size_t pt=0; pt < quad.size(); ++pt)
{
......@@ -81,8 +90,7 @@ class NewPFEAssembler : public LocalOperatorAssembler < GridType, TrialLocalFE,
const double integrationElement = geometry.integrationElement(quadPos);
// Value of the phase field at this quadrature point
typename FunctionType::RangeType phaseValue;
phaseField_.evaluateLocal(element, quadPos, phaseValue);
PhaseFieldRangeType phaseValue = localPhaseField_(quadPos);
// ncompspec !!
//Dune::FieldVector<double,2> phaseFieldVector((1.0+phaseValue)*0.5); /* this line for order parameter u \in [-1,1] */
......
......@@ -4,7 +4,6 @@
#include <dune/common/parallel/mpihelper.hh>
#include <dune/fufem/assemblers/assembler.hh>
#include <dune/fufem/functions/basisgridfunction.hh>
#include <dune/fufem/functionspacebases/dunefunctionsbasis.hh>
#include <dune/fufem/mechanics/cubictensor.hh>
#include <dune/fufem/mechanics/isotropictensor.hh>
......@@ -13,6 +12,7 @@
#include <dune/fufem/test/common.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>
#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
template <int n>
struct eigenStrainComputer {
......@@ -73,7 +73,9 @@ void runTest() {
using Basis = DuneFunctionsBasis<Dune::Functions::LagrangeBasis<GridView, 1> >;
using LocalBasis = typename Basis::LocalFiniteElement;
using Function = BasisGridFunction<Basis, Vector>;
using GlobalCoordinate = typename Grid::template Codim<0>::Geometry::GlobalCoordinate;
using PhaseFieldRange = typename Vector::value_type;
using Function = PhaseFieldRange(GlobalCoordinate);
using MaterialAssembler =
NewPFEAssembler<Grid, LocalBasis, LocalBasis, Function, Alloy>;
......@@ -81,7 +83,7 @@ void runTest() {
Basis basis(grid->leafGridView());
Vector coefficients(basis.size());
coefficients = 1.0;
Function oneFunction(basis, coefficients);
auto oneFunction = Dune::Functions::makeDiscreteGlobalBasisFunction<PhaseFieldRange>(basis.functionsBasis(), coefficients);
Assembler<Basis, Basis> assembler(basis, basis);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment