kelly_error_estimator.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 namespace libMesh
38 {
39 
40 
41 
42 void
44  ErrorVector&,
45  bool)
46 {
47  // Hang onto the system - we may need it for variable names later.
48  my_system = &system;
49 
50  // We'll need gradients and normal vectors for flux jump computation
51  fe_fine->get_dphi();
52  fe_fine->get_normals();
53  fe_coarse->get_dphi();
54 }
55 
56 
57 
58 void
60 {
61  Real error = 1.e-30;
62  unsigned int n_qp = fe_fine->n_quadrature_points();
63  unsigned int n_fine_dofs = Ufine.size();
64  unsigned int n_coarse_dofs = Ucoarse.size();
65 
66  std::vector<std::vector<RealGradient> > dphi_coarse = fe_coarse->get_dphi();
67  std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi();
68  std::vector<Point> face_normals = fe_fine->get_normals();
69  std::vector<Real> JxW_face = fe_fine->get_JxW();
70 
71  for (unsigned int qp=0; qp != n_qp; ++qp)
72  {
73  // Calculate solution gradients on fine and coarse elements
74  // at this quadrature point
75  Gradient grad_fine, grad_coarse;
76  for (unsigned int i=0; i != n_coarse_dofs; ++i)
77  grad_coarse.add_scaled (dphi_coarse[i][qp], Ucoarse(i));
78 
79  for (unsigned int i=0; i != n_fine_dofs; ++i)
80  grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i));
81 
82  // Find the jump in the normal derivative
83  // at this quadrature point
84  const Number jump = (grad_fine - grad_coarse)*face_normals[qp];
85  const Real jump2 = TensorTools::norm_sq(jump);
86 
87  // Accumulate the jump integral
88  error += JxW_face[qp] * jump2;
89  }
90 
91  // Add the h-weighted jump integral to each error term
92  fine_error =
93  error * fine_elem->hmax() * error_norm.weight(var);
94  coarse_error =
95  error * coarse_elem->hmax() * error_norm.weight(var);
96 }
97 
98 
99 bool
101 {
102  const std::string &var_name = my_system->variable_name(var);
103  const unsigned int n_fine_dofs = Ufine.size();
104 
105  std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi();
106  std::vector<Point> face_normals = fe_fine->get_normals();
107  std::vector<Real> JxW_face = fe_fine->get_JxW();
108  std::vector<Point> qface_point = fe_fine->get_xyz();
109 
110  // The reinitialization also recomputes the locations of
111  // the quadrature points on the side. By checking if the
112  // first quadrature point on the side is on a flux boundary
113  // for a particular variable, we will determine if the whole
114  // element is on a flux boundary (assuming quadrature points
115  // are strictly contained in the side).
116  if (this->_bc_function(*my_system, qface_point[0], var_name).first)
117  {
118  const Real h = fine_elem->hmax();
119 
120  // The number of quadrature points
121  const unsigned int n_qp = fe_fine->n_quadrature_points();
122 
123  // The error contribution from this face
124  Real error = 1.e-30;
125 
126  // loop over the integration points on the face.
127  for (unsigned int qp=0; qp<n_qp; qp++)
128  {
129  // Value of the imposed flux BC at this quadrature point.
130  const std::pair<bool,Real> flux_bc =
131  this->_bc_function(*my_system, qface_point[qp], var_name);
132 
133  // Be sure the BC function still thinks we're on the
134  // flux boundary.
135  libmesh_assert_equal_to (flux_bc.first, true);
136 
137  // The solution gradient from each element
138  Gradient grad_fine;
139 
140  // Compute the solution gradient on element e
141  for (unsigned int i=0; i != n_fine_dofs; i++)
142  grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i));
143 
144  // The difference between the desired BC and the approximate solution.
145  const Number jump = flux_bc.second - grad_fine*face_normals[qp];
146 
147  // The flux jump squared. If using complex numbers,
148  // TensorTools::norm_sq(z) returns |z|^2, where |z| is the modulus of z.
149  const Real jump2 = TensorTools::norm_sq(jump);
150 
151  // Integrate the error on the face. The error is
152  // scaled by an additional power of h, where h is
153  // the maximum side length for the element. This
154  // arises in the definition of the indicator.
155  error += JxW_face[qp]*jump2;
156 
157  } // End quadrature point loop
158 
159  fine_error = error*h*error_norm.weight(var);
160 
161  return true;
162  } // end if side on flux boundary
163  return false;
164 }
165 
166 
167 
168 void
169 KellyErrorEstimator::attach_flux_bc_function (std::pair<bool,Real> fptr(const System& system,
170  const Point& p,
171  const std::string& var_name))
172 {
173  _bc_function = fptr;
174 
175 // We may be turning boundary side integration on or off
176  if (fptr)
178  else
179  integrate_boundary_sides = false;
180 }
181 
182 } // namespace libMesh

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

Hosted By:
SourceForge.net Logo