Amandus: Simulations based on multilevel Schwarz methods
laplace/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 __laplace_noforce_h
10 #define __laplace_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 LaplaceIntegrators
24 {
30 template <int dim>
31 class NoForceRHS : public AmandusIntegrator<dim>
32 {
33 public:
34  NoForceRHS();
35 
36  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
37  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
38  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
39  IntegrationInfo<dim>& info2) const;
40 };
41 
47 template <int dim>
49 {
50 public:
52 
53  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
54  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
55  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
56  IntegrationInfo<dim>& info2) const;
57 };
58 
59 //----------------------------------------------------------------------//
60 
61 template <int dim>
63 {
64  this->use_boundary = false;
65  this->use_face = false;
66 }
67 
68 template <int dim>
69 void
70 NoForceRHS<dim>::cell(DoFInfo<dim>&, IntegrationInfo<dim>&) const
71 {
72 }
73 
74 template <int dim>
75 void
76 NoForceRHS<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
77 {
78 }
79 
80 template <int dim>
81 void
82 NoForceRHS<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
83  IntegrationInfo<dim>&) const
84 {
85 }
86 
87 //----------------------------------------------------------------------//
88 
89 template <int dim>
91 {
92  this->input_vector_names.push_back("Newton iterate");
93 }
94 
95 template <int dim>
96 void
97 NoForceResidual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
98 {
99  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
100  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
101 
102  std::vector<std::vector<double>> rhs(1,
103  std::vector<double>(info.fe_values(0).n_quadrature_points));
104  Laplace::cell_residual(dinfo.vector(0).block(0), info.fe_values(0), info.gradients[0][0]);
105 }
106 
107 template <int dim>
108 void
109 NoForceResidual<dim>::boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
110 {
111  std::vector<double> null(info.fe_values(0).n_quadrature_points, 0.);
112 
113  const unsigned int deg = info.fe_values(0).get_fe().tensor_degree();
114  Laplace::nitsche_residual(dinfo.vector(0).block(0),
115  info.fe_values(0),
116  info.values[0][0],
117  info.gradients[0][0],
118  null,
119  Laplace::compute_penalty(dinfo, dinfo, deg, deg));
120 }
121 
122 template <int dim>
123 void
124 NoForceResidual<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
125  IntegrationInfo<dim>& info2) const
126 {
127  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
128  Laplace::ip_residual(dinfo1.vector(0).block(0),
129  dinfo2.vector(0).block(0),
130  info1.fe_values(0),
131  info2.fe_values(0),
132  info1.values[0][0],
133  info1.gradients[0][0],
134  info2.values[0][0],
135  info2.gradients[0][0],
136  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg));
137 }
138 }
139 
140 #endif
Definition: laplace/noforce.h:48
Definition: laplace/eigen.h:19
Definition: laplace/noforce.h:31
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