$extrastylesheet
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 // C++ includes 00019 00020 00021 // Local includes 00022 #include "libmesh/dof_map.h" 00023 #include "libmesh/elem.h" 00024 #include "libmesh/exact_solution.h" 00025 #include "libmesh/equation_systems.h" 00026 #include "libmesh/fe_base.h" 00027 #include "libmesh/function_base.h" 00028 #include "libmesh/mesh_base.h" 00029 #include "libmesh/mesh_function.h" 00030 #include "libmesh/numeric_vector.h" 00031 #include "libmesh/parallel.h" 00032 #include "libmesh/quadrature.h" 00033 #include "libmesh/wrapped_function.h" 00034 #include "libmesh/fe_interface.h" 00035 #include "libmesh/raw_accessor.h" 00036 #include "libmesh/tensor_tools.h" 00037 00038 namespace libMesh 00039 { 00040 00041 ExactSolution::ExactSolution(const EquationSystems& es) : 00042 _equation_systems(es), 00043 _equation_systems_fine(NULL), 00044 _extra_order(0) 00045 { 00046 // Initialize the _errors data structure which holds all 00047 // the eventual values of the error. 00048 for (unsigned int sys=0; sys<_equation_systems.n_systems(); ++sys) 00049 { 00050 // Reference to the system 00051 const System& system = _equation_systems.get_system(sys); 00052 00053 // The name of the system 00054 const std::string& sys_name = system.name(); 00055 00056 // The SystemErrorMap to be inserted 00057 ExactSolution::SystemErrorMap sem; 00058 00059 for (unsigned int var=0; var<system.n_vars(); ++var) 00060 { 00061 // The name of this variable 00062 const std::string& var_name = system.variable_name(var); 00063 sem[var_name] = std::vector<Real>(5, 0.); 00064 } 00065 00066 _errors[sys_name] = sem; 00067 } 00068 } 00069 00070 00071 ExactSolution::~ExactSolution() 00072 { 00073 // delete will clean up any cloned functors and no-op on any NULL 00074 // pointers 00075 00076 for (unsigned int i=0; i != _exact_values.size(); ++i) 00077 delete (_exact_values[i]); 00078 00079 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00080 delete (_exact_derivs[i]); 00081 00082 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00083 delete (_exact_hessians[i]); 00084 } 00085 00086 00087 void ExactSolution::attach_reference_solution (const EquationSystems* es_fine) 00088 { 00089 libmesh_assert(es_fine); 00090 _equation_systems_fine = es_fine; 00091 00092 // If we're using a fine grid solution, we're not using exact values 00093 _exact_values.clear(); 00094 _exact_derivs.clear(); 00095 _exact_hessians.clear(); 00096 } 00097 00098 00099 void ExactSolution::attach_exact_value (Number fptr(const Point& p, 00100 const Parameters& parameters, 00101 const std::string& sys_name, 00102 const std::string& unknown_name)) 00103 { 00104 libmesh_assert(fptr); 00105 00106 // Clear out any previous _exact_values entries, then add a new 00107 // entry for each system. 00108 _exact_values.clear(); 00109 00110 for (unsigned int sys=0; sys<_equation_systems.n_systems(); ++sys) 00111 { 00112 const System& system = _equation_systems.get_system(sys); 00113 _exact_values.push_back 00114 (new WrappedFunction<Number> 00115 (system, fptr, &_equation_systems.parameters)); 00116 } 00117 00118 // If we're using exact values, we're not using a fine grid solution 00119 _equation_systems_fine = NULL; 00120 } 00121 00122 00123 void ExactSolution::attach_exact_values (std::vector<FunctionBase<Number> *> f) 00124 { 00125 // Clear out any previous _exact_values entries, then add a new 00126 // entry for each system. 00127 for (unsigned int i=0; i != _exact_values.size(); ++i) 00128 delete (_exact_values[i]); 00129 00130 _exact_values.clear(); 00131 _exact_values.resize(f.size(), NULL); 00132 00133 // We use clone() to get non-sliced copies of FunctionBase 00134 // subclasses, but we don't currently put the resulting UniquePtrs 00135 // into an STL container. 00136 for (unsigned int i=0; i != f.size(); ++i) 00137 if (f[i]) 00138 _exact_values[i] = f[i]->clone().release(); 00139 } 00140 00141 00142 void ExactSolution::attach_exact_value (unsigned int sys_num, 00143 FunctionBase<Number> * f) 00144 { 00145 if (_exact_values.size() <= sys_num) 00146 _exact_values.resize(sys_num+1, NULL); 00147 00148 if (f) 00149 _exact_values[sys_num] = f->clone().release(); 00150 } 00151 00152 00153 void ExactSolution::attach_exact_deriv (Gradient gptr(const Point& p, 00154 const Parameters& parameters, 00155 const std::string& sys_name, 00156 const std::string& unknown_name)) 00157 { 00158 libmesh_assert(gptr); 00159 00160 // Clear out any previous _exact_derivs entries, then add a new 00161 // entry for each system. 00162 _exact_derivs.clear(); 00163 00164 for (unsigned int sys=0; sys<_equation_systems.n_systems(); ++sys) 00165 { 00166 const System& system = _equation_systems.get_system(sys); 00167 _exact_derivs.push_back 00168 (new WrappedFunction<Gradient> 00169 (system, gptr, &_equation_systems.parameters)); 00170 } 00171 00172 // If we're using exact values, we're not using a fine grid solution 00173 _equation_systems_fine = NULL; 00174 } 00175 00176 00177 void ExactSolution::attach_exact_derivs (std::vector<FunctionBase<Gradient> *> g) 00178 { 00179 // Clear out any previous _exact_derivs entries, then add a new 00180 // entry for each system. 00181 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00182 delete (_exact_derivs[i]); 00183 00184 _exact_derivs.clear(); 00185 _exact_derivs.resize(g.size(), NULL); 00186 00187 // We use clone() to get non-sliced copies of FunctionBase 00188 // subclasses, but we don't currently put the resulting UniquePtrs 00189 // into an STL container. 00190 for (unsigned int i=0; i != g.size(); ++i) 00191 if (g[i]) 00192 _exact_derivs[i] = g[i]->clone().release(); 00193 } 00194 00195 00196 void ExactSolution::attach_exact_deriv (unsigned int sys_num, 00197 FunctionBase<Gradient>* g) 00198 { 00199 if (_exact_derivs.size() <= sys_num) 00200 _exact_derivs.resize(sys_num+1, NULL); 00201 00202 if (g) 00203 _exact_derivs[sys_num] = g->clone().release(); 00204 } 00205 00206 00207 void ExactSolution::attach_exact_hessian (Tensor hptr(const Point& p, 00208 const Parameters& parameters, 00209 const std::string& sys_name, 00210 const std::string& unknown_name)) 00211 { 00212 libmesh_assert(hptr); 00213 00214 // Clear out any previous _exact_hessians entries, then add a new 00215 // entry for each system. 00216 _exact_hessians.clear(); 00217 00218 for (unsigned int sys=0; sys<_equation_systems.n_systems(); ++sys) 00219 { 00220 const System& system = _equation_systems.get_system(sys); 00221 _exact_hessians.push_back 00222 (new WrappedFunction<Tensor> 00223 (system, hptr, &_equation_systems.parameters)); 00224 } 00225 00226 // If we're using exact values, we're not using a fine grid solution 00227 _equation_systems_fine = NULL; 00228 } 00229 00230 00231 void ExactSolution::attach_exact_hessians (std::vector<FunctionBase<Tensor> *> h) 00232 { 00233 // Clear out any previous _exact_hessians entries, then add a new 00234 // entry for each system. 00235 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00236 delete (_exact_hessians[i]); 00237 00238 _exact_hessians.clear(); 00239 _exact_hessians.resize(h.size(), NULL); 00240 00241 // We use clone() to get non-sliced copies of FunctionBase 00242 // subclasses, but we don't currently put the resulting UniquePtrs 00243 // into an STL container. 00244 for (unsigned int i=0; i != h.size(); ++i) 00245 if (h[i]) 00246 _exact_hessians[i] = h[i]->clone().release(); 00247 } 00248 00249 00250 void ExactSolution::attach_exact_hessian (unsigned int sys_num, 00251 FunctionBase<Tensor>* h) 00252 { 00253 if (_exact_hessians.size() <= sys_num) 00254 _exact_hessians.resize(sys_num+1, NULL); 00255 00256 if (h) 00257 _exact_hessians[sys_num] = h->clone().release(); 00258 } 00259 00260 00261 std::vector<Real>& ExactSolution::_check_inputs(const std::string& sys_name, 00262 const std::string& unknown_name) 00263 { 00264 // If no exact solution function or fine grid solution has been 00265 // attached, we now just compute the solution norm (i.e. the 00266 // difference from an "exact solution" of zero 00267 00268 // Make sure the requested sys_name exists. 00269 std::map<std::string, SystemErrorMap>::iterator sys_iter = 00270 _errors.find(sys_name); 00271 00272 if (sys_iter == _errors.end()) 00273 libmesh_error_msg("Sorry, couldn't find the requested system '" << sys_name << "'."); 00274 00275 // Make sure the requested unknown_name exists. 00276 SystemErrorMap::iterator var_iter = (*sys_iter).second.find(unknown_name); 00277 00278 if (var_iter == (*sys_iter).second.end()) 00279 libmesh_error_msg("Sorry, couldn't find the requested variable '" << unknown_name << "'."); 00280 00281 // Return a reference to the proper error entry 00282 return (*var_iter).second; 00283 } 00284 00285 00286 00287 void ExactSolution::compute_error(const std::string& sys_name, 00288 const std::string& unknown_name) 00289 { 00290 // Check the inputs for validity, and get a reference 00291 // to the proper location to store the error 00292 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00293 unknown_name); 00294 00295 libmesh_assert( _equation_systems.has_system(sys_name) ); 00296 const System& sys = _equation_systems.get_system<System>( sys_name ); 00297 00298 libmesh_assert( sys.has_variable( unknown_name ) ); 00299 switch( FEInterface::field_type(sys.variable_type( unknown_name )) ) 00300 { 00301 case TYPE_SCALAR: 00302 { 00303 this->_compute_error<Real>(sys_name, 00304 unknown_name, 00305 error_vals); 00306 break; 00307 } 00308 case TYPE_VECTOR: 00309 { 00310 this->_compute_error<RealGradient>(sys_name, 00311 unknown_name, 00312 error_vals); 00313 break; 00314 } 00315 default: 00316 libmesh_error_msg("Invalid variable type!"); 00317 } 00318 } 00319 00320 00321 00322 00323 00324 Real ExactSolution::error_norm(const std::string& sys_name, 00325 const std::string& unknown_name, 00326 const FEMNormType& norm) 00327 { 00328 // Check the inputs for validity, and get a reference 00329 // to the proper location to store the error 00330 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00331 unknown_name); 00332 00333 libmesh_assert(_equation_systems.has_system(sys_name)); 00334 libmesh_assert(_equation_systems.get_system(sys_name).has_variable( unknown_name )); 00335 const FEType& fe_type = _equation_systems.get_system(sys_name).variable_type(unknown_name); 00336 00337 switch (norm) 00338 { 00339 case L2: 00340 return std::sqrt(error_vals[0]); 00341 case H1: 00342 return std::sqrt(error_vals[0] + error_vals[1]); 00343 case H2: 00344 return std::sqrt(error_vals[0] + error_vals[1] + error_vals[2]); 00345 case HCURL: 00346 { 00347 if(FEInterface::field_type(fe_type) == TYPE_SCALAR) 00348 libmesh_error_msg("Cannot compute HCurl error norm of scalar-valued variables!"); 00349 else 00350 return std::sqrt(error_vals[0] + error_vals[5]); 00351 } 00352 case HDIV: 00353 { 00354 if(FEInterface::field_type(fe_type) == TYPE_SCALAR) 00355 libmesh_error_msg("Cannot compute HDiv error norm of scalar-valued variables!"); 00356 else 00357 return std::sqrt(error_vals[0] + error_vals[6]); 00358 } 00359 case H1_SEMINORM: 00360 return std::sqrt(error_vals[1]); 00361 case H2_SEMINORM: 00362 return std::sqrt(error_vals[2]); 00363 case HCURL_SEMINORM: 00364 { 00365 if(FEInterface::field_type(fe_type) == TYPE_SCALAR) 00366 libmesh_error_msg("Cannot compute HCurl error seminorm of scalar-valued variables!"); 00367 else 00368 return std::sqrt(error_vals[5]); 00369 } 00370 case HDIV_SEMINORM: 00371 { 00372 if(FEInterface::field_type(fe_type) == TYPE_SCALAR) 00373 libmesh_error_msg("Cannot compute HDiv error seminorm of scalar-valued variables!"); 00374 else 00375 return std::sqrt(error_vals[6]); 00376 } 00377 case L1: 00378 return error_vals[3]; 00379 case L_INF: 00380 return error_vals[4]; 00381 00382 default: 00383 libmesh_error_msg("Currently only Sobolev norms/seminorms are supported!"); 00384 } 00385 } 00386 00387 00388 00389 00390 00391 00392 00393 Real ExactSolution::l2_error(const std::string& sys_name, 00394 const std::string& unknown_name) 00395 { 00396 00397 // Check the inputs for validity, and get a reference 00398 // to the proper location to store the error 00399 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00400 unknown_name); 00401 00402 // Return the square root of the first component of the 00403 // computed error. 00404 return std::sqrt(error_vals[0]); 00405 } 00406 00407 00408 00409 00410 00411 00412 00413 Real ExactSolution::l1_error(const std::string& sys_name, 00414 const std::string& unknown_name) 00415 { 00416 00417 // Check the inputs for validity, and get a reference 00418 // to the proper location to store the error 00419 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00420 unknown_name); 00421 00422 // Return the square root of the first component of the 00423 // computed error. 00424 return error_vals[3]; 00425 } 00426 00427 00428 00429 00430 00431 00432 00433 Real ExactSolution::l_inf_error(const std::string& sys_name, 00434 const std::string& unknown_name) 00435 { 00436 00437 // Check the inputs for validity, and get a reference 00438 // to the proper location to store the error 00439 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00440 unknown_name); 00441 00442 // Return the square root of the first component of the 00443 // computed error. 00444 return error_vals[4]; 00445 } 00446 00447 00448 00449 00450 00451 00452 00453 Real ExactSolution::h1_error(const std::string& sys_name, 00454 const std::string& unknown_name) 00455 { 00456 // If the user has supplied no exact derivative function, we 00457 // just integrate the H1 norm of the solution; i.e. its 00458 // difference from an "exact solution" of zero. 00459 00460 // Check the inputs for validity, and get a reference 00461 // to the proper location to store the error 00462 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00463 unknown_name); 00464 00465 // Return the square root of the sum of the computed errors. 00466 return std::sqrt(error_vals[0] + error_vals[1]); 00467 } 00468 00469 00470 Real ExactSolution::hcurl_error(const std::string& sys_name, 00471 const std::string& unknown_name) 00472 { 00473 return this->error_norm(sys_name,unknown_name,HCURL); 00474 } 00475 00476 00477 Real ExactSolution::hdiv_error(const std::string& sys_name, 00478 const std::string& unknown_name) 00479 { 00480 return this->error_norm(sys_name,unknown_name,HDIV); 00481 } 00482 00483 00484 00485 Real ExactSolution::h2_error(const std::string& sys_name, 00486 const std::string& unknown_name) 00487 { 00488 // If the user has supplied no exact derivative functions, we 00489 // just integrate the H2 norm of the solution; i.e. its 00490 // difference from an "exact solution" of zero. 00491 00492 // Check the inputs for validity, and get a reference 00493 // to the proper location to store the error 00494 std::vector<Real>& error_vals = this->_check_inputs(sys_name, 00495 unknown_name); 00496 00497 // Return the square root of the sum of the computed errors. 00498 return std::sqrt(error_vals[0] + error_vals[1] + error_vals[2]); 00499 } 00500 00501 00502 00503 00504 00505 00506 00507 template< typename OutputShape> 00508 void ExactSolution::_compute_error(const std::string& sys_name, 00509 const std::string& unknown_name, 00510 std::vector<Real>& error_vals) 00511 { 00512 // Make sure we aren't "overconfigured" 00513 libmesh_assert (!(_exact_values.size() && _equation_systems_fine)); 00514 00515 // We need a commmunicator. 00516 const Parallel::Communicator &communicator(_equation_systems.comm()); 00517 00518 // This function must be run on all processors at once 00519 libmesh_parallel_only(communicator); 00520 00521 // Get a reference to the system whose error is being computed. 00522 // If we have a fine grid, however, we'll integrate on that instead 00523 // for more accuracy. 00524 const System& computed_system = _equation_systems_fine ? 00525 _equation_systems_fine->get_system(sys_name) : 00526 _equation_systems.get_system (sys_name); 00527 00528 const Real time = _equation_systems.get_system(sys_name).time; 00529 00530 const unsigned int sys_num = computed_system.number(); 00531 const unsigned int var = computed_system.variable_number(unknown_name); 00532 const unsigned int var_component = 00533 computed_system.variable_scalar_number(var, 0); 00534 00535 // Prepare a global solution and a MeshFunction of the coarse system if we need one 00536 UniquePtr<MeshFunction> coarse_values; 00537 UniquePtr<NumericVector<Number> > comparison_soln = NumericVector<Number>::build(_equation_systems.comm()); 00538 if (_equation_systems_fine) 00539 { 00540 const System& comparison_system 00541 = _equation_systems.get_system(sys_name); 00542 00543 std::vector<Number> global_soln; 00544 comparison_system.update_global_solution(global_soln); 00545 comparison_soln->init(comparison_system.solution->size(), true, SERIAL); 00546 (*comparison_soln) = global_soln; 00547 00548 coarse_values = UniquePtr<MeshFunction> 00549 (new MeshFunction(_equation_systems, 00550 *comparison_soln, 00551 comparison_system.get_dof_map(), 00552 comparison_system.variable_number(unknown_name))); 00553 coarse_values->init(); 00554 } 00555 00556 // Initialize any functors we're going to use 00557 for (unsigned int i=0; i != _exact_values.size(); ++i) 00558 if (_exact_values[i]) 00559 _exact_values[i]->init(); 00560 00561 for (unsigned int i=0; i != _exact_derivs.size(); ++i) 00562 if (_exact_derivs[i]) 00563 _exact_derivs[i]->init(); 00564 00565 for (unsigned int i=0; i != _exact_hessians.size(); ++i) 00566 if (_exact_hessians[i]) 00567 _exact_hessians[i]->init(); 00568 00569 // Get a reference to the dofmap and mesh for that system 00570 const DofMap& computed_dof_map = computed_system.get_dof_map(); 00571 00572 const MeshBase& _mesh = computed_system.get_mesh(); 00573 00574 // Grab which element dimensions are present in the mesh 00575 const std::set<unsigned char>& elem_dims = _mesh.elem_dimensions(); 00576 00577 // Zero the error before summation 00578 // 0 - sum of square of function error (L2) 00579 // 1 - sum of square of gradient error (H1 semi) 00580 // 2 - sum of square of Hessian error (H2 semi) 00581 // 3 - sum of sqrt(square of function error) (L1) 00582 // 4 - max of sqrt(square of function error) (Linfty) 00583 // 5 - sum of square of curl error (HCurl semi) 00584 // 6 - sum of square of div error (HDiv semi) 00585 error_vals = std::vector<Real>(7, 0.); 00586 00587 // Construct Quadrature rule based on default quadrature order 00588 const FEType& fe_type = computed_dof_map.variable_type(var); 00589 00590 unsigned int n_vec_dim = FEInterface::n_vec_dim( _mesh, fe_type ); 00591 00592 // FIXME: MeshFunction needs to be updated to support vector-valued 00593 // elements before we can use a reference solution. 00594 if( (n_vec_dim > 1) && _equation_systems_fine ) 00595 { 00596 libMesh::err << "Error calculation using reference solution not yet\n" 00597 << "supported for vector-valued elements." 00598 << std::endl; 00599 libmesh_not_implemented(); 00600 } 00601 00602 00603 // Allow space for dims 0-3, even if we don't use them all 00604 std::vector<FEGenericBase<OutputShape>*> fe_ptrs(4,NULL); 00605 std::vector<QBase*> q_rules(4,NULL); 00606 00607 // Prepare finite elements for each dimension present in the mesh 00608 for( std::set<unsigned char>::const_iterator d_it = elem_dims.begin(); 00609 d_it != elem_dims.end(); ++d_it ) 00610 { 00611 q_rules[*d_it] = 00612 fe_type.default_quadrature_rule (*d_it, _extra_order).release(); 00613 00614 // Construct finite element object 00615 00616 fe_ptrs[*d_it] = FEGenericBase<OutputShape>::build(*d_it, fe_type).release(); 00617 00618 // Attach quadrature rule to FE object 00619 fe_ptrs[*d_it]->attach_quadrature_rule (q_rules[*d_it]); 00620 } 00621 00622 // The global degree of freedom indices associated 00623 // with the local degrees of freedom. 00624 std::vector<dof_id_type> dof_indices; 00625 00626 00627 // 00628 // Begin the loop over the elements 00629 // 00630 // TODO: this ought to be threaded (and using subordinate 00631 // MeshFunction objects in each thread rather than a single 00632 // master) 00633 MeshBase::const_element_iterator el = _mesh.active_local_elements_begin(); 00634 const MeshBase::const_element_iterator end_el = _mesh.active_local_elements_end(); 00635 00636 for ( ; el != end_el; ++el) 00637 { 00638 // Store a pointer to the element we are currently 00639 // working on. This allows for nicer syntax later. 00640 const Elem* elem = *el; 00641 const unsigned int dim = elem->dim(); 00642 00643 const subdomain_id_type elem_subid = elem->subdomain_id(); 00644 00645 // If the variable is not active on this subdomain, don't bother 00646 if(!computed_system.variable(var).active_on_subdomain(elem_subid)) 00647 continue; 00648 00649 /* If the variable is active, then we're going to restrict the 00650 MeshFunction evaluations to the current element subdomain. 00651 This is for cases such as mixed dimension meshes where we want 00652 to restrict the calculation to one particular domain. */ 00653 std::set<subdomain_id_type> subdomain_id; 00654 subdomain_id.insert(elem_subid); 00655 00656 FEGenericBase<OutputShape>* fe = fe_ptrs[dim]; 00657 QBase* qrule = q_rules[dim]; 00658 libmesh_assert(fe); 00659 libmesh_assert(qrule); 00660 00661 // The Jacobian*weight at the quadrature points. 00662 const std::vector<Real>& JxW = fe->get_JxW(); 00663 00664 // The value of the shape functions at the quadrature points 00665 // i.e. phi(i) = phi_values[i][qp] 00666 const std::vector<std::vector<OutputShape> >& phi_values = fe->get_phi(); 00667 00668 // The value of the shape function gradients at the quadrature points 00669 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >& 00670 dphi_values = fe->get_dphi(); 00671 00672 // The value of the shape function curls at the quadrature points 00673 // Only computed for vector-valued elements 00674 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> >* curl_values = NULL; 00675 00676 // The value of the shape function divergences at the quadrature points 00677 // Only computed for vector-valued elements 00678 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence> >* div_values = NULL; 00679 00680 if( FEInterface::field_type(fe_type) == TYPE_VECTOR ) 00681 { 00682 curl_values = &fe->get_curl_phi(); 00683 div_values = &fe->get_div_phi(); 00684 } 00685 00686 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00687 // The value of the shape function second derivatives at the quadrature points 00688 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >& 00689 d2phi_values = fe->get_d2phi(); 00690 #endif 00691 00692 // The XYZ locations (in physical space) of the quadrature points 00693 const std::vector<Point>& q_point = fe->get_xyz(); 00694 00695 // reinitialize the element-specific data 00696 // for the current element 00697 fe->reinit (elem); 00698 00699 // Get the local to global degree of freedom maps 00700 computed_dof_map.dof_indices (elem, dof_indices, var); 00701 00702 // The number of quadrature points 00703 const unsigned int n_qp = qrule->n_points(); 00704 00705 // The number of shape functions 00706 const unsigned int n_sf = 00707 cast_int<unsigned int>(dof_indices.size()); 00708 00709 // 00710 // Begin the loop over the Quadrature points. 00711 // 00712 for (unsigned int qp=0; qp<n_qp; qp++) 00713 { 00714 // Real u_h = 0.; 00715 // RealGradient grad_u_h; 00716 00717 typename FEGenericBase<OutputShape>::OutputNumber u_h = 0.; 00718 00719 typename FEGenericBase<OutputShape>::OutputNumberGradient grad_u_h; 00720 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00721 typename FEGenericBase<OutputShape>::OutputNumberTensor grad2_u_h; 00722 #endif 00723 typename FEGenericBase<OutputShape>::OutputNumber curl_u_h = 0.0; 00724 typename FEGenericBase<OutputShape>::OutputNumberDivergence div_u_h = 0.0; 00725 00726 // Compute solution values at the current 00727 // quadrature point. This reqiures a sum 00728 // over all the shape functions evaluated 00729 // at the quadrature point. 00730 for (unsigned int i=0; i<n_sf; i++) 00731 { 00732 // Values from current solution. 00733 u_h += phi_values[i][qp]*computed_system.current_solution (dof_indices[i]); 00734 grad_u_h += dphi_values[i][qp]*computed_system.current_solution (dof_indices[i]); 00735 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00736 grad2_u_h += d2phi_values[i][qp]*computed_system.current_solution (dof_indices[i]); 00737 #endif 00738 if( FEInterface::field_type(fe_type) == TYPE_VECTOR ) 00739 { 00740 curl_u_h += (*curl_values)[i][qp]*computed_system.current_solution (dof_indices[i]); 00741 div_u_h += (*div_values)[i][qp]*computed_system.current_solution (dof_indices[i]); 00742 } 00743 } 00744 00745 // Compute the value of the error at this quadrature point 00746 typename FEGenericBase<OutputShape>::OutputNumber exact_val = 0; 00747 RawAccessor<typename FEGenericBase<OutputShape>::OutputNumber> exact_val_accessor( exact_val, dim ); 00748 if (_exact_values.size() > sys_num && _exact_values[sys_num]) 00749 { 00750 for( unsigned int c = 0; c < n_vec_dim; c++) 00751 exact_val_accessor(c) = 00752 _exact_values[sys_num]-> 00753 component(var_component+c, q_point[qp], time); 00754 } 00755 else if (_equation_systems_fine) 00756 { 00757 // FIXME: Needs to be updated for vector-valued elements 00758 DenseVector<Number> output(1); 00759 (*coarse_values)(q_point[qp],time,output,&subdomain_id); 00760 exact_val = output(0); 00761 } 00762 const typename FEGenericBase<OutputShape>::OutputNumber val_error = u_h - exact_val; 00763 00764 // Add the squares of the error to each contribution 00765 Real error_sq = TensorTools::norm_sq(val_error); 00766 error_vals[0] += JxW[qp]*error_sq; 00767 00768 Real norm = sqrt(error_sq); 00769 error_vals[3] += JxW[qp]*norm; 00770 00771 if(error_vals[4]<norm) { error_vals[4] = norm; } 00772 00773 // Compute the value of the error in the gradient at this 00774 // quadrature point 00775 typename FEGenericBase<OutputShape>::OutputNumberGradient exact_grad; 00776 RawAccessor<typename FEGenericBase<OutputShape>::OutputNumberGradient> exact_grad_accessor( exact_grad, _mesh.spatial_dimension() ); 00777 if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num]) 00778 { 00779 for( unsigned int c = 0; c < n_vec_dim; c++) 00780 for( unsigned int d = 0; d < _mesh.spatial_dimension(); d++ ) 00781 exact_grad_accessor(d + c*_mesh.spatial_dimension() ) = 00782 _exact_derivs[sys_num]-> 00783 component(var_component+c, q_point[qp], time)(d); 00784 } 00785 else if (_equation_systems_fine) 00786 { 00787 // FIXME: Needs to be updated for vector-valued elements 00788 std::vector<Gradient> output(1); 00789 coarse_values->gradient(q_point[qp],time,output,&subdomain_id); 00790 exact_grad = output[0]; 00791 } 00792 00793 const typename FEGenericBase<OutputShape>::OutputNumberGradient grad_error = grad_u_h - exact_grad; 00794 00795 error_vals[1] += JxW[qp]*grad_error.size_sq(); 00796 00797 00798 if( FEInterface::field_type(fe_type) == TYPE_VECTOR ) 00799 { 00800 // Compute the value of the error in the curl at this 00801 // quadrature point 00802 typename FEGenericBase<OutputShape>::OutputNumber exact_curl = 0.0; 00803 if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num]) 00804 { 00805 exact_curl = TensorTools::curl_from_grad( exact_grad ); 00806 } 00807 else if (_equation_systems_fine) 00808 { 00809 // FIXME: Need to implement curl for MeshFunction and support reference 00810 // solution for vector-valued elements 00811 } 00812 00813 const typename FEGenericBase<OutputShape>::OutputNumber curl_error = curl_u_h - exact_curl; 00814 00815 error_vals[5] += JxW[qp]*TensorTools::norm_sq(curl_error); 00816 00817 // Compute the value of the error in the divergence at this 00818 // quadrature point 00819 typename FEGenericBase<OutputShape>::OutputNumberDivergence exact_div = 0.0; 00820 if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num]) 00821 { 00822 exact_div = TensorTools::div_from_grad( exact_grad ); 00823 } 00824 else if (_equation_systems_fine) 00825 { 00826 // FIXME: Need to implement div for MeshFunction and support reference 00827 // solution for vector-valued elements 00828 } 00829 00830 const typename FEGenericBase<OutputShape>::OutputNumberDivergence div_error = div_u_h - exact_div; 00831 00832 error_vals[6] += JxW[qp]*TensorTools::norm_sq(div_error); 00833 } 00834 00835 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00836 // Compute the value of the error in the hessian at this 00837 // quadrature point 00838 typename FEGenericBase<OutputShape>::OutputNumberTensor exact_hess; 00839 RawAccessor<typename FEGenericBase<OutputShape>::OutputNumberTensor> exact_hess_accessor( exact_hess, dim ); 00840 if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num]) 00841 { 00842 //FIXME: This needs to be implemented to support rank 3 tensors 00843 // which can't happen until type_n_tensor is fully implemented 00844 // and a RawAccessor<TypeNTensor> is fully implemented 00845 if( FEInterface::field_type(fe_type) == TYPE_VECTOR ) 00846 libmesh_not_implemented(); 00847 00848 for( unsigned int c = 0; c < n_vec_dim; c++) 00849 for( unsigned int d = 0; d < dim; d++ ) 00850 for( unsigned int e =0; e < dim; e++ ) 00851 exact_hess_accessor(d + e*dim + c*dim*dim) = 00852 _exact_hessians[sys_num]-> 00853 component(var_component+c, q_point[qp], time)(d,e); 00854 } 00855 else if (_equation_systems_fine) 00856 { 00857 // FIXME: Needs to be updated for vector-valued elements 00858 std::vector<Tensor> output(1); 00859 coarse_values->hessian(q_point[qp],time,output,&subdomain_id); 00860 exact_hess = output[0]; 00861 } 00862 00863 const typename FEGenericBase<OutputShape>::OutputNumberTensor grad2_error = grad2_u_h - exact_hess; 00864 00865 // FIXME: PB: Is this what we want for rank 3 tensors? 00866 error_vals[2] += JxW[qp]*grad2_error.size_sq(); 00867 #endif 00868 00869 } // end qp loop 00870 } // end element loop 00871 00872 // Clean up the FE and QBase pointers we created 00873 for( std::set<unsigned char>::const_iterator d_it = elem_dims.begin(); 00874 d_it != elem_dims.end(); ++d_it ) 00875 { 00876 delete fe_ptrs[*d_it]; 00877 delete q_rules[*d_it]; 00878 } 00879 00880 // Add up the error values on all processors, except for the L-infty 00881 // norm, for which the maximum is computed. 00882 Real l_infty_norm = error_vals[4]; 00883 communicator.max(l_infty_norm); 00884 communicator.sum(error_vals); 00885 error_vals[4] = l_infty_norm; 00886 } 00887 00888 // Explicit instantiations of templated member functions 00889 template void ExactSolution::_compute_error<Real>(const std::string&, const std::string&, std::vector<Real>&); 00890 template void ExactSolution::_compute_error<RealGradient>(const std::string&, const std::string&, std::vector<Real>&); 00891 00892 } // namespace libMesh