exact_error_estimator.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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/exact_error_estimator.h" 00028 #include "libmesh/dof_map.h" 00029 #include "libmesh/equation_systems.h" 00030 #include "libmesh/error_vector.h" 00031 #include "libmesh/fe_base.h" 00032 #include "libmesh/libmesh_logging.h" 00033 #include "libmesh/elem.h" 00034 #include "libmesh/mesh_base.h" 00035 #include "libmesh/mesh_function.h" 00036 #include "libmesh/numeric_vector.h" 00037 #include "libmesh/quadrature.h" 00038 #include "libmesh/system.h" 00039 #include "libmesh/tensor_tools.h" 00040 00041 namespace libMesh 00042 { 00043 00044 //----------------------------------------------------------------- 00045 // ErrorEstimator implementations 00046 void ExactErrorEstimator::attach_exact_value (Number fptr(const Point& p, 00047 const Parameters& parameters, 00048 const std::string& sys_name, 00049 const std::string& unknown_name)) 00050 { 00051 libmesh_assert(fptr); 00052 _exact_value = fptr; 00053 00054 // We're not using a fine grid solution 00055 _equation_systems_fine = NULL; 00056 00057 // We're not using user-provided functors 00058 this->clear_functors(); 00059 } 00060 00061 00062 void ExactErrorEstimator::attach_exact_values (std::vector<FunctionBase<Number> *> f) 00063 { 00064 // Clear out any previous _exact_values entries, then add a new 00065 // entry for each system. 00066 for (unsigned int i=0; i != _exact_values.size(); ++i) 00067 delete (_exact_values[i]); 00068 00069 _exact_values.clear(); 00070 _exact_values.resize(f.size(), NULL); 00071 00072 // We use clone() to get non-sliced copies of FunctionBase 00073 // subclasses, but we can't put the resulting AutoPtrs into an STL 00074 // container. 00075 for (unsigned int i=0; i != f.size(); ++i) 00076 if (f[i]) 00077 _exact_values[i] = f[i]->clone().release(); 00078 } 00079 00080 00081 void ExactErrorEstimator::attach_exact_value (unsigned int sys_num, 00082 FunctionBase<Number> * f) 00083 { 00084 if (_exact_values.size() <= sys_num) 00085 _exact_values.resize(sys_num+1, NULL); 00086 00087 if (f) 00088 _exact_values[sys_num] = f->clone().release(); 00089 } 00090 00091 00092 void ExactErrorEstimator::attach_exact_deriv (Gradient gptr(const Point& p, 00093 const Parameters& parameters, 00094 const std::string& sys_name, 00095 const std::string& unknown_name)) 00096 { 00097 libmesh_assert(gptr); 00098 _exact_deriv = gptr; 00099 00100 // We're not using a fine grid solution 00101 _equation_systems_fine = NULL; 00102 00103 // We're not using user-provided functors 00104 this->clear_functors(); 00105 } 00106 00107 00108 void ExactErrorEstimator::attach_exact_derivs (std::vector<FunctionBase<Gradient> *> g) 00109 { 00110 // Clear out any previous _exact_derivs entries, then add a new 00111 // entry for each system. 00112 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00113 delete (_exact_derivs[i]); 00114 00115 _exact_derivs.clear(); 00116 _exact_derivs.resize(g.size(), NULL); 00117 00118 // We use clone() to get non-sliced copies of FunctionBase 00119 // subclasses, but we can't put the resulting AutoPtrs into an STL 00120 // container. 00121 for (unsigned int i=0; i != g.size(); ++i) 00122 if (g[i]) 00123 _exact_derivs[i] = g[i]->clone().release(); 00124 } 00125 00126 00127 void ExactErrorEstimator::attach_exact_deriv (unsigned int sys_num, 00128 FunctionBase<Gradient>* g) 00129 { 00130 if (_exact_derivs.size() <= sys_num) 00131 _exact_derivs.resize(sys_num+1, NULL); 00132 00133 if (g) 00134 _exact_derivs[sys_num] = g->clone().release(); 00135 } 00136 00137 00138 00139 00140 void ExactErrorEstimator::attach_exact_hessian (Tensor hptr(const Point& p, 00141 const Parameters& parameters, 00142 const std::string& sys_name, 00143 const std::string& unknown_name)) 00144 { 00145 libmesh_assert(hptr); 00146 _exact_hessian = hptr; 00147 00148 // We're not using a fine grid solution 00149 _equation_systems_fine = NULL; 00150 00151 // We're not using user-provided functors 00152 this->clear_functors(); 00153 } 00154 00155 00156 void ExactErrorEstimator::attach_exact_hessians (std::vector<FunctionBase<Tensor> *> h) 00157 { 00158 // Clear out any previous _exact_hessians entries, then add a new 00159 // entry for each system. 00160 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00161 delete (_exact_hessians[i]); 00162 00163 _exact_hessians.clear(); 00164 _exact_hessians.resize(h.size(), NULL); 00165 00166 // We use clone() to get non-sliced copies of FunctionBase 00167 // subclasses, but we can't put the resulting AutoPtrs into an STL 00168 // container. 00169 for (unsigned int i=0; i != h.size(); ++i) 00170 if (h[i]) 00171 _exact_hessians[i] = h[i]->clone().release(); 00172 } 00173 00174 00175 void ExactErrorEstimator::attach_exact_hessian (unsigned int sys_num, 00176 FunctionBase<Tensor>* h) 00177 { 00178 if (_exact_hessians.size() <= sys_num) 00179 _exact_hessians.resize(sys_num+1, NULL); 00180 00181 if (h) 00182 _exact_hessians[sys_num] = h->clone().release(); 00183 } 00184 00185 00186 void ExactErrorEstimator::attach_reference_solution (EquationSystems* es_fine) 00187 { 00188 libmesh_assert(es_fine); 00189 _equation_systems_fine = es_fine; 00190 00191 // If we're using a fine grid solution, we're not using exact value 00192 // function pointers or functors. 00193 _exact_value = NULL; 00194 _exact_deriv = NULL; 00195 _exact_hessian = NULL; 00196 00197 this->clear_functors(); 00198 } 00199 00200 #ifdef LIBMESH_ENABLE_AMR 00201 void ExactErrorEstimator::estimate_error (const System& system, 00202 ErrorVector& error_per_cell, 00203 const NumericVector<Number>* solution_vector, 00204 bool estimate_parent_error) 00205 #else 00206 void ExactErrorEstimator::estimate_error (const System& system, 00207 ErrorVector& error_per_cell, 00208 const NumericVector<Number>* solution_vector, 00209 bool /* estimate_parent_error */ ) 00210 #endif 00211 { 00212 // The current mesh 00213 const MeshBase& mesh = system.get_mesh(); 00214 00215 // The dimensionality of the mesh 00216 const unsigned int dim = mesh.mesh_dimension(); 00217 00218 // The number of variables in the system 00219 const unsigned int n_vars = system.n_vars(); 00220 00221 // The DofMap for this system 00222 const DofMap& dof_map = system.get_dof_map(); 00223 00224 // Resize the error_per_cell vector to be 00225 // the number of elements, initialize it to 0. 00226 error_per_cell.resize (mesh.max_elem_id()); 00227 std::fill (error_per_cell.begin(), error_per_cell.end(), 0.); 00228 00229 // Prepare current_local_solution to localize a non-standard 00230 // solution vector if necessary 00231 if (solution_vector && solution_vector != system.solution.get()) 00232 { 00233 NumericVector<Number>* newsol = 00234 const_cast<NumericVector<Number>*>(solution_vector); 00235 System &sys = const_cast<System&>(system); 00236 newsol->swap(*sys.solution); 00237 sys.update(); 00238 } 00239 00240 // Loop over all the variables in the system 00241 for (unsigned int var=0; var<n_vars; var++) 00242 { 00243 // Possibly skip this variable 00244 if (error_norm.weight(var) == 0.0) continue; 00245 00246 // The (string) name of this variable 00247 const std::string& var_name = system.variable_name(var); 00248 00249 // The type of finite element to use for this variable 00250 const FEType& fe_type = dof_map.variable_type (var); 00251 00252 AutoPtr<FEBase> fe (FEBase::build (dim, fe_type)); 00253 00254 // Build an appropriate Gaussian quadrature rule 00255 AutoPtr<QBase> qrule = 00256 fe_type.default_quadrature_rule (dim, 00257 _extra_order); 00258 00259 fe->attach_quadrature_rule (qrule.get()); 00260 00261 // Prepare a global solution and a MeshFunction of the fine system if we need one 00262 AutoPtr<MeshFunction> fine_values; 00263 AutoPtr<NumericVector<Number> > fine_soln = NumericVector<Number>::build(); 00264 if (_equation_systems_fine) 00265 { 00266 const System& fine_system = _equation_systems_fine->get_system(system.name()); 00267 00268 std::vector<Number> global_soln; 00269 // FIXME - we're assuming that the fine system solution gets 00270 // used even when a different vector is used for the coarse 00271 // system 00272 fine_system.update_global_solution(global_soln); 00273 fine_soln->init (global_soln.size(), true, SERIAL); 00274 (*fine_soln) = global_soln; 00275 00276 fine_values = AutoPtr<MeshFunction> 00277 (new MeshFunction(*_equation_systems_fine, 00278 *fine_soln, 00279 fine_system.get_dof_map(), 00280 fine_system.variable_number(var_name))); 00281 fine_values->init(); 00282 } else { 00283 // Initialize functors if we're using them 00284 for (unsigned int i=0; i != _exact_values.size(); ++i) 00285 if (_exact_values[i]) 00286 _exact_values[i]->init(); 00287 00288 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00289 if (_exact_derivs[i]) 00290 _exact_derivs[i]->init(); 00291 00292 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00293 if (_exact_hessians[i]) 00294 _exact_hessians[i]->init(); 00295 } 00296 00297 // Request the data we'll need to compute with 00298 fe->get_JxW(); 00299 fe->get_phi(); 00300 fe->get_dphi(); 00301 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00302 fe->get_d2phi(); 00303 #endif 00304 fe->get_xyz(); 00305 00306 // If we compute on parent elements, we'll want to do so only 00307 // once on each, so we need to keep track of which we've done. 00308 std::vector<bool> computed_var_on_parent; 00309 00310 #ifdef LIBMESH_ENABLE_AMR 00311 if (estimate_parent_error) 00312 computed_var_on_parent.resize(error_per_cell.size(), false); 00313 #endif 00314 00315 // TODO: this ought to be threaded (and using subordinate 00316 // MeshFunction objects in each thread rather than a single 00317 // master) 00318 00319 // Iterate over all the active elements in the mesh 00320 // that live on this processor. 00321 MeshBase::const_element_iterator 00322 elem_it = mesh.active_local_elements_begin(); 00323 const MeshBase::const_element_iterator 00324 elem_end = mesh.active_local_elements_end(); 00325 00326 for (;elem_it != elem_end; ++elem_it) 00327 { 00328 // e is necessarily an active element on the local processor 00329 const Elem* elem = *elem_it; 00330 const dof_id_type e_id = elem->id(); 00331 00332 #ifdef LIBMESH_ENABLE_AMR 00333 // See if the parent of element e has been examined yet; 00334 // if not, we may want to compute the estimator on it 00335 const Elem* parent = elem->parent(); 00336 00337 // We only can compute and only need to compute on 00338 // parents with all active children 00339 bool compute_on_parent = true; 00340 if (!parent || !estimate_parent_error) 00341 compute_on_parent = false; 00342 else 00343 for (unsigned int c=0; c != parent->n_children(); ++c) 00344 if (!parent->child(c)->active()) 00345 compute_on_parent = false; 00346 00347 if (compute_on_parent && 00348 !computed_var_on_parent[parent->id()]) 00349 { 00350 computed_var_on_parent[parent->id()] = true; 00351 00352 // Compute a projection onto the parent 00353 DenseVector<Number> Uparent; 00354 FEBase::coarsened_dof_values(*(system.current_local_solution), 00355 dof_map, parent, Uparent, 00356 var, false); 00357 00358 error_per_cell[parent->id()] += 00359 static_cast<ErrorVectorReal> 00360 (find_squared_element_error(system, var_name, 00361 parent, Uparent, 00362 fe.get(), 00363 fine_values.get())); 00364 } 00365 #endif 00366 00367 // Get the local to global degree of freedom maps 00368 std::vector<dof_id_type> dof_indices; 00369 dof_map.dof_indices (elem, dof_indices, var); 00370 const unsigned int n_dofs = 00371 libmesh_cast_int<unsigned int>(dof_indices.size()); 00372 DenseVector<Number> Uelem(n_dofs); 00373 for (unsigned int i=0; i != n_dofs; ++i) 00374 Uelem(i) = system.current_solution(dof_indices[i]); 00375 00376 error_per_cell[e_id] += 00377 static_cast<ErrorVectorReal> 00378 (find_squared_element_error(system, var_name, elem, 00379 Uelem, fe.get(), 00380 fine_values.get())); 00381 00382 } // End loop over active local elements 00383 } // End loop over variables 00384 00385 00386 00387 // Each processor has now computed the error contribuions 00388 // for its local elements. We need to sum the vector 00389 // and then take the square-root of each component. Note 00390 // that we only need to sum if we are running on multiple 00391 // processors, and we only need to take the square-root 00392 // if the value is nonzero. There will in general be many 00393 // zeros for the inactive elements. 00394 00395 // First sum the vector of estimated error values 00396 this->reduce_error(error_per_cell); 00397 00398 // Compute the square-root of each component. 00399 START_LOG("std::sqrt()", "ExactErrorEstimator"); 00400 for (dof_id_type i=0; i<error_per_cell.size(); i++) 00401 { 00402 00403 if (error_per_cell[i] != 0.) 00404 { 00405 libmesh_assert_greater (error_per_cell[i], 0.); 00406 error_per_cell[i] = std::sqrt(error_per_cell[i]); 00407 } 00408 00409 00410 } 00411 STOP_LOG("std::sqrt()", "ExactErrorEstimator"); 00412 00413 // If we used a non-standard solution before, now is the time to fix 00414 // the current_local_solution 00415 if (solution_vector && solution_vector != system.solution.get()) 00416 { 00417 NumericVector<Number>* newsol = 00418 const_cast<NumericVector<Number>*>(solution_vector); 00419 System &sys = const_cast<System&>(system); 00420 newsol->swap(*sys.solution); 00421 sys.update(); 00422 } 00423 } 00424 00425 00426 00427 Real ExactErrorEstimator::find_squared_element_error(const System& system, 00428 const std::string& var_name, 00429 const Elem *elem, 00430 const DenseVector<Number> &Uelem, 00431 FEBase *fe, 00432 MeshFunction *fine_values) const 00433 { 00434 // The (string) name of this system 00435 const std::string& sys_name = system.name(); 00436 const unsigned int sys_num = system.number(); 00437 00438 const unsigned int var = system.variable_number(var_name); 00439 const unsigned int var_component = 00440 system.variable_scalar_number(var, 0); 00441 00442 const Parameters& parameters = system.get_equation_systems().parameters; 00443 00444 // reinitialize the element-specific data 00445 // for the current element 00446 fe->reinit (elem); 00447 00448 // Get the data we need to compute with 00449 const std::vector<Real> & JxW = fe->get_JxW(); 00450 const std::vector<std::vector<Real> >& phi_values = fe->get_phi(); 00451 const std::vector<std::vector<RealGradient> >& dphi_values = fe->get_dphi(); 00452 const std::vector<Point>& q_point = fe->get_xyz(); 00453 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00454 const std::vector<std::vector<RealTensor> >& d2phi_values = fe->get_d2phi(); 00455 #endif 00456 00457 // The number of shape functions 00458 const unsigned int n_sf = 00459 libmesh_cast_int<unsigned int>(Uelem.size()); 00460 00461 // The number of quadrature points 00462 const unsigned int n_qp = 00463 libmesh_cast_int<unsigned int>(JxW.size()); 00464 00465 Real error_val = 0; 00466 00467 // Begin the loop over the Quadrature points. 00468 // 00469 for (unsigned int qp=0; qp<n_qp; qp++) 00470 { 00471 // Real u_h = 0.; 00472 // RealGradient grad_u_h; 00473 00474 Number u_h = 0.; 00475 00476 Gradient grad_u_h; 00477 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00478 Tensor grad2_u_h; 00479 #endif 00480 00481 // Compute solution values at the current 00482 // quadrature point. This reqiures a sum 00483 // over all the shape functions evaluated 00484 // at the quadrature point. 00485 for (unsigned int i=0; i<n_sf; i++) 00486 { 00487 // Values from current solution. 00488 u_h += phi_values[i][qp]*Uelem(i); 00489 grad_u_h += dphi_values[i][qp]*Uelem(i); 00490 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00491 grad2_u_h += d2phi_values[i][qp]*Uelem(i); 00492 #endif 00493 } 00494 00495 // Compute the value of the error at this quadrature point 00496 if (error_norm.type(var) == L2 || 00497 error_norm.type(var) == H1 || 00498 error_norm.type(var) == H2) 00499 { 00500 Number val_error = u_h; 00501 if (_exact_value) 00502 val_error -= _exact_value(q_point[qp],parameters,sys_name,var_name); 00503 else if (_exact_values.size() > sys_num && _exact_values[sys_num]) 00504 val_error -= _exact_values[sys_num]-> 00505 component(var_component, q_point[qp], system.time); 00506 else if (_equation_systems_fine) 00507 val_error -= (*fine_values)(q_point[qp]); 00508 00509 // Add the squares of the error to each contribution 00510 error_val += JxW[qp]*TensorTools::norm_sq(val_error); 00511 } 00512 00513 // Compute the value of the error in the gradient at this 00514 // quadrature point 00515 if (error_norm.type(var) == H1 || 00516 error_norm.type(var) == H1_SEMINORM || 00517 error_norm.type(var) == H2) 00518 { 00519 Gradient grad_error = grad_u_h; 00520 if(_exact_deriv) 00521 grad_error -= _exact_deriv(q_point[qp],parameters,sys_name,var_name); 00522 else if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num]) 00523 grad_error -= _exact_derivs[sys_num]-> 00524 component(var_component, q_point[qp], system.time); 00525 else if(_equation_systems_fine) 00526 grad_error -= fine_values->gradient(q_point[qp]); 00527 00528 error_val += JxW[qp]*grad_error.size_sq(); 00529 } 00530 00531 00532 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00533 // Compute the value of the error in the hessian at this 00534 // quadrature point 00535 if ((error_norm.type(var) == H2_SEMINORM || 00536 error_norm.type(var) == H2)) 00537 { 00538 Tensor grad2_error = grad2_u_h; 00539 if(_exact_hessian) 00540 grad2_error -= _exact_hessian(q_point[qp],parameters,sys_name,var_name); 00541 else if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num]) 00542 grad2_error -= _exact_hessians[sys_num]-> 00543 component(var_component, q_point[qp], system.time); 00544 else if (_equation_systems_fine) 00545 grad2_error -= fine_values->hessian(q_point[qp]); 00546 00547 error_val += JxW[qp]*grad2_error.size_sq(); 00548 } 00549 #endif 00550 00551 } // end qp loop 00552 00553 libmesh_assert_greater_equal (error_val, 0.); 00554 00555 return error_val; 00556 } 00557 00558 00559 00560 void ExactErrorEstimator::clear_functors() 00561 { 00562 // delete will clean up any cloned functors and no-op on any NULL 00563 // pointers 00564 00565 for (unsigned int i=0; i != _exact_values.size(); ++i) 00566 delete (_exact_values[i]); 00567 00568 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00569 delete (_exact_derivs[i]); 00570 00571 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00572 delete (_exact_hessians[i]); 00573 } 00574 00575 00576 00577 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: