$extrastylesheet
discontinuity_measure.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 // C++ includes
00020 #include <algorithm> // for std::fill
00021 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
00022 #include <cmath>    // for sqrt
00023 
00024 
00025 // Local Includes
00026 #include "libmesh/libmesh_common.h"
00027 #include "libmesh/discontinuity_measure.h"
00028 #include "libmesh/error_vector.h"
00029 #include "libmesh/fe_base.h"
00030 #include "libmesh/libmesh_logging.h"
00031 #include "libmesh/elem.h"
00032 #include "libmesh/system.h"
00033 
00034 #include "libmesh/dense_vector.h"
00035 #include "libmesh/tensor_tools.h"
00036 
00037 
00038 namespace libMesh
00039 {
00040 
00041 void
00042 DiscontinuityMeasure::initialize(const System& system,
00043                                  ErrorVector&,
00044                                  bool)
00045 {
00046   // Hang onto the system - we may need it for variable names later.
00047   my_system = &system;
00048 
00049   // We'll need values for jump computation
00050   fe_fine->get_phi();
00051   fe_coarse->get_phi();
00052 }
00053 
00054 
00055 
00056 void
00057 DiscontinuityMeasure::internal_side_integration ()
00058 {
00059   Real error = 1.e-30;
00060   unsigned int n_qp = fe_fine->n_quadrature_points();
00061   unsigned int n_fine_dofs = Ufine.size();
00062   unsigned int n_coarse_dofs = Ucoarse.size();
00063 
00064   std::vector<std::vector<Real> > phi_coarse = fe_coarse->get_phi();
00065   std::vector<std::vector<Real> > phi_fine = fe_fine->get_phi();
00066   std::vector<Real> JxW_face = fe_fine->get_JxW();
00067 
00068   for (unsigned int qp=0; qp != n_qp; ++qp)
00069     {
00070       // Calculate solution values on fine and coarse elements
00071       // at this quadrature point
00072       Number u_fine = 0., u_coarse = 0.;
00073       for (unsigned int i=0; i != n_coarse_dofs; ++i)
00074         u_coarse += phi_coarse[i][qp] * Ucoarse(i);
00075 
00076       for (unsigned int i=0; i != n_fine_dofs; ++i)
00077         u_fine += phi_fine[i][qp] * Ufine(i);
00078 
00079       // Find the jump in the value
00080       // at this quadrature point
00081       const Number jump = u_fine - u_coarse;
00082       const Real jump2 = TensorTools::norm_sq(jump);
00083       // Accumulate the jump integral
00084       error += JxW_face[qp] * jump2;
00085     }
00086 
00087   // Add the h-weighted jump integral to each error term
00088   fine_error =
00089     error * fine_elem->hmax() * error_norm.weight(var);
00090   coarse_error =
00091     error * coarse_elem->hmax() * error_norm.weight(var);
00092 }
00093 
00094 
00095 bool
00096 DiscontinuityMeasure::boundary_side_integration ()
00097 {
00098   const std::string &var_name = my_system->variable_name(var);
00099 
00100   std::vector<std::vector<Real> > phi_fine = fe_fine->get_phi();
00101   std::vector<Real> JxW_face = fe_fine->get_JxW();
00102   std::vector<Point> qface_point = fe_fine->get_xyz();
00103 
00104   // The reinitialization also recomputes the locations of
00105   // the quadrature points on the side.  By checking if the
00106   // first quadrature point on the side is on an essential boundary
00107   // for a particular variable, we will determine if the whole
00108   // element is on an essential boundary (assuming quadrature points
00109   // are strictly contained in the side).
00110   if (this->_bc_function(*my_system, qface_point[0], var_name).first)
00111     {
00112       const Real h = fine_elem->hmax();
00113 
00114       // The number of quadrature points
00115       const unsigned int n_qp = fe_fine->n_quadrature_points();
00116 
00117       // The error contribution from this face
00118       Real error = 1.e-30;
00119 
00120       // loop over the integration points on the face.
00121       for (unsigned int qp=0; qp<n_qp; qp++)
00122         {
00123           // Value of the imposed essential BC at this quadrature point.
00124           const std::pair<bool,Real> essential_bc =
00125             this->_bc_function(*my_system, qface_point[qp], var_name);
00126 
00127           // Be sure the BC function still thinks we're on the
00128           // essential boundary.
00129           libmesh_assert_equal_to (essential_bc.first, true);
00130 
00131           // The solution gradient from each element
00132           Number u_fine = 0.;
00133 
00134           // Compute the solution gradient on element e
00135           for (unsigned int i=0; i != Ufine.size(); i++)
00136             u_fine += phi_fine[i][qp] * Ufine(i);
00137 
00138           // The difference between the desired BC and the approximate solution.
00139           const Number jump = essential_bc.second - u_fine;
00140 
00141           // The flux jump squared.  If using complex numbers,
00142           // norm_sq(z) returns |z|^2, where |z| is the modulus of z.
00143           const Real jump2 = TensorTools::norm_sq(jump);
00144 
00145           // Integrate the error on the face.  The error is
00146           // scaled by an additional power of h, where h is
00147           // the maximum side length for the element.  This
00148           // arises in the definition of the indicator.
00149           error += JxW_face[qp]*jump2;
00150 
00151         } // End quadrature point loop
00152 
00153       fine_error = error*h*error_norm.weight(var);
00154 
00155       return true;
00156     } // end if side on flux boundary
00157   return false;
00158 }
00159 
00160 
00161 
00162 void
00163 DiscontinuityMeasure::attach_essential_bc_function
00164 (std::pair<bool,Real> fptr(const System& system,
00165                            const Point& p,
00166                            const std::string& var_name))
00167 {
00168   _bc_function = fptr;
00169 
00170   // We may be turning boundary side integration on or off
00171   if (fptr)
00172     integrate_boundary_sides = true;
00173   else
00174     integrate_boundary_sides = false;
00175 }
00176 
00177 } // namespace libMesh