Reference documentation for deal.II version 9.1.0-pre
l2.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_integrators_l2_h
17 #define dealii_integrators_l2_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/quadrature.h>
24 
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/fe/mapping.h>
27 
28 #include <deal.II/lac/full_matrix.h>
29 
30 #include <deal.II/meshworker/dof_info.h>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 namespace LocalIntegrators
35 {
43  namespace L2
44  {
56  template <int dim>
57  void
59  const FEValuesBase<dim> &fe,
60  const double factor = 1.)
61  {
62  const unsigned int n_dofs = fe.dofs_per_cell;
63  const unsigned int n_components = fe.get_fe().n_components();
64 
65  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
66  {
67  const double dx = fe.JxW(k) * factor;
68  for (unsigned int i = 0; i < n_dofs; ++i)
69  {
70  double Mii = 0.0;
71  for (unsigned int d = 0; d < n_components; ++d)
72  Mii += dx * fe.shape_value_component(i, k, d) *
73  fe.shape_value_component(i, k, d);
74 
75  M(i, i) += Mii;
76 
77  for (unsigned int j = i + 1; j < n_dofs; ++j)
78  {
79  double Mij = 0.0;
80  for (unsigned int d = 0; d < n_components; ++d)
81  Mij += dx * fe.shape_value_component(j, k, d) *
82  fe.shape_value_component(i, k, d);
83 
84  M(i, j) += Mij;
85  M(j, i) += Mij;
86  }
87  }
88  }
89  }
90 
106  template <int dim>
107  void
109  const FEValuesBase<dim> & fe,
110  const std::vector<double> &weights)
111  {
112  const unsigned int n_dofs = fe.dofs_per_cell;
113  const unsigned int n_components = fe.get_fe().n_components();
114  AssertDimension(M.m(), n_dofs);
115  AssertDimension(M.n(), n_dofs);
116  AssertDimension(weights.size(), fe.n_quadrature_points);
117 
118  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
119  {
120  const double dx = fe.JxW(k) * weights[k];
121  for (unsigned int i = 0; i < n_dofs; ++i)
122  {
123  double Mii = 0.0;
124  for (unsigned int d = 0; d < n_components; ++d)
125  Mii += dx * fe.shape_value_component(i, k, d) *
126  fe.shape_value_component(i, k, d);
127 
128  M(i, i) += Mii;
129 
130  for (unsigned int j = i + 1; j < n_dofs; ++j)
131  {
132  double Mij = 0.0;
133  for (unsigned int d = 0; d < n_components; ++d)
134  Mij += dx * fe.shape_value_component(j, k, d) *
135  fe.shape_value_component(i, k, d);
136 
137  M(i, j) += Mij;
138  M(j, i) += Mij;
139  }
140  }
141  }
142  }
143 
152  template <int dim, typename number>
153  void
154  L2(Vector<number> & result,
155  const FEValuesBase<dim> & fe,
156  const std::vector<double> &input,
157  const double factor = 1.)
158  {
159  const unsigned int n_dofs = fe.dofs_per_cell;
160  AssertDimension(result.size(), n_dofs);
161  AssertDimension(fe.get_fe().n_components(), 1);
162  AssertDimension(input.size(), fe.n_quadrature_points);
163 
164  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
165  for (unsigned int i = 0; i < n_dofs; ++i)
166  result(i) += fe.JxW(k) * factor * input[k] * fe.shape_value(i, k);
167  }
168 
177  template <int dim, typename number>
178  void
179  L2(Vector<number> & result,
180  const FEValuesBase<dim> & fe,
181  const VectorSlice<const std::vector<std::vector<double>>> &input,
182  const double factor = 1.)
183  {
184  const unsigned int n_dofs = fe.dofs_per_cell;
185  const unsigned int n_components = input.size();
186 
187  AssertDimension(result.size(), n_dofs);
188  AssertDimension(input.size(), fe.get_fe().n_components());
189 
190  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
191  for (unsigned int i = 0; i < n_dofs; ++i)
192  for (unsigned int d = 0; d < n_components; ++d)
193  result(i) += fe.JxW(k) * factor *
194  fe.shape_value_component(i, k, d) * input[d][k];
195  }
196 
209  template <int dim>
210  void
212  FullMatrix<double> & M12,
213  FullMatrix<double> & M21,
214  FullMatrix<double> & M22,
215  const FEValuesBase<dim> &fe1,
216  const FEValuesBase<dim> &fe2,
217  const double factor1 = 1.,
218  const double factor2 = 1.)
219  {
220  const unsigned int n1_dofs = fe1.dofs_per_cell;
221  const unsigned int n2_dofs = fe2.dofs_per_cell;
222  const unsigned int n_components = fe1.get_fe().n_components();
223 
224  Assert(n1_dofs == n2_dofs, ExcNotImplemented());
225  (void)n2_dofs;
226  AssertDimension(n_components, fe2.get_fe().n_components());
227  AssertDimension(M11.m(), n1_dofs);
228  AssertDimension(M12.m(), n1_dofs);
229  AssertDimension(M21.m(), n2_dofs);
230  AssertDimension(M22.m(), n2_dofs);
231  AssertDimension(M11.n(), n1_dofs);
232  AssertDimension(M12.n(), n2_dofs);
233  AssertDimension(M21.n(), n1_dofs);
234  AssertDimension(M22.n(), n2_dofs);
235 
236  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
237  {
238  const double dx = fe1.JxW(k);
239 
240  for (unsigned int i = 0; i < n1_dofs; ++i)
241  for (unsigned int j = 0; j < n1_dofs; ++j)
242  for (unsigned int d = 0; d < n_components; ++d)
243  {
244  const double u1 =
245  factor1 * fe1.shape_value_component(j, k, d);
246  const double u2 =
247  -factor2 * fe2.shape_value_component(j, k, d);
248  const double v1 =
249  factor1 * fe1.shape_value_component(i, k, d);
250  const double v2 =
251  -factor2 * fe2.shape_value_component(i, k, d);
252 
253  M11(i, j) += dx * u1 * v1;
254  M12(i, j) += dx * u2 * v1;
255  M21(i, j) += dx * u1 * v2;
256  M22(i, j) += dx * u2 * v2;
257  }
258  }
259  }
260  } // namespace L2
261 } // namespace LocalIntegrators
262 
263 DEAL_II_NAMESPACE_CLOSE
264 
265 #endif
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
const unsigned int dofs_per_cell
Definition: fe_values.h:2033
size_type size() const
void jump_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const double factor1=1., const double factor2=1.)
Definition: l2.h:211
Library of integrals over cells and faces.
Definition: advection.h:34
void weighted_mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &weights)
Definition: l2.h:108
size_type n() const
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: l2.h:58
#define Assert(cond, exc)
Definition: exceptions.h:1227
const FiniteElement< dim, spacedim > & get_fe() const
size_type m() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
static::ExceptionBase & ExcNotImplemented()
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
Definition: l2.h:154
double JxW(const unsigned int quadrature_point) const