Amandus: Simulations based on multilevel Schwarz methods
darcy/noforce.h
Go to the documentation of this file.
1 /**********************************************************************
2  * Copyright (C) 2011 - 2014 by the authors
3  * Distributed under the MIT License
4  *
5  * See the files AUTHORS and LICENSE in the project root directory
6  *
7  **********************************************************************/
8 
9 #ifndef __darcy_noforce_h
10 #define __darcy_noforce_h
11 
12 #include <amandus/integrator.h>
13 #include <deal.II/base/polynomial.h>
14 #include <deal.II/base/tensor.h>
15 #include <deal.II/integrators/divergence.h>
16 #include <deal.II/integrators/l2.h>
17 #include <deal.II/integrators/laplace.h>
18 
19 using namespace dealii;
20 using namespace LocalIntegrators;
21 using namespace MeshWorker;
22 
23 namespace DarcyIntegrators
24 {
32 template <int dim>
34 {
35 public:
36  DarcyNoForceResidual();
37  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
38  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
39  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
40  IntegrationInfo<dim>& info2) const;
41 };
42 
43 //----------------------------------------------------------------------//
44 
45 template <int dim>
46 inline DarcyNoForceResidual<dim>::DarcyNoForceResidual()
47 {
48  this->input_vector_names.push_back("Newton iterate");
49 }
50 
51 template <int dim>
52 inline void
53 DarcyNoForceResidual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
54 {
55  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
56  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
57 
58  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), make_slice(info.values[0], 0, dim));
59  Divergence::gradient_residual(
60  dinfo.vector(0).block(0), info.fe_values(0), info.values[0][dim], -1.);
62  dinfo.vector(0).block(1), info.fe_values(1), make_slice(info.gradients[0], 0, dim), 1.);
63 }
64 
65 template <int dim>
66 inline void
67 DarcyNoForceResidual<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
68 {
69 }
70 
71 template <int dim>
72 inline void
73 DarcyNoForceResidual<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
74  IntegrationInfo<dim>&) const
75 {
76 }
77 }
78 
79 #endif
Definition: darcy/noforce.h:33
Definition: estimator.h:33
void cell_residual(dealii::Vector< number > &result, const dealii::FEValuesBase< dim > &fe, const dealii::VectorSlice< const std::vector< std::vector< dealii::Tensor< 1, dim >>>> &input, double lambda=0., double mu=1.)
Definition: elasticity/integrators.h:23
Definition: integrator.h:29