discontinuity_measure.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 #include <algorithm> // for std::fill
21 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
22 #include <cmath> // for sqrt
23 
24 
25 // Local Includes
26 #include "libmesh/libmesh_common.h"
28 #include "libmesh/error_vector.h"
29 #include "libmesh/fe_base.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/system.h"
33 
34 #include "libmesh/dense_vector.h"
35 #include "libmesh/tensor_tools.h"
36 
37 
38 namespace libMesh
39 {
40 
41 void
43  ErrorVector&,
44  bool)
45 {
46  // Hang onto the system - we may need it for variable names later.
47  my_system = &system;
48 
49  // We'll need values for jump computation
50  fe_fine->get_phi();
51  fe_coarse->get_phi();
52 }
53 
54 
55 
56 void
58 {
59  Real error = 1.e-30;
60  unsigned int n_qp = fe_fine->n_quadrature_points();
61  unsigned int n_fine_dofs = Ufine.size();
62  unsigned int n_coarse_dofs = Ucoarse.size();
63 
64  std::vector<std::vector<Real> > phi_coarse = fe_coarse->get_phi();
65  std::vector<std::vector<Real> > phi_fine = fe_fine->get_phi();
66  std::vector<Real> JxW_face = fe_fine->get_JxW();
67 
68  for (unsigned int qp=0; qp != n_qp; ++qp)
69  {
70  // Calculate solution values on fine and coarse elements
71  // at this quadrature point
72  Number u_fine = 0., u_coarse = 0.;
73  for (unsigned int i=0; i != n_coarse_dofs; ++i)
74  u_coarse += phi_coarse[i][qp] * Ucoarse(i);
75 
76  for (unsigned int i=0; i != n_fine_dofs; ++i)
77  u_fine += phi_fine[i][qp] * Ufine(i);
78 
79  // Find the jump in the value
80  // at this quadrature point
81  const Number jump = u_fine - u_coarse;
82  const Real jump2 = TensorTools::norm_sq(jump);
83  // Accumulate the jump integral
84  error += JxW_face[qp] * jump2;
85  }
86 
87  // Add the h-weighted jump integral to each error term
88  fine_error =
89  error * fine_elem->hmax() * error_norm.weight(var);
90  coarse_error =
91  error * coarse_elem->hmax() * error_norm.weight(var);
92 }
93 
94 
95 bool
97 {
98  const std::string &var_name = my_system->variable_name(var);
99 
100  std::vector<std::vector<Real> > phi_fine = fe_fine->get_phi();
101  std::vector<Real> JxW_face = fe_fine->get_JxW();
102  std::vector<Point> qface_point = fe_fine->get_xyz();
103 
104  // The reinitialization also recomputes the locations of
105  // the quadrature points on the side. By checking if the
106  // first quadrature point on the side is on an essential boundary
107  // for a particular variable, we will determine if the whole
108  // element is on an essential boundary (assuming quadrature points
109  // are strictly contained in the side).
110  if (this->_bc_function(*my_system, qface_point[0], var_name).first)
111  {
112  const Real h = fine_elem->hmax();
113 
114  // The number of quadrature points
115  const unsigned int n_qp = fe_fine->n_quadrature_points();
116 
117  // The error contribution from this face
118  Real error = 1.e-30;
119 
120  // loop over the integration points on the face.
121  for (unsigned int qp=0; qp<n_qp; qp++)
122  {
123  // Value of the imposed essential BC at this quadrature point.
124  const std::pair<bool,Real> essential_bc =
125  this->_bc_function(*my_system, qface_point[qp], var_name);
126 
127  // Be sure the BC function still thinks we're on the
128  // essential boundary.
129  libmesh_assert_equal_to (essential_bc.first, true);
130 
131  // The solution gradient from each element
132  Number u_fine = 0.;
133 
134  // Compute the solution gradient on element e
135  for (unsigned int i=0; i != Ufine.size(); i++)
136  u_fine += phi_fine[i][qp] * Ufine(i);
137 
138  // The difference between the desired BC and the approximate solution.
139  const Number jump = essential_bc.second - u_fine;
140 
141  // The flux jump squared. If using complex numbers,
142  // norm_sq(z) returns |z|^2, where |z| is the modulus of z.
143  const Real jump2 = TensorTools::norm_sq(jump);
144 
145  // Integrate the error on the face. The error is
146  // scaled by an additional power of h, where h is
147  // the maximum side length for the element. This
148  // arises in the definition of the indicator.
149  error += JxW_face[qp]*jump2;
150 
151  } // End quadrature point loop
152 
153  fine_error = error*h*error_norm.weight(var);
154 
155  return true;
156  } // end if side on flux boundary
157  return false;
158 }
159 
160 
161 
162 void
164  (std::pair<bool,Real> fptr(const System& system,
165  const Point& p,
166  const std::string& var_name))
167 {
168  _bc_function = fptr;
169 
170 // We may be turning boundary side integration on or off
171  if (fptr)
172  integrate_boundary_sides = true;
173  else
174  integrate_boundary_sides = false;
175 }
176 
177 } // namespace libMesh

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:04 UTC

Hosted By:
SourceForge.net Logo