2#ifndef DUNE_PDELAB_LOCALOPERATOR_DARCYFEM_HH
3#define DUNE_PDELAB_LOCALOPERATOR_DARCYFEM_HH
5#include <dune/common/fvector.hh>
6#include <dune/common/shared_ptr.hh>
8#include <dune/geometry/referenceelements.hh>
23template<
typename P,
typename T,
typename X>
26 Dune::PDELab::GridFunctionTraits<
27 typename T::Traits::GridViewType,
28 typename T::Traits::FiniteElementType::Traits::LocalBasisType
29 ::Traits::RangeFieldType,
30 T::Traits::FiniteElementType::Traits::LocalBasisType::Traits
33 typename T::Traits::FiniteElementType::Traits
34 ::LocalBasisType::Traits::RangeFieldType,
35 T::Traits::FiniteElementType::Traits::LocalBasisType::Traits
37 DarcyVelocityFromHeadFEM<P,T,X> >
40 using LBTraits =
typename GFS::Traits::FiniteElementType::
41 Traits::LocalBasisType::Traits;
45 using LView =
typename X::template LocalView<LFSCache>;
49 typename GFS::Traits::GridViewType,
50 typename LBTraits::RangeFieldType,
53 typename LBTraits::RangeFieldType,
54 LBTraits::dimDomain> >;
68 : pgfs(stackobject_to_shared_ptr(gfs))
69 , pxg(stackobject_to_shared_ptr(x_))
70 , pp(stackobject_to_shared_ptr(
p))
86 std::vector<typename Traits::RangeFieldType> xl(lfs.
size());
87 lview.bind(lfs_cache);
92 auto geo = e.geometry();
95 const typename Traits::ElementType::Geometry::JacobianInverseTransposed
96 JgeoIT(geo.jacobianInverseTransposed(x));
99 std::vector<typename LBTraits::JacobianType> J(lfs.
size());
100 lfs.finiteElement().localBasis().evaluateJacobian(x,J);
104 for(
unsigned int i = 0; i < lfs.
size(); ++i) {
107 JgeoIT.umv(J[i][0], gradphi);
110 minusgrad.axpy(-xl[i], gradphi);
114 auto inside_cell_center_local = referenceElement(geo).position(0,0);
115 typename P::Traits::PermTensorType A(pp->A(e,inside_cell_center_local));
122 return pgfs->gridView();
126 std::shared_ptr<const GFS> pgfs;
127 std::shared_ptr<X> pxg;
128 std::shared_ptr<const P> pp;
130 mutable LFSCache lfs_cache;
const P & p
Definition: constraints.hh:148
Dune::FieldVector< GV::Grid::ctype, GV::dimension > DomainType
domain type in dim-size coordinates
Definition: function.hh:50
R RangeType
range type
Definition: function.hh:62
GV::Traits::template Codim< 0 >::Entity ElementType
codim 0 entity
Definition: function.hh:119
GV GridViewType
The type of the grid view the function lives on.
Definition: function.hh:116
traits class holding the function signature, same as in local function
Definition: function.hh:183
leaf of a function tree
Definition: function.hh:302
void update()
Definition: lfsindexcache.hh:304
Definition: lfsindexcache.hh:979
Traits::IndexContainer::size_type size() const
number of degrees of freedom contained in this lfs node
Definition: localfunctionspace.hh:223
Provide Darcy velocity as a vector-valued grid function.
Definition: darcyfem.hh:38
const Traits::GridViewType & getGridView() const
get a reference to the GridView
Definition: darcyfem.hh:120
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Definition: darcyfem.hh:77
DarcyVelocityFromHeadFEM(const P &p, const GFS &gfs, X &x_)
Construct a DarcyVelocityFromHeadFEM.
Definition: darcyfem.hh:67
Dune::PDELab::GridFunctionTraits< typename GFS::Traits::GridViewType, typename LBTraits::RangeFieldType, LBTraits::dimDomain, Dune::FieldVector< typename LBTraits::RangeFieldType, LBTraits::dimDomain > > Traits
Definition: darcyfem.hh:54