$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 #include <limits> // for std::numeric_limits::max 00020 #include <math.h> // for sqrt 00021 00022 00023 // Local Includes 00024 #include "libmesh/hp_coarsentest.h" 00025 #include "libmesh/dense_matrix.h" 00026 #include "libmesh/dense_vector.h" 00027 #include "libmesh/dof_map.h" 00028 #include "libmesh/fe_base.h" 00029 #include "libmesh/fe_interface.h" 00030 #include "libmesh/libmesh_logging.h" 00031 #include "libmesh/elem.h" 00032 #include "libmesh/error_vector.h" 00033 #include "libmesh/mesh_base.h" 00034 #include "libmesh/quadrature.h" 00035 #include "libmesh/system.h" 00036 #include "libmesh/tensor_value.h" 00037 00038 #ifdef LIBMESH_ENABLE_AMR 00039 00040 namespace libMesh 00041 { 00042 00043 //----------------------------------------------------------------- 00044 // HPCoarsenTest implementations 00045 00046 void HPCoarsenTest::add_projection(const System &system, 00047 const Elem *elem, 00048 unsigned int var) 00049 { 00050 // If we have children, we need to add their projections instead 00051 if (!elem->active()) 00052 { 00053 libmesh_assert(!elem->subactive()); 00054 for (unsigned int c = 0; c != elem->n_children(); ++c) 00055 this->add_projection(system, elem->child(c), var); 00056 return; 00057 } 00058 00059 // The DofMap for this system 00060 const DofMap& dof_map = system.get_dof_map(); 00061 00062 // The type of finite element to use for this variable 00063 const FEType& fe_type = dof_map.variable_type (var); 00064 00065 const FEContinuity cont = fe->get_continuity(); 00066 00067 fe->reinit(elem); 00068 00069 dof_map.dof_indices(elem, dof_indices, var); 00070 00071 const unsigned int n_dofs = 00072 cast_int<unsigned int>(dof_indices.size()); 00073 00074 FEInterface::inverse_map (system.get_mesh().mesh_dimension(), 00075 fe_type, coarse, *xyz_values, coarse_qpoints); 00076 00077 fe_coarse->reinit(coarse, &coarse_qpoints); 00078 00079 const unsigned int n_coarse_dofs = 00080 cast_int<unsigned int>(phi_coarse->size()); 00081 00082 if (Uc.size() == 0) 00083 { 00084 Ke.resize(n_coarse_dofs, n_coarse_dofs); 00085 Ke.zero(); 00086 Fe.resize(n_coarse_dofs); 00087 Fe.zero(); 00088 Uc.resize(n_coarse_dofs); 00089 Uc.zero(); 00090 } 00091 libmesh_assert_equal_to (Uc.size(), phi_coarse->size()); 00092 00093 // Loop over the quadrature points 00094 for (unsigned int qp=0; qp<qrule->n_points(); qp++) 00095 { 00096 // The solution value at the quadrature point 00097 Number val = libMesh::zero; 00098 Gradient grad; 00099 Tensor hess; 00100 00101 for (unsigned int i=0; i != n_dofs; i++) 00102 { 00103 dof_id_type dof_num = dof_indices[i]; 00104 val += (*phi)[i][qp] * 00105 system.current_solution(dof_num); 00106 if (cont == C_ZERO || cont == C_ONE) 00107 grad.add_scaled((*dphi)[i][qp],system.current_solution(dof_num)); 00108 // grad += (*dphi)[i][qp] * 00109 // system.current_solution(dof_num); 00110 if (cont == C_ONE) 00111 hess.add_scaled((*d2phi)[i][qp], system.current_solution(dof_num)); 00112 // hess += (*d2phi)[i][qp] * 00113 // system.current_solution(dof_num); 00114 } 00115 00116 // The projection matrix and vector 00117 for (unsigned int i=0; i != Fe.size(); ++i) 00118 { 00119 Fe(i) += (*JxW)[qp] * 00120 (*phi_coarse)[i][qp]*val; 00121 if (cont == C_ZERO || cont == C_ONE) 00122 Fe(i) += (*JxW)[qp] * 00123 (grad*(*dphi_coarse)[i][qp]); 00124 if (cont == C_ONE) 00125 Fe(i) += (*JxW)[qp] * 00126 hess.contract((*d2phi_coarse)[i][qp]); 00127 // Fe(i) += (*JxW)[qp] * 00128 // (*d2phi_coarse)[i][qp].contract(hess); 00129 00130 for (unsigned int j=0; j != Fe.size(); ++j) 00131 { 00132 Ke(i,j) += (*JxW)[qp] * 00133 (*phi_coarse)[i][qp]*(*phi_coarse)[j][qp]; 00134 if (cont == C_ZERO || cont == C_ONE) 00135 Ke(i,j) += (*JxW)[qp] * 00136 (*dphi_coarse)[i][qp]*(*dphi_coarse)[j][qp]; 00137 if (cont == C_ONE) 00138 Ke(i,j) += (*JxW)[qp] * 00139 ((*d2phi_coarse)[i][qp].contract((*d2phi_coarse)[j][qp])); 00140 } 00141 } 00142 } 00143 } 00144 00145 void HPCoarsenTest::select_refinement (System &system) 00146 { 00147 START_LOG("select_refinement()", "HPCoarsenTest"); 00148 00149 // The current mesh 00150 MeshBase& mesh = system.get_mesh(); 00151 00152 // The dimensionality of the mesh 00153 const unsigned int dim = mesh.mesh_dimension(); 00154 00155 // The number of variables in the system 00156 const unsigned int n_vars = system.n_vars(); 00157 00158 // The DofMap for this system 00159 const DofMap& dof_map = system.get_dof_map(); 00160 00161 // The system number (for doing bad hackery) 00162 const unsigned int sys_num = system.number(); 00163 00164 // Check for a valid component_scale 00165 if (!component_scale.empty()) 00166 { 00167 if (component_scale.size() != n_vars) 00168 libmesh_error_msg("ERROR: component_scale is the wrong size:\n" \ 00169 << " component_scale.size()=" \ 00170 << component_scale.size() \ 00171 << "\n n_vars=" \ 00172 << n_vars); 00173 } 00174 else 00175 { 00176 // No specified scaling. Scale all variables by one. 00177 component_scale.resize (n_vars, 1.0); 00178 } 00179 00180 // Resize the error_per_cell vectors to handle 00181 // the number of elements, initialize them to 0. 00182 std::vector<ErrorVectorReal> h_error_per_cell(mesh.n_elem(), 0.); 00183 std::vector<ErrorVectorReal> p_error_per_cell(mesh.n_elem(), 0.); 00184 00185 // Loop over all the variables in the system 00186 for (unsigned int var=0; var<n_vars; var++) 00187 { 00188 // Possibly skip this variable 00189 if (!component_scale.empty()) 00190 if (component_scale[var] == 0.0) continue; 00191 00192 // The type of finite element to use for this variable 00193 const FEType& fe_type = dof_map.variable_type (var); 00194 00195 // Finite element objects for a fine (and probably a coarse) 00196 // element will be needed 00197 fe = FEBase::build (dim, fe_type); 00198 fe_coarse = FEBase::build (dim, fe_type); 00199 00200 // Any cached coarse element results have expired 00201 coarse = NULL; 00202 unsigned int cached_coarse_p_level = 0; 00203 00204 const FEContinuity cont = fe->get_continuity(); 00205 libmesh_assert (cont == DISCONTINUOUS || cont == C_ZERO || 00206 cont == C_ONE); 00207 00208 // Build an appropriate quadrature rule 00209 qrule = fe_type.default_quadrature_rule(dim); 00210 00211 // Tell the refined finite element about the quadrature 00212 // rule. The coarse finite element need not know about it 00213 fe->attach_quadrature_rule (qrule.get()); 00214 00215 // We will always do the integration 00216 // on the fine elements. Get their Jacobian values, etc.. 00217 JxW = &(fe->get_JxW()); 00218 xyz_values = &(fe->get_xyz()); 00219 00220 // The shape functions 00221 phi = &(fe->get_phi()); 00222 phi_coarse = &(fe_coarse->get_phi()); 00223 00224 // The shape function derivatives 00225 if (cont == C_ZERO || cont == C_ONE) 00226 { 00227 dphi = &(fe->get_dphi()); 00228 dphi_coarse = &(fe_coarse->get_dphi()); 00229 } 00230 00231 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00232 // The shape function second derivatives 00233 if (cont == C_ONE) 00234 { 00235 d2phi = &(fe->get_d2phi()); 00236 d2phi_coarse = &(fe_coarse->get_d2phi()); 00237 } 00238 #endif // defined (LIBMESH_ENABLE_SECOND_DERIVATIVES) 00239 00240 // Iterate over all the active elements in the mesh 00241 // that live on this processor. 00242 00243 MeshBase::const_element_iterator elem_it = 00244 mesh.active_local_elements_begin(); 00245 const MeshBase::const_element_iterator elem_end = 00246 mesh.active_local_elements_end(); 00247 00248 for (; elem_it != elem_end; ++elem_it) 00249 { 00250 const Elem* elem = *elem_it; 00251 00252 // We're only checking elements that are already flagged for h 00253 // refinement 00254 if (elem->refinement_flag() != Elem::REFINE) 00255 continue; 00256 00257 const dof_id_type e_id = elem->id(); 00258 00259 // Find the projection onto the parent element, 00260 // if necessary 00261 if (elem->parent() && 00262 (coarse != elem->parent() || 00263 cached_coarse_p_level != elem->p_level())) 00264 { 00265 Uc.resize(0); 00266 00267 coarse = elem->parent(); 00268 cached_coarse_p_level = elem->p_level(); 00269 00270 unsigned int old_parent_level = coarse->p_level(); 00271 (const_cast<Elem *>(coarse))->hack_p_level(elem->p_level()); 00272 00273 this->add_projection(system, coarse, var); 00274 00275 (const_cast<Elem *>(coarse))->hack_p_level(old_parent_level); 00276 00277 // Solve the h-coarsening projection problem 00278 Ke.cholesky_solve(Fe, Uc); 00279 } 00280 00281 fe->reinit(elem); 00282 00283 // Get the DOF indices for the fine element 00284 dof_map.dof_indices (elem, dof_indices, var); 00285 00286 // The number of quadrature points 00287 const unsigned int n_qp = qrule->n_points(); 00288 00289 // The number of DOFS on the fine element 00290 const unsigned int n_dofs = 00291 cast_int<unsigned int>(dof_indices.size()); 00292 00293 // The number of nodes on the fine element 00294 const unsigned int n_nodes = elem->n_nodes(); 00295 00296 // The average element value (used as an ugly hack 00297 // when we have nothing p-coarsened to compare to) 00298 // Real average_val = 0.; 00299 Number average_val = 0.; 00300 00301 // Calculate this variable's contribution to the p 00302 // refinement error 00303 00304 if (elem->p_level() == 0) 00305 { 00306 unsigned int n_vertices = 0; 00307 for (unsigned int n = 0; n != n_nodes; ++n) 00308 if (elem->is_vertex(n)) 00309 { 00310 n_vertices++; 00311 const Node * const node = elem->get_node(n); 00312 average_val += system.current_solution 00313 (node->dof_number(sys_num,var,0)); 00314 } 00315 average_val /= n_vertices; 00316 } 00317 else 00318 { 00319 unsigned int old_elem_level = elem->p_level(); 00320 (const_cast<Elem *>(elem))->hack_p_level(old_elem_level - 1); 00321 00322 fe_coarse->reinit(elem, &(qrule->get_points())); 00323 00324 const unsigned int n_coarse_dofs = 00325 cast_int<unsigned int>(phi_coarse->size()); 00326 00327 (const_cast<Elem *>(elem))->hack_p_level(old_elem_level); 00328 00329 Ke.resize(n_coarse_dofs, n_coarse_dofs); 00330 Ke.zero(); 00331 Fe.resize(n_coarse_dofs); 00332 Fe.zero(); 00333 00334 // Loop over the quadrature points 00335 for (unsigned int qp=0; qp<qrule->n_points(); qp++) 00336 { 00337 // The solution value at the quadrature point 00338 Number val = libMesh::zero; 00339 Gradient grad; 00340 Tensor hess; 00341 00342 for (unsigned int i=0; i != n_dofs; i++) 00343 { 00344 dof_id_type dof_num = dof_indices[i]; 00345 val += (*phi)[i][qp] * 00346 system.current_solution(dof_num); 00347 if (cont == C_ZERO || cont == C_ONE) 00348 grad.add_scaled((*dphi)[i][qp], system.current_solution(dof_num)); 00349 // grad += (*dphi)[i][qp] * 00350 // system.current_solution(dof_num); 00351 if (cont == C_ONE) 00352 hess.add_scaled((*d2phi)[i][qp], system.current_solution(dof_num)); 00353 // hess += (*d2phi)[i][qp] * 00354 // system.current_solution(dof_num); 00355 } 00356 00357 // The projection matrix and vector 00358 for (unsigned int i=0; i != Fe.size(); ++i) 00359 { 00360 Fe(i) += (*JxW)[qp] * 00361 (*phi_coarse)[i][qp]*val; 00362 if (cont == C_ZERO || cont == C_ONE) 00363 Fe(i) += (*JxW)[qp] * 00364 grad * (*dphi_coarse)[i][qp]; 00365 if (cont == C_ONE) 00366 Fe(i) += (*JxW)[qp] * 00367 hess.contract((*d2phi_coarse)[i][qp]); 00368 00369 for (unsigned int j=0; j != Fe.size(); ++j) 00370 { 00371 Ke(i,j) += (*JxW)[qp] * 00372 (*phi_coarse)[i][qp]*(*phi_coarse)[j][qp]; 00373 if (cont == C_ZERO || cont == C_ONE) 00374 Ke(i,j) += (*JxW)[qp] * 00375 (*dphi_coarse)[i][qp]*(*dphi_coarse)[j][qp]; 00376 if (cont == C_ONE) 00377 Ke(i,j) += (*JxW)[qp] * 00378 ((*d2phi_coarse)[i][qp].contract((*d2phi_coarse)[j][qp])); 00379 } 00380 } 00381 } 00382 00383 // Solve the p-coarsening projection problem 00384 Ke.cholesky_solve(Fe, Up); 00385 } 00386 00387 // loop over the integration points on the fine element 00388 for (unsigned int qp=0; qp<n_qp; qp++) 00389 { 00390 Number value_error = 0.; 00391 Gradient grad_error; 00392 Tensor hessian_error; 00393 for (unsigned int i=0; i<n_dofs; i++) 00394 { 00395 const dof_id_type dof_num = dof_indices[i]; 00396 value_error += (*phi)[i][qp] * 00397 system.current_solution(dof_num); 00398 if (cont == C_ZERO || cont == C_ONE) 00399 grad_error.add_scaled((*dphi)[i][qp], system.current_solution(dof_num)); 00400 // grad_error += (*dphi)[i][qp] * 00401 // system.current_solution(dof_num); 00402 if (cont == C_ONE) 00403 hessian_error.add_scaled((*d2phi)[i][qp], system.current_solution(dof_num)); 00404 // hessian_error += (*d2phi)[i][qp] * 00405 // system.current_solution(dof_num); 00406 } 00407 if (elem->p_level() == 0) 00408 { 00409 value_error -= average_val; 00410 } 00411 else 00412 { 00413 for (unsigned int i=0; i<Up.size(); i++) 00414 { 00415 value_error -= (*phi_coarse)[i][qp] * Up(i); 00416 if (cont == C_ZERO || cont == C_ONE) 00417 grad_error.subtract_scaled((*dphi_coarse)[i][qp], Up(i)); 00418 // grad_error -= (*dphi_coarse)[i][qp] * Up(i); 00419 if (cont == C_ONE) 00420 hessian_error.subtract_scaled((*d2phi_coarse)[i][qp], Up(i)); 00421 // hessian_error -= (*d2phi_coarse)[i][qp] * Up(i); 00422 } 00423 } 00424 00425 p_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00426 (component_scale[var] * 00427 (*JxW)[qp] * TensorTools::norm_sq(value_error)); 00428 if (cont == C_ZERO || cont == C_ONE) 00429 p_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00430 (component_scale[var] * 00431 (*JxW)[qp] * grad_error.size_sq()); 00432 if (cont == C_ONE) 00433 p_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00434 (component_scale[var] * 00435 (*JxW)[qp] * hessian_error.size_sq()); 00436 } 00437 00438 // Calculate this variable's contribution to the h 00439 // refinement error 00440 00441 if (!elem->parent()) 00442 { 00443 // For now, we'll always start with an h refinement 00444 h_error_per_cell[e_id] = 00445 std::numeric_limits<ErrorVectorReal>::max() / 2; 00446 } 00447 else 00448 { 00449 FEInterface::inverse_map (dim, fe_type, coarse, 00450 *xyz_values, coarse_qpoints); 00451 00452 unsigned int old_parent_level = coarse->p_level(); 00453 (const_cast<Elem *>(coarse))->hack_p_level(elem->p_level()); 00454 00455 fe_coarse->reinit(coarse, &coarse_qpoints); 00456 00457 (const_cast<Elem *>(coarse))->hack_p_level(old_parent_level); 00458 00459 // The number of DOFS on the coarse element 00460 unsigned int n_coarse_dofs = 00461 cast_int<unsigned int>(phi_coarse->size()); 00462 00463 // Loop over the quadrature points 00464 for (unsigned int qp=0; qp<n_qp; qp++) 00465 { 00466 // The solution difference at the quadrature point 00467 Number value_error = libMesh::zero; 00468 Gradient grad_error; 00469 Tensor hessian_error; 00470 00471 for (unsigned int i=0; i != n_dofs; ++i) 00472 { 00473 const dof_id_type dof_num = dof_indices[i]; 00474 value_error += (*phi)[i][qp] * 00475 system.current_solution(dof_num); 00476 if (cont == C_ZERO || cont == C_ONE) 00477 grad_error.add_scaled((*dphi)[i][qp], system.current_solution(dof_num)); 00478 // grad_error += (*dphi)[i][qp] * 00479 // system.current_solution(dof_num); 00480 if (cont == C_ONE) 00481 hessian_error.add_scaled((*d2phi)[i][qp], system.current_solution(dof_num)); 00482 // hessian_error += (*d2phi)[i][qp] * 00483 // system.current_solution(dof_num); 00484 } 00485 00486 for (unsigned int i=0; i != n_coarse_dofs; ++i) 00487 { 00488 value_error -= (*phi_coarse)[i][qp] * Uc(i); 00489 if (cont == C_ZERO || cont == C_ONE) 00490 // grad_error -= (*dphi_coarse)[i][qp] * Uc(i); 00491 grad_error.subtract_scaled((*dphi_coarse)[i][qp], Uc(i)); 00492 if (cont == C_ONE) 00493 hessian_error.subtract_scaled((*d2phi_coarse)[i][qp], Uc(i)); 00494 // hessian_error -= (*d2phi_coarse)[i][qp] * Uc(i); 00495 } 00496 00497 h_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00498 (component_scale[var] * 00499 (*JxW)[qp] * TensorTools::norm_sq(value_error)); 00500 if (cont == C_ZERO || cont == C_ONE) 00501 h_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00502 (component_scale[var] * 00503 (*JxW)[qp] * grad_error.size_sq()); 00504 if (cont == C_ONE) 00505 h_error_per_cell[e_id] += static_cast<ErrorVectorReal> 00506 (component_scale[var] * 00507 (*JxW)[qp] * hessian_error.size_sq()); 00508 } 00509 00510 } 00511 } 00512 } 00513 00514 // Now that we've got our approximations for p_error and h_error, let's see 00515 // if we want to switch any h refinement flags to p refinement 00516 00517 // Iterate over all the active elements in the mesh 00518 // that live on this processor. 00519 00520 MeshBase::element_iterator elem_it = 00521 mesh.active_local_elements_begin(); 00522 const MeshBase::element_iterator elem_end = 00523 mesh.active_local_elements_end(); 00524 00525 for (; elem_it != elem_end; ++elem_it) 00526 { 00527 Elem* elem = *elem_it; 00528 00529 // We're only checking elements that are already flagged for h 00530 // refinement 00531 if (elem->refinement_flag() != Elem::REFINE) 00532 continue; 00533 00534 const dof_id_type e_id = elem->id(); 00535 00536 unsigned int dofs_per_elem = 0, dofs_per_p_elem = 0; 00537 00538 // Loop over all the variables in the system 00539 for (unsigned int var=0; var<n_vars; var++) 00540 { 00541 // The type of finite element to use for this variable 00542 const FEType& fe_type = dof_map.variable_type (var); 00543 00544 // FIXME: we're overestimating the number of DOFs added by h 00545 // refinement 00546 FEType elem_fe_type = fe_type; 00547 elem_fe_type.order = 00548 static_cast<Order>(fe_type.order + elem->p_level()); 00549 dofs_per_elem += 00550 FEInterface::n_dofs(dim, elem_fe_type, elem->type()); 00551 00552 elem_fe_type.order = 00553 static_cast<Order>(fe_type.order + elem->p_level() + 1); 00554 dofs_per_p_elem += 00555 FEInterface::n_dofs(dim, elem_fe_type, elem->type()); 00556 } 00557 00558 const unsigned int new_h_dofs = dofs_per_elem * 00559 (elem->n_children() - 1); 00560 00561 const unsigned int new_p_dofs = dofs_per_p_elem - 00562 dofs_per_elem; 00563 00564 /* 00565 libMesh::err << "Cell " << e_id << ": h = " << elem->hmax() 00566 << ", p = " << elem->p_level() + 1 << "," << std::endl 00567 << " h_error = " << h_error_per_cell[e_id] 00568 << ", p_error = " << p_error_per_cell[e_id] << std::endl 00569 << " new_h_dofs = " << new_h_dofs 00570 << ", new_p_dofs = " << new_p_dofs << std::endl; 00571 */ 00572 const Real p_value = 00573 std::sqrt(p_error_per_cell[e_id]) * p_weight / new_p_dofs; 00574 const Real h_value = 00575 std::sqrt(h_error_per_cell[e_id]) / 00576 static_cast<Real>(new_h_dofs); 00577 if (p_value > h_value) 00578 { 00579 elem->set_p_refinement_flag(Elem::REFINE); 00580 elem->set_refinement_flag(Elem::DO_NOTHING); 00581 } 00582 } 00583 00584 STOP_LOG("select_refinement()", "HPCoarsenTest"); 00585 } 00586 00587 } // namespace libMesh 00588 00589 #endif // #ifdef LIBMESH_ENABLE_AMR