uniform_refinement_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 // C++ includes 00019 #include <algorithm> // for std::fill 00020 #include <sstream> 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/dof_map.h" 00027 #include "libmesh/elem.h" 00028 #include "libmesh/equation_systems.h" 00029 #include "libmesh/error_vector.h" 00030 #include "libmesh/fe.h" 00031 #include "libmesh/libmesh_common.h" 00032 #include "libmesh/libmesh_logging.h" 00033 #include "libmesh/mesh_base.h" 00034 #include "libmesh/mesh_refinement.h" 00035 #include "libmesh/numeric_vector.h" 00036 #include "libmesh/quadrature.h" 00037 #include "libmesh/system.h" 00038 #include "libmesh/uniform_refinement_estimator.h" 00039 #include "libmesh/partitioner.h" 00040 #include "libmesh/tensor_tools.h" 00041 00042 #ifdef LIBMESH_ENABLE_AMR 00043 00044 namespace libMesh 00045 { 00046 00047 //----------------------------------------------------------------- 00048 // ErrorEstimator implementations 00049 void UniformRefinementEstimator::estimate_error (const System& _system, 00050 ErrorVector& error_per_cell, 00051 const NumericVector<Number>* solution_vector, 00052 bool estimate_parent_error) 00053 { 00054 START_LOG("estimate_error()", "UniformRefinementEstimator"); 00055 std::map<const System*, const NumericVector<Number>* > solution_vectors; 00056 solution_vectors[&_system] = solution_vector; 00057 this->_estimate_error (NULL, &_system, &error_per_cell, NULL, NULL, 00058 &solution_vectors, estimate_parent_error); 00059 STOP_LOG("estimate_error()", "UniformRefinementEstimator"); 00060 } 00061 00062 void UniformRefinementEstimator::estimate_errors (const EquationSystems& _es, 00063 ErrorVector& error_per_cell, 00064 const std::map<const System*, SystemNorm>& error_norms, 00065 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00066 bool estimate_parent_error) 00067 { 00068 START_LOG("estimate_errors()", "UniformRefinementEstimator"); 00069 this->_estimate_error (&_es, NULL, &error_per_cell, NULL, 00070 &error_norms, solution_vectors, 00071 estimate_parent_error); 00072 STOP_LOG("estimate_errors()", "UniformRefinementEstimator"); 00073 } 00074 00075 void UniformRefinementEstimator::estimate_errors (const EquationSystems& _es, 00076 ErrorMap& errors_per_cell, 00077 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00078 bool estimate_parent_error) 00079 { 00080 START_LOG("estimate_errors()", "UniformRefinementEstimator"); 00081 this->_estimate_error (&_es, NULL, NULL, &errors_per_cell, NULL, 00082 solution_vectors, estimate_parent_error); 00083 STOP_LOG("estimate_errors()", "UniformRefinementEstimator"); 00084 } 00085 00086 void UniformRefinementEstimator::_estimate_error (const EquationSystems* _es, 00087 const System* _system, 00088 ErrorVector* error_per_cell, 00089 ErrorMap* errors_per_cell, 00090 const std::map<const System*, SystemNorm > *_error_norms, 00091 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00092 bool) 00093 { 00094 // Get a vector of the Systems we're going to work on, 00095 // and set up a error_norms map if necessary 00096 std::vector<System *> system_list; 00097 AutoPtr<std::map<const System*, SystemNorm > > error_norms = 00098 AutoPtr<std::map<const System*, SystemNorm > > 00099 (new std::map<const System*, SystemNorm>); 00100 00101 if (_es) 00102 { 00103 libmesh_assert(!_system); 00104 libmesh_assert(_es->n_systems()); 00105 _system = &(_es->get_system(0)); 00106 libmesh_assert_equal_to (&(_system->get_equation_systems()), _es); 00107 00108 libmesh_assert(_es->n_systems()); 00109 for (unsigned int i=0; i != _es->n_systems(); ++i) 00110 // We have to break the rules here, because we can't refine a const System 00111 system_list.push_back(const_cast<System *>(&(_es->get_system(i)))); 00112 00113 // If we're computing one vector, we need to know how to scale 00114 // each variable's contributions to it. 00115 if (_error_norms) 00116 { 00117 libmesh_assert(!errors_per_cell); 00118 } 00119 else 00120 // If we're computing many vectors, we just need to know which 00121 // variables to skip 00122 { 00123 libmesh_assert (errors_per_cell); 00124 00125 _error_norms = error_norms.get(); 00126 00127 for (unsigned int i=0; i!= _es->n_systems(); ++i) 00128 { 00129 const System &sys = _es->get_system(i); 00130 unsigned int n_vars = sys.n_vars(); 00131 00132 std::vector<Real> weights(n_vars, 0.0); 00133 for (unsigned int v = 0; v != n_vars; ++v) 00134 { 00135 if (errors_per_cell->find(std::make_pair(&sys, v)) == 00136 errors_per_cell->end()) 00137 continue; 00138 00139 weights[v] = 1.0; 00140 } 00141 (*error_norms)[&sys] = 00142 SystemNorm(std::vector<FEMNormType>(n_vars, error_norm.type(0)), 00143 weights); 00144 } 00145 } 00146 } 00147 else 00148 { 00149 libmesh_assert(_system); 00150 // We have to break the rules here, because we can't refine a const System 00151 system_list.push_back(const_cast<System *>(_system)); 00152 00153 libmesh_assert(!_error_norms); 00154 (*error_norms)[_system] = error_norm; 00155 _error_norms = error_norms.get(); 00156 } 00157 00158 // An EquationSystems reference will be convenient. 00159 // We have to break the rules here, because we can't refine a const System 00160 EquationSystems& es = 00161 const_cast<EquationSystems &>(_system->get_equation_systems()); 00162 00163 // The current mesh 00164 MeshBase& mesh = es.get_mesh(); 00165 00166 // The dimensionality of the mesh 00167 const unsigned int dim = mesh.mesh_dimension(); 00168 00169 // Resize the error_per_cell vectors to be 00170 // the number of elements, initialize them to 0. 00171 if (error_per_cell) 00172 { 00173 error_per_cell->clear(); 00174 error_per_cell->resize (mesh.max_elem_id(), 0.); 00175 } 00176 else 00177 { 00178 libmesh_assert(errors_per_cell); 00179 for (ErrorMap::iterator i = errors_per_cell->begin(); 00180 i != errors_per_cell->end(); ++i) 00181 { 00182 ErrorVector *e = i->second; 00183 e->clear(); 00184 e->resize(mesh.max_elem_id(), 0.); 00185 } 00186 } 00187 00188 // We'll want to back up all coarse grid vectors 00189 std::vector<std::map<std::string, NumericVector<Number> *> > 00190 coarse_vectors(system_list.size()); 00191 std::vector<NumericVector<Number> *> 00192 coarse_solutions(system_list.size()); 00193 std::vector<NumericVector<Number> *> 00194 coarse_local_solutions(system_list.size()); 00195 // And make copies of projected solutions 00196 std::vector<NumericVector<Number> *> 00197 projected_solutions(system_list.size()); 00198 00199 // And we'll need to temporarily change solution projection settings 00200 std::vector<bool> old_projection_settings(system_list.size()); 00201 00202 // And it'll be best to avoid any repartitioning 00203 AutoPtr<Partitioner> old_partitioner = mesh.partitioner(); 00204 mesh.partitioner().reset(NULL); 00205 00206 for (unsigned int i=0; i != system_list.size(); ++i) 00207 { 00208 System &system = *system_list[i]; 00209 00210 // Check for valid error_norms 00211 libmesh_assert (_error_norms->find(&system) != 00212 _error_norms->end()); 00213 00214 // Back up the solution vector 00215 coarse_solutions[i] = system.solution->clone().release(); 00216 coarse_local_solutions[i] = 00217 system.current_local_solution->clone().release(); 00218 00219 // Back up all other coarse grid vectors 00220 for (System::vectors_iterator vec = system.vectors_begin(); vec != 00221 system.vectors_end(); ++vec) 00222 { 00223 // The (string) name of this vector 00224 const std::string& var_name = vec->first; 00225 00226 coarse_vectors[i][var_name] = vec->second->clone().release(); 00227 } 00228 00229 // Use a non-standard solution vector if necessary 00230 if (solution_vectors && 00231 solution_vectors->find(&system) != solution_vectors->end() && 00232 solution_vectors->find(&system)->second && 00233 solution_vectors->find(&system)->second != system.solution.get()) 00234 { 00235 NumericVector<Number>* newsol = 00236 const_cast<NumericVector<Number>*> 00237 (solution_vectors->find(&system)->second); 00238 newsol->swap(*system.solution); 00239 system.update(); 00240 } 00241 00242 // Make sure the solution is projected when we refine the mesh 00243 old_projection_settings[i] = system.project_solution_on_reinit(); 00244 system.project_solution_on_reinit() = true; 00245 } 00246 00247 // Find the number of coarse mesh elements, to make it possible 00248 // to find correct coarse elem ids later 00249 const dof_id_type max_coarse_elem_id = mesh.max_elem_id(); 00250 #ifndef NDEBUG 00251 // n_coarse_elem is only used in an assertion later so 00252 // avoid declaring it unless asserts are active. 00253 const dof_id_type n_coarse_elem = mesh.n_elem(); 00254 #endif 00255 00256 // Uniformly refine the mesh 00257 MeshRefinement mesh_refinement(mesh); 00258 00259 libmesh_assert (number_h_refinements > 0 || number_p_refinements > 0); 00260 00261 // FIXME: this may break if there is more than one System 00262 // on this mesh but estimate_error was still called instead of 00263 // estimate_errors 00264 for (unsigned int i = 0; i != number_h_refinements; ++i) 00265 { 00266 mesh_refinement.uniformly_refine(1); 00267 es.reinit(); 00268 } 00269 00270 for (unsigned int i = 0; i != number_p_refinements; ++i) 00271 { 00272 mesh_refinement.uniformly_p_refine(1); 00273 es.reinit(); 00274 } 00275 00276 for (unsigned int i=0; i != system_list.size(); ++i) 00277 { 00278 System &system = *system_list[i]; 00279 00280 // Copy the projected coarse grid solutions, which will be 00281 // overwritten by solve() 00282 // projected_solutions[i] = system.solution->clone().release(); 00283 projected_solutions[i] = NumericVector<Number>::build().release(); 00284 projected_solutions[i]->init(system.solution->size(), true, SERIAL); 00285 system.solution->localize(*projected_solutions[i], 00286 system.get_dof_map().get_send_list()); 00287 } 00288 00289 // Are we doing a forward or an adjoint solve? 00290 bool solve_adjoint = false; 00291 if (solution_vectors) 00292 { 00293 System *sys = system_list[0]; 00294 libmesh_assert (solution_vectors->find(sys) != 00295 solution_vectors->end()); 00296 const NumericVector<Number> *vec = solution_vectors->find(sys)->second; 00297 for (unsigned int j=0; j != sys->qoi.size(); ++j) 00298 { 00299 std::ostringstream adjoint_name; 00300 adjoint_name << "adjoint_solution" << j; 00301 00302 if (vec == sys->request_vector(adjoint_name.str())) 00303 { 00304 solve_adjoint = true; 00305 break; 00306 } 00307 } 00308 } 00309 00310 // Get the uniformly refined solution. 00311 00312 if (_es) 00313 { 00314 // Even if we had a decent preconditioner, valid matrix etc. before 00315 // refinement, we don't any more. 00316 for (unsigned int i=0; i != es.n_systems(); ++i) 00317 es.get_system(i).disable_cache(); 00318 00319 // No specified vectors == forward solve 00320 if (!solution_vectors) 00321 es.solve(); 00322 else 00323 { 00324 libmesh_assert_equal_to (solution_vectors->size(), es.n_systems()); 00325 libmesh_assert (solution_vectors->find(system_list[0]) != 00326 solution_vectors->end()); 00327 libmesh_assert(solve_adjoint || 00328 (solution_vectors->find(system_list[0])->second == 00329 system_list[0]->solution.get()) || 00330 !solution_vectors->find(system_list[0])->second); 00331 00332 #ifdef DEBUG 00333 for (unsigned int i=0; i != system_list.size(); ++i) 00334 { 00335 System *sys = system_list[i]; 00336 libmesh_assert (solution_vectors->find(sys) != 00337 solution_vectors->end()); 00338 const NumericVector<Number> *vec = solution_vectors->find(sys)->second; 00339 if (solve_adjoint) 00340 { 00341 bool found_vec = false; 00342 for (unsigned int j=0; j != sys->qoi.size(); ++j) 00343 { 00344 std::ostringstream adjoint_name; 00345 adjoint_name << "adjoint_solution" << j; 00346 00347 if (vec == sys->request_vector(adjoint_name.str())) 00348 { 00349 found_vec = true; 00350 break; 00351 } 00352 } 00353 libmesh_assert(found_vec); 00354 } 00355 else 00356 libmesh_assert(vec == sys->solution.get() || !vec); 00357 } 00358 #endif 00359 00360 if (solve_adjoint) 00361 { 00362 std::vector<unsigned int> adjs(system_list.size(), 00363 libMesh::invalid_uint); 00364 // Set up proper initial guesses 00365 for (unsigned int i=0; i != system_list.size(); ++i) 00366 { 00367 System *sys = system_list[i]; 00368 libmesh_assert (solution_vectors->find(sys) != 00369 solution_vectors->end()); 00370 const NumericVector<Number> *vec = solution_vectors->find(sys)->second; 00371 for (unsigned int j=0; j != sys->qoi.size(); ++j) 00372 { 00373 std::ostringstream adjoint_name; 00374 adjoint_name << "adjoint_solution" << j; 00375 00376 if (vec == sys->request_vector(adjoint_name.str())) 00377 { 00378 adjs[i] = j; 00379 break; 00380 } 00381 } 00382 libmesh_assert_not_equal_to (adjs[i], libMesh::invalid_uint); 00383 system_list[i]->get_adjoint_solution(adjs[i]) = 00384 *system_list[i]->solution; 00385 } 00386 00387 es.adjoint_solve(); 00388 00389 // Put the adjoint_solution into solution for 00390 // comparisons 00391 for (unsigned int i=0; i != system_list.size(); ++i) 00392 { 00393 system_list[i]->get_adjoint_solution(adjs[i]).swap(*system_list[i]->solution); 00394 system_list[i]->update(); 00395 } 00396 } 00397 else 00398 es.solve(); 00399 } 00400 } 00401 else 00402 { 00403 System *sys = system_list[0]; 00404 00405 // Even if we had a decent preconditioner, valid matrix etc. before 00406 // refinement, we don't any more. 00407 sys->disable_cache(); 00408 00409 // No specified vectors == forward solve 00410 if (!solution_vectors) 00411 sys->solve(); 00412 else 00413 { 00414 libmesh_assert (solution_vectors->find(sys) != 00415 solution_vectors->end()); 00416 00417 const NumericVector<Number> *vec = solution_vectors->find(sys)->second; 00418 00419 libmesh_assert(solve_adjoint || 00420 (solution_vectors->find(sys)->second == 00421 sys->solution.get()) || 00422 !solution_vectors->find(sys)->second); 00423 00424 if (solve_adjoint) 00425 { 00426 unsigned int adj = libMesh::invalid_uint; 00427 for (unsigned int j=0; j != sys->qoi.size(); ++j) 00428 { 00429 std::ostringstream adjoint_name; 00430 adjoint_name << "adjoint_solution" << j; 00431 00432 if (vec == sys->request_vector(adjoint_name.str())) 00433 { 00434 adj = j; 00435 break; 00436 } 00437 } 00438 libmesh_assert_not_equal_to (adj, libMesh::invalid_uint); 00439 00440 // Set up proper initial guess 00441 sys->get_adjoint_solution(adj) = *sys->solution; 00442 sys->adjoint_solve(); 00443 // Put the adjoint_solution into solution for 00444 // comparisons 00445 sys->get_adjoint_solution(adj).swap(*sys->solution); 00446 sys->update(); 00447 } 00448 else 00449 sys->solve(); 00450 } 00451 } 00452 00453 // Get the error in the uniformly refined solution(s). 00454 00455 for (unsigned int sysnum=0; sysnum != system_list.size(); ++sysnum) 00456 { 00457 System &system = *system_list[sysnum]; 00458 00459 unsigned int n_vars = system.n_vars(); 00460 00461 DofMap &dof_map = system.get_dof_map(); 00462 00463 const SystemNorm &system_i_norm = 00464 _error_norms->find(&system)->second; 00465 00466 NumericVector<Number> *projected_solution = projected_solutions[sysnum]; 00467 00468 // Loop over all the variables in the system 00469 for (unsigned int var=0; var<n_vars; var++) 00470 { 00471 // Get the error vector to fill for this system and variable 00472 ErrorVector *err_vec = error_per_cell; 00473 if (!err_vec) 00474 { 00475 libmesh_assert(errors_per_cell); 00476 err_vec = 00477 (*errors_per_cell)[std::make_pair(&system,var)]; 00478 } 00479 00480 // The type of finite element to use for this variable 00481 const FEType& fe_type = dof_map.variable_type (var); 00482 00483 // Finite element object for each fine element 00484 AutoPtr<FEBase> fe (FEBase::build (dim, fe_type)); 00485 00486 // Build and attach an appropriate quadrature rule 00487 AutoPtr<QBase> qrule = fe_type.default_quadrature_rule(dim); 00488 fe->attach_quadrature_rule (qrule.get()); 00489 00490 const std::vector<Real>& JxW = fe->get_JxW(); 00491 const std::vector<std::vector<Real> >& phi = fe->get_phi(); 00492 const std::vector<std::vector<RealGradient> >& dphi = 00493 fe->get_dphi(); 00494 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00495 const std::vector<std::vector<RealTensor> >& d2phi = 00496 fe->get_d2phi(); 00497 #endif 00498 00499 // The global DOF indices for the fine element 00500 std::vector<dof_id_type> dof_indices; 00501 00502 // Iterate over all the active elements in the fine mesh 00503 // that live on this processor. 00504 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 00505 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 00506 00507 for (; elem_it != elem_end; ++elem_it) 00508 { 00509 // e is necessarily an active element on the local processor 00510 const Elem* elem = *elem_it; 00511 00512 // Find the element id for the corresponding coarse grid element 00513 const Elem* coarse = elem; 00514 dof_id_type e_id = coarse->id(); 00515 while (e_id >= max_coarse_elem_id) 00516 { 00517 libmesh_assert (coarse->parent()); 00518 coarse = coarse->parent(); 00519 e_id = coarse->id(); 00520 } 00521 00522 double L2normsq = 0., H1seminormsq = 0., H2seminormsq = 0.; 00523 00524 // reinitialize the element-specific data 00525 // for the current element 00526 fe->reinit (elem); 00527 00528 // Get the local to global degree of freedom maps 00529 dof_map.dof_indices (elem, dof_indices, var); 00530 00531 // The number of quadrature points 00532 const unsigned int n_qp = qrule->n_points(); 00533 00534 // The number of shape functions 00535 const unsigned int n_sf = 00536 libmesh_cast_int<unsigned int>(dof_indices.size()); 00537 00538 // 00539 // Begin the loop over the Quadrature points. 00540 // 00541 for (unsigned int qp=0; qp<n_qp; qp++) 00542 { 00543 Number u_fine = 0., u_coarse = 0.; 00544 00545 Gradient grad_u_fine, grad_u_coarse; 00546 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00547 Tensor grad2_u_fine, grad2_u_coarse; 00548 #endif 00549 00550 // Compute solution values at the current 00551 // quadrature point. This reqiures a sum 00552 // over all the shape functions evaluated 00553 // at the quadrature point. 00554 for (unsigned int i=0; i<n_sf; i++) 00555 { 00556 u_fine += phi[i][qp]*system.current_solution (dof_indices[i]); 00557 u_coarse += phi[i][qp]*(*projected_solution) (dof_indices[i]); 00558 grad_u_fine += dphi[i][qp]*system.current_solution (dof_indices[i]); 00559 grad_u_coarse += dphi[i][qp]*(*projected_solution) (dof_indices[i]); 00560 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00561 grad2_u_fine += d2phi[i][qp]*system.current_solution (dof_indices[i]); 00562 grad2_u_coarse += d2phi[i][qp]*(*projected_solution) (dof_indices[i]); 00563 #endif 00564 } 00565 00566 // Compute the value of the error at this quadrature point 00567 const Number val_error = u_fine - u_coarse; 00568 00569 // Add the squares of the error to each contribution 00570 if (system_i_norm.type(var) == L2 || 00571 system_i_norm.type(var) == H1 || 00572 system_i_norm.type(var) == H2) 00573 { 00574 L2normsq += JxW[qp] * system_i_norm.weight_sq(var) * 00575 TensorTools::norm_sq(val_error); 00576 libmesh_assert_greater_equal (L2normsq, 0.); 00577 } 00578 00579 00580 // Compute the value of the error in the gradient at this 00581 // quadrature point 00582 if (system_i_norm.type(var) == H1 || 00583 system_i_norm.type(var) == H2 || 00584 system_i_norm.type(var) == H1_SEMINORM) 00585 { 00586 Gradient grad_error = grad_u_fine - grad_u_coarse; 00587 00588 H1seminormsq += JxW[qp] * system_i_norm.weight_sq(var) * 00589 grad_error.size_sq(); 00590 libmesh_assert_greater_equal (H1seminormsq, 0.); 00591 } 00592 00593 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00594 // Compute the value of the error in the hessian at this 00595 // quadrature point 00596 if (system_i_norm.type(var) == H2 || 00597 system_i_norm.type(var) == H2_SEMINORM) 00598 { 00599 Tensor grad2_error = grad2_u_fine - grad2_u_coarse; 00600 00601 H2seminormsq += JxW[qp] * system_i_norm.weight_sq(var) * 00602 grad2_error.size_sq(); 00603 libmesh_assert_greater_equal (H2seminormsq, 0.); 00604 } 00605 #endif 00606 } // end qp loop 00607 00608 if (system_i_norm.type(var) == L2 || 00609 system_i_norm.type(var) == H1 || 00610 system_i_norm.type(var) == H2) 00611 (*err_vec)[e_id] += 00612 static_cast<ErrorVectorReal>(L2normsq); 00613 if (system_i_norm.type(var) == H1 || 00614 system_i_norm.type(var) == H2 || 00615 system_i_norm.type(var) == H1_SEMINORM) 00616 (*err_vec)[e_id] += 00617 static_cast<ErrorVectorReal>(H1seminormsq); 00618 00619 if (system_i_norm.type(var) == H2 || 00620 system_i_norm.type(var) == H2_SEMINORM) 00621 (*err_vec)[e_id] += 00622 static_cast<ErrorVectorReal>(H2seminormsq); 00623 } // End loop over active local elements 00624 } // End loop over variables 00625 00626 // Don't bother projecting the solution; we'll restore from backup 00627 // after coarsening 00628 system.project_solution_on_reinit() = false; 00629 } 00630 00631 00632 // Uniformly coarsen the mesh, without projecting the solution 00633 libmesh_assert (number_h_refinements > 0 || number_p_refinements > 0); 00634 00635 for (unsigned int i = 0; i != number_h_refinements; ++i) 00636 { 00637 mesh_refinement.uniformly_coarsen(1); 00638 // FIXME - should the reinits here be necessary? - RHS 00639 es.reinit(); 00640 } 00641 00642 for (unsigned int i = 0; i != number_p_refinements; ++i) 00643 { 00644 mesh_refinement.uniformly_p_coarsen(1); 00645 es.reinit(); 00646 } 00647 00648 // We should be back where we started 00649 libmesh_assert_equal_to (n_coarse_elem, mesh.n_elem()); 00650 00651 // Each processor has now computed the error contribuions 00652 // for its local elements. We need to sum the vector 00653 // and then take the square-root of each component. Note 00654 // that we only need to sum if we are running on multiple 00655 // processors, and we only need to take the square-root 00656 // if the value is nonzero. There will in general be many 00657 // zeros for the inactive elements. 00658 00659 if (error_per_cell) 00660 { 00661 // First sum the vector of estimated error values 00662 this->reduce_error(*error_per_cell); 00663 00664 // Compute the square-root of each component. 00665 START_LOG("std::sqrt()", "UniformRefinementEstimator"); 00666 for (unsigned int i=0; i<error_per_cell->size(); i++) 00667 if ((*error_per_cell)[i] != 0.) 00668 (*error_per_cell)[i] = std::sqrt((*error_per_cell)[i]); 00669 STOP_LOG("std::sqrt()", "UniformRefinementEstimator"); 00670 } 00671 else 00672 { 00673 for (ErrorMap::iterator it = errors_per_cell->begin(); 00674 it != errors_per_cell->end(); ++it) 00675 { 00676 ErrorVector *e = it->second; 00677 // First sum the vector of estimated error values 00678 this->reduce_error(*e); 00679 00680 // Compute the square-root of each component. 00681 START_LOG("std::sqrt()", "UniformRefinementEstimator"); 00682 for (unsigned int i=0; i<e->size(); i++) 00683 if ((*e)[i] != 0.) 00684 (*e)[i] = std::sqrt((*e)[i]); 00685 STOP_LOG("std::sqrt()", "UniformRefinementEstimator"); 00686 } 00687 } 00688 00689 // Restore old solutions and clean up the heap 00690 for (unsigned int i=0; i != system_list.size(); ++i) 00691 { 00692 System &system = *system_list[i]; 00693 00694 system.project_solution_on_reinit() = old_projection_settings[i]; 00695 00696 // Restore the coarse solution vectors and delete their copies 00697 *system.solution = *coarse_solutions[i]; 00698 delete coarse_solutions[i]; 00699 *system.current_local_solution = *coarse_local_solutions[i]; 00700 delete coarse_local_solutions[i]; 00701 delete projected_solutions[i]; 00702 00703 for (System::vectors_iterator vec = system.vectors_begin(); vec != 00704 system.vectors_end(); ++vec) 00705 { 00706 // The (string) name of this vector 00707 const std::string& var_name = vec->first; 00708 00709 system.get_vector(var_name) = *coarse_vectors[i][var_name]; 00710 00711 coarse_vectors[i][var_name]->clear(); 00712 delete coarse_vectors[i][var_name]; 00713 } 00714 } 00715 00716 // Restore old partitioner settings 00717 mesh.partitioner() = old_partitioner; 00718 } 00719 00720 } // namespace libMesh 00721 00722 #endif // #ifdef LIBMESH_ENABLE_AMR
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:49 UTC
Hosted By: