Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/quadrature.h>
23 #include <deal.II/lac/full_matrix.h>
24 #include <deal.II/fe/mapping.h>
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/meshworker/dof_info.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace LocalIntegrators
31 {
39  namespace L2
40  {
52  template <int dim>
53  void mass_matrix (
55  const FEValuesBase<dim> &fe,
56  const double factor = 1.)
57  {
58  const unsigned int n_dofs = fe.dofs_per_cell;
59  const unsigned int n_components = fe.get_fe().n_components();
60 
61  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
62  {
63  const double dx = fe.JxW(k) * factor;
64  for (unsigned int i=0; i<n_dofs; ++i)
65  {
66  double Mii = 0.0;
67  for (unsigned int d=0; d<n_components; ++d)
68  Mii += dx
69  * fe.shape_value_component(i,k,d)
70  * fe.shape_value_component(i,k,d);
71 
72  M(i,i) += Mii;
73 
74  for (unsigned int j=i+1; j<n_dofs; ++j)
75  {
76  double Mij = 0.0;
77  for (unsigned int d=0; d<n_components; ++d)
78  Mij += dx
79  * fe.shape_value_component(j,k,d)
80  * fe.shape_value_component(i,k,d);
81 
82  M(i,j) += Mij;
83  M(j,i) += Mij;
84  }
85  }
86  }
87  }
88 
104  template <int dim>
107  const FEValuesBase<dim> &fe,
108  const std::vector<double> &weights)
109  {
110  const unsigned int n_dofs = fe.dofs_per_cell;
111  const unsigned int n_components = fe.get_fe().n_components();
112  AssertDimension(M.m(), n_dofs);
113  AssertDimension(M.n(), n_dofs);
114  AssertDimension(weights.size(), fe.n_quadrature_points);
115 
116  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
117  {
118  const double dx = fe.JxW(k) * weights[k];
119  for (unsigned int i=0; i<n_dofs; ++i)
120  {
121  double Mii = 0.0;
122  for (unsigned int d=0; d<n_components; ++d)
123  Mii += dx
124  * fe.shape_value_component(i,k,d)
125  * fe.shape_value_component(i,k,d);
126 
127  M(i,i) += Mii;
128 
129  for (unsigned int j=i+1; j<n_dofs; ++j)
130  {
131  double Mij = 0.0;
132  for (unsigned int d=0; d<n_components; ++d)
133  Mij += dx
134  * 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 L2 (
154  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 L2 (
179  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 * fe.shape_value_component(i,k,d) * input[d][k];
194  }
195 
208  template <int dim>
209  void jump_matrix (
210  FullMatrix<double> &M11,
211  FullMatrix<double> &M12,
212  FullMatrix<double> &M21,
213  FullMatrix<double> &M22,
214  const FEValuesBase<dim> &fe1,
215  const FEValuesBase<dim> &fe2,
216  const double factor1 = 1.,
217  const double factor2 = 1.)
218  {
219  const unsigned int n1_dofs = fe1.dofs_per_cell;
220  const unsigned int n2_dofs = fe2.dofs_per_cell;
221  const unsigned int n_components = fe1.get_fe().n_components();
222 
223  Assert(n1_dofs == n2_dofs, ExcNotImplemented());
224  (void) n2_dofs;
225  AssertDimension(n_components, fe2.get_fe().n_components());
226  AssertDimension(M11.m(), n1_dofs);
227  AssertDimension(M12.m(), n1_dofs);
228  AssertDimension(M21.m(), n2_dofs);
229  AssertDimension(M22.m(), n2_dofs);
230  AssertDimension(M11.n(), n1_dofs);
231  AssertDimension(M12.n(), n2_dofs);
232  AssertDimension(M21.n(), n1_dofs);
233  AssertDimension(M22.n(), n2_dofs);
234 
235  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
236  {
237  const double dx = fe1.JxW(k);
238 
239  for (unsigned int i=0; i<n1_dofs; ++i)
240  for (unsigned int j=0; j<n1_dofs; ++j)
241  for (unsigned int d=0; d<n_components; ++d)
242  {
243  const double u1 = factor1*fe1.shape_value_component(j,k,d);
244  const double u2 =-factor2*fe2.shape_value_component(j,k,d);
245  const double v1 = factor1*fe1.shape_value_component(i,k,d);
246  const double v2 =-factor2*fe2.shape_value_component(i,k,d);
247 
248  M11(i,j) += dx * u1*v1;
249  M12(i,j) += dx * u2*v1;
250  M21(i,j) += dx * u1*v2;
251  M22(i,j) += dx * u2*v2;
252  }
253  }
254  }
255  }
256 }
257 
258 DEAL_II_NAMESPACE_CLOSE
259 
260 #endif
size_type m() const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
const unsigned int dofs_per_cell
Definition: fe_values.h:1847
const FiniteElement< dim, spacedim > & get_fe() 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:209
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Library of integrals over cells and faces.
Definition: advection.h:30
void weighted_mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &weights)
Definition: l2.h:105
size_type n() const
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: l2.h:53
#define Assert(cond, exc)
Definition: exceptions.h:1142
const unsigned int n_quadrature_points
Definition: fe_values.h:1840
double JxW(const unsigned int quadrature_point) const
size_type size() 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:153