2 #ifndef DUNE_PDELAB_CONVECTIONDIFFUSION_HH
3 #define DUNE_PDELAB_CONVECTIONDIFFUSION_HH
7 #include<dune/common/exceptions.hh>
8 #include<dune/common/fvector.hh>
9 #include<dune/geometry/type.hh>
11 #include<dune/geometry/referenceelements.hh>
12 #include<dune/geometry/quadraturerules.hh>
38 template<
typename GV,
typename RF>
55 typedef Dune::FieldVector<DomainFieldType,dimDomain>
DomainType;
64 typedef Dune::FieldVector<RF,GV::dimensionworld>
RangeType;
67 typedef Dune::FieldMatrix<RangeFieldType,dimDomain,dimDomain>
PermTensorType;
70 typedef typename GV::Traits::template Codim<0>::Entity
ElementType;
75 template<
class T,
class Imp>
82 typename Traits::RangeFieldType
83 f (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x,
84 typename Traits::RangeFieldType u)
const
86 return asImp().f(e,x,u);
90 typename Traits::RangeFieldType
91 w (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x,
92 typename Traits::RangeFieldType u)
const
94 return asImp().w(e,x,u);
98 typename Traits::RangeFieldType
99 v (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x,
100 typename Traits::RangeFieldType u)
const
102 return asImp().v(e,x,u);
106 typename Traits::PermTensorType
107 D (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x)
const
109 return asImp().D(e,x);
113 typename Traits::RangeType
114 q (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x,
115 typename Traits::RangeFieldType u)
const
117 return asImp().q(e,x,u);
122 const I & intersection,
123 const Dune::FieldVector<typename I::ctype, I::dimension-1> & coord
126 return asImp().isDirichlet( intersection, coord );
130 typename Traits::RangeFieldType
131 g (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x)
const
133 return asImp().g(e,x);
139 typename Traits::RangeFieldType
140 j (
const typename Traits::ElementType&
e,
const typename Traits::DomainType& x,
141 typename Traits::RangeFieldType u)
const
143 return asImp().j(e,x);
147 Imp& asImp () {
return static_cast<Imp &
> (*this);}
148 const Imp& asImp ()
const {
return static_cast<const Imp &
>(*this);}
160 const typename T::Traits::GridViewType gv;
171 const I & intersection,
172 const Dune::FieldVector<typename I::ctype, I::dimension-1> & coord
175 return t.isDirichlet( intersection, coord );
186 :
public GridFunctionBase<Dune::PDELab::GridFunctionTraits<typename T::Traits::GridViewType,
187 typename T::Traits::RangeFieldType,
188 1,Dune::FieldVector<typename T::Traits::RangeFieldType,1> >
189 ,DirichletBoundaryCondition_CD<T> >
193 typename T::Traits::RangeFieldType,
194 1,Dune::FieldVector<typename T::Traits::RangeFieldType,1> >
Traits;
242 : param(param_), intorder(intorder_)
246 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename R>
247 void alpha_volume (
const EG&
eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, R& r)
const
250 typedef typename LFSU::Traits::FiniteElementType::
251 Traits::LocalBasisType::Traits::DomainFieldType DF;
252 typedef typename LFSU::Traits::FiniteElementType::
253 Traits::LocalBasisType::Traits::RangeFieldType RF;
254 typedef typename LFSU::Traits::FiniteElementType::
255 Traits::LocalBasisType::Traits::JacobianType JacobianType;
256 typedef typename LFSU::Traits::FiniteElementType::
257 Traits::LocalBasisType::Traits::RangeType RangeType;
259 typedef typename LFSU::Traits::SizeType size_type;
262 const int dim = EG::Geometry::dimension;
265 Dune::GeometryType gt = eg.geometry().type();
266 const Dune::QuadratureRule<DF,dim>& rule = Dune::QuadratureRules<DF,dim>::rule(gt,intorder);
269 typename T::Traits::PermTensorType tensor;
270 Dune::FieldVector<DF,dim> localcenter = Dune::ReferenceElements<DF,dim>::general(gt).position(0,0);
271 tensor = param.D(eg.entity(),localcenter);
274 std::vector<typename T::Traits::RangeFieldType>
w(lfsu.size());
275 for (size_type i=0; i<lfsu.size(); i++)
276 w[i] = param.w(eg.entity(),localcenter,x(lfsu,i));
279 for (
typename Dune::QuadratureRule<DF,dim>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
282 std::vector<RangeType> phi(lfsu.size());
283 lfsu.finiteElement().localBasis().evaluateFunction(it->position(),phi);
287 for (size_type i=0; i<lfsu.size(); i++)
291 typename T::Traits::RangeFieldType f = param.f(eg.entity(),it->position(),u);
294 typename T::Traits::RangeType q = param.q(eg.entity(),it->position(),u);
297 std::vector<JacobianType> js(lfsu.size());
298 lfsu.finiteElement().localBasis().evaluateJacobian(it->position(),js);
301 const typename EG::Geometry::JacobianInverseTransposed jac
302 = eg.geometry().jacobianInverseTransposed(it->position());
303 std::vector<Dune::FieldVector<RF,dim> > gradphi(lfsu.size());
304 for (size_type i=0; i<lfsu.size(); i++)
307 jac.umv(js[i][0],gradphi[i]);
311 Dune::FieldVector<RF,dim> vgradu(0.0);
312 for (size_type i=0; i<lfsu.size(); i++)
313 vgradu.axpy(
w[i],gradphi[i]);
314 vgradu *= param.v(eg.entity(),it->position(),u);
317 Dune::FieldVector<RF,dim> Dvgradu(0.0);
318 tensor.umv(vgradu,Dvgradu);
321 RF factor = it->weight() * eg.geometry().integrationElement(it->position());
322 for (size_type i=0; i<lfsu.size(); i++)
323 r.accumulate(lfsu,i,( Dvgradu*gradphi[i] - q*gradphi[i] - f*phi[i] )*factor);
328 template<
typename IG,
typename LFSU,
typename X,
typename LFSV,
typename R>
330 const LFSU& lfsu_s,
const X& x_s,
const LFSV& lfsv_s,
334 typedef typename LFSV::Traits::FiniteElementType::
335 Traits::LocalBasisType::Traits::DomainFieldType DF;
336 typedef typename LFSV::Traits::FiniteElementType::
337 Traits::LocalBasisType::Traits::RangeFieldType RF;
338 typedef typename LFSV::Traits::FiniteElementType::
339 Traits::LocalBasisType::Traits::RangeType RangeType;
341 typedef typename LFSV::Traits::SizeType size_type;
344 const int dim = IG::dimension;
347 Dune::GeometryType gtface = ig.geometryInInside().type();
348 const Dune::QuadratureRule<DF,dim-1>& rule = Dune::QuadratureRules<DF,dim-1>::rule(gtface,intorder);
351 Dune::FieldVector<DF,dim-1> facecenterlocal = Dune::ReferenceElements<DF,dim-1>::general(gtface).position(0,0);
352 Dune::FieldVector<DF,dim> facecenterinelement = ig.geometryInInside().global( facecenterlocal );
353 std::vector<typename T::Traits::RangeFieldType>
w(lfsu_s.size());
354 auto cell_inside = ig.inside();
355 for (size_type i=0; i<lfsu_s.size(); i++)
356 w[i] = param.w(cell_inside,facecenterinelement,x_s(lfsu_s,i));
359 for (
typename Dune::QuadratureRule<DF,dim-1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
363 if( param.isDirichlet( ig.intersection(), it->position() ) )
367 Dune::FieldVector<DF,dim> local = ig.geometryInInside().global(it->position());
370 std::vector<RangeType> phi(lfsv_s.size());
371 lfsv_s.finiteElement().localBasis().evaluateFunction(local,phi);
375 for (size_type i=0; i<lfsu_s.size(); i++)
379 typename T::Traits::RangeFieldType j;
380 j = param.j(cell_inside,local,u);
383 RF factor = it->weight()*ig.geometry().integrationElement(it->position());
384 for (size_type i=0; i<lfsv_s.size(); i++)
385 r_s.accumulate(lfsu_s,i,j*phi[i]*factor);
GV::Traits::template Codim< 0 >::Entity ElementType
grid types
Definition: convectiondiffusion.hh:70
Dune::FieldVector< RF, GV::dimensionworld > RangeType
range type
Definition: convectiondiffusion.hh:64
GV GridViewType
the grid view
Definition: convectiondiffusion.hh:43
void alpha_boundary(const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, R &r_s) const
Definition: convectiondiffusion.hh:329
GV::Intersection IntersectionType
Definition: convectiondiffusion.hh:71
Definition: convectiondiffusion.hh:157
BCTypeParam_CD(const typename T::Traits::GridViewType &gv_, const T &t_)
Definition: convectiondiffusion.hh:164
dimension of the domain
Definition: convectiondiffusion.hh:48
base class for parameter class
Definition: convectiondiffusion.hh:76
RF RangeFieldType
Export type for range field.
Definition: convectiondiffusion.hh:61
Traits::RangeFieldType f(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeFieldType u) const
source/reaction term
Definition: convectiondiffusion.hh:83
const E & e
Definition: interpolate.hh:172
Implement jacobian_apply_volume() based on alpha_volume()
Definition: defaultimp.hh:321
Dune::FieldMatrix< RangeFieldType, dimDomain, dimDomain > PermTensorType
permeability tensor type
Definition: convectiondiffusion.hh:67
void alpha_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
Definition: convectiondiffusion.hh:247
Dune::FieldVector< DomainFieldType, dimDomain-1 > IntersectionDomainType
domain type
Definition: convectiondiffusion.hh:58
GV GridViewType
The type of the grid view the function lives on.
Definition: function.hh:114
Definition: convectiondiffusion.hh:235
bool isDirichlet(const I &intersection, const Dune::FieldVector< typename I::ctype, I::dimension-1 > &coord) const
Definition: convectiondiffusion.hh:170
const IG & ig
Definition: constraints.hh:147
Definition: constraintsparameters.hh:24
GV::Grid::ctype DomainFieldType
Export type for domain field.
Definition: convectiondiffusion.hh:52
leaf of a function tree
Definition: function.hh:576
Default flags for all local operators.
Definition: flags.hh:18
void setTime(double t)
set time in parameter class
Definition: convectiondiffusion.hh:390
Dune::PDELab::GridFunctionTraits< typename T::Traits::GridViewType, typename T::Traits::RangeFieldType, 1, Dune::FieldVector< typename T::Traits::RangeFieldType, 1 > > Traits
Definition: convectiondiffusion.hh:194
Definition: convectiondiffusion.hh:238
sparsity pattern generator
Definition: pattern.hh:13
Traits::RangeType q(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeFieldType u) const
nonlinear flux vector
Definition: convectiondiffusion.hh:114
Traits::PermTensorType D(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
tensor diffusion coefficient
Definition: convectiondiffusion.hh:107
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Evaluate the GridFunction at given position.
Definition: convectiondiffusion.hh:200
Dune::FieldVector< GV::Grid::ctype, GV::dimension > DomainType
domain type in dim-size coordinates
Definition: function.hh:48
Definition: convectiondiffusion.hh:239
T Traits
Definition: convectiondiffusion.hh:79
DirichletBoundaryCondition_CD(const typename Traits::GridViewType &g_, const T &t_)
constructor
Definition: convectiondiffusion.hh:197
traits class holding the function signature, same as in local function
Definition: function.hh:175
Definition: convectiondiffusion.hh:224
Definition: adaptivity.hh:27
Traits::RangeFieldType w(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeFieldType u) const
nonlinearity under gradient
Definition: convectiondiffusion.hh:91
ConvectionDiffusion(T ¶m_, int intorder_=2)
Definition: convectiondiffusion.hh:241
Implement jacobian_volume() based on alpha_volume()
Definition: defaultimp.hh:32
Implement jacobian_boundary() based on alpha_boundary()
Definition: defaultimp.hh:251
Traits::RangeFieldType g(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
Dirichlet boundary condition value.
Definition: convectiondiffusion.hh:131
static const int dim
Definition: adaptivity.hh:83
GV::Traits::template Codim< 0 >::Entity ElementType
codim 0 entity
Definition: function.hh:117
Traits::RangeFieldType v(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeFieldType u) const
scalar nonlinearity in diffusion coefficient
Definition: convectiondiffusion.hh:99
R RangeType
range type
Definition: function.hh:60
const Traits::GridViewType & getGridView() const
Definition: convectiondiffusion.hh:207
Traits::RangeFieldType j(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeFieldType u) const
Neumann boundary condition.
Definition: convectiondiffusion.hh:140
VTKWriter & w
Definition: function.hh:1101
Definition: convectiondiffusion.hh:185
Dune::FieldVector< DomainFieldType, dimDomain > DomainType
domain type
Definition: convectiondiffusion.hh:55
const EG & eg
Definition: constraints.hh:280
Default class for additional methods in instationary local operators.
Definition: idefault.hh:89
Implement jacobian_apply_boundary() based on alpha_boundary()
Definition: defaultimp.hh:538
bool isDirichlet(const I &intersection, const Dune::FieldVector< typename I::ctype, I::dimension-1 > &coord) const
Definition: convectiondiffusion.hh:121
traits class for two phase parameter class
Definition: convectiondiffusion.hh:40