weighted_patch_recovery_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 std::sqrt std::pow std::abs 00023 00024 00025 // Local Includes 00026 #include "libmesh/dense_matrix.h" 00027 #include "libmesh/dense_vector.h" 00028 #include "libmesh/dof_map.h" 00029 #include "libmesh/elem.h" 00030 #include "libmesh/error_vector.h" 00031 #include "libmesh/fe_base.h" 00032 #include "libmesh/fem_context.h" 00033 #include "libmesh/libmesh_logging.h" 00034 #include "libmesh/mesh_base.h" 00035 #include "libmesh/numeric_vector.h" 00036 #include "libmesh/quadrature_grid.h" 00037 #include "libmesh/system.h" 00038 #include "libmesh/tensor_value.h" 00039 #include "libmesh/threads.h" 00040 #include "libmesh/weighted_patch_recovery_error_estimator.h" 00041 00042 namespace libMesh 00043 { 00044 void WeightedPatchRecoveryErrorEstimator::estimate_error (const System& system, 00045 ErrorVector& error_per_cell, 00046 const NumericVector<Number>* solution_vector, 00047 bool) 00048 { 00049 START_LOG("estimate_error()", "WeightedPatchRecoveryErrorEstimator"); 00050 00051 // The current mesh 00052 const MeshBase& mesh = system.get_mesh(); 00053 00054 // Resize the error_per_cell vector to be 00055 // the number of elements, initialize it to 0. 00056 error_per_cell.resize (mesh.max_elem_id()); 00057 std::fill (error_per_cell.begin(), error_per_cell.end(), 0.); 00058 00059 // Prepare current_local_solution to localize a non-standard 00060 // solution vector if necessary 00061 if (solution_vector && solution_vector != system.solution.get()) 00062 { 00063 NumericVector<Number>* newsol = 00064 const_cast<NumericVector<Number>*>(solution_vector); 00065 System &sys = const_cast<System&>(system); 00066 newsol->swap(*sys.solution); 00067 sys.update(); 00068 } 00069 00070 //------------------------------------------------------------ 00071 // Iterate over all the active elements in the mesh 00072 // that live on this processor. 00073 Threads::parallel_for (ConstElemRange(mesh.active_local_elements_begin(), 00074 mesh.active_local_elements_end(), 00075 200), 00076 EstimateError(system, 00077 *this, 00078 error_per_cell) 00079 ); 00080 00081 // Each processor has now computed the error contributions 00082 // for its local elements, and error_per_cell contains 0 for all the 00083 // non-local elements. Summing the vector will provide the true 00084 // value for each element, local or remote 00085 this->reduce_error(error_per_cell); 00086 00087 // If we used a non-standard solution before, now is the time to fix 00088 // the current_local_solution 00089 if (solution_vector && solution_vector != system.solution.get()) 00090 { 00091 NumericVector<Number>* newsol = 00092 const_cast<NumericVector<Number>*>(solution_vector); 00093 System &sys = const_cast<System&>(system); 00094 newsol->swap(*sys.solution); 00095 sys.update(); 00096 } 00097 00098 STOP_LOG("estimate_error()", "WeightedPatchRecoveryErrorEstimator"); 00099 } 00100 00101 00102 00103 void WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(const ConstElemRange &range) const 00104 { 00105 // The current mesh 00106 const MeshBase& mesh = system.get_mesh(); 00107 00108 // The dimensionality of the mesh 00109 const unsigned int dim = mesh.mesh_dimension(); 00110 00111 // The number of variables in the system 00112 const unsigned int n_vars = system.n_vars(); 00113 00114 // The DofMap for this system 00115 const DofMap& dof_map = system.get_dof_map(); 00116 00117 //------------------------------------------------------------ 00118 // Iterate over all the elements in the range. 00119 for (ConstElemRange::const_iterator elem_it=range.begin(); elem_it!=range.end(); ++elem_it) 00120 { 00121 // elem is necessarily an active element on the local processor 00122 const Elem* elem = *elem_it; 00123 00124 // We'll need an index into the error vector 00125 const dof_id_type e_id=elem->id(); 00126 00127 // We are going to build a patch containing the current element 00128 // and its neighbors on the local processor 00129 Patch patch; 00130 00131 // If we are reusing patches and the current element 00132 // already has an estimate associated with it, move on the 00133 // next element 00134 if(this->error_estimator.patch_reuse && error_per_cell[e_id] != 0) 00135 continue; 00136 00137 // If we are not reusing patches or havent built one containing this element, we build one 00138 00139 // Use user specified patch size and growth strategy 00140 patch.build_around_element (elem, error_estimator.target_patch_size, 00141 error_estimator.patch_growth_strategy); 00142 00143 // Declare a new_error_per_cell vector to hold error estimates 00144 // from each element in this patch, or one estimate if we are 00145 // not reusing patches since we will only be computing error for 00146 // one cell 00147 std::vector<Real> new_error_per_cell(1, 0.); 00148 if(this->error_estimator.patch_reuse) 00149 new_error_per_cell.resize(patch.size(), 0.); 00150 00151 //------------------------------------------------------------ 00152 // Process each variable in the system using the current patch 00153 for (unsigned int var=0; var<n_vars; var++) 00154 { 00155 #ifndef DEBUG 00156 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00157 libmesh_assert (error_estimator.error_norm.type(var) == L2 || 00158 error_estimator.error_norm.type(var) == H1_SEMINORM || 00159 error_estimator.error_norm.type(var) == H2_SEMINORM || 00160 error_estimator.error_norm.type(var) == H1_X_SEMINORM || 00161 error_estimator.error_norm.type(var) == H1_Y_SEMINORM || 00162 error_estimator.error_norm.type(var) == H1_Z_SEMINORM || 00163 error_estimator.error_norm.type(var) == L_INF || 00164 error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00165 error_estimator.error_norm.type(var) == W2_INF_SEMINORM); 00166 #else 00167 libmesh_assert (error_estimator.error_norm.type(var) == L2 || 00168 error_estimator.error_norm.type(var) == L_INF || 00169 error_estimator.error_norm.type(var) == H1_SEMINORM || 00170 error_estimator.error_norm.type(var) == H1_X_SEMINORM || 00171 error_estimator.error_norm.type(var) == H1_Y_SEMINORM || 00172 error_estimator.error_norm.type(var) == H1_Z_SEMINORM || 00173 error_estimator.error_norm.type(var) == W1_INF_SEMINORM); 00174 #endif 00175 if (var > 0) 00176 // We can't mix L_inf and L_2 norms 00177 libmesh_assert (((error_estimator.error_norm.type(var) == L2 || 00178 error_estimator.error_norm.type(var) == H1_SEMINORM || 00179 error_estimator.error_norm.type(var) == H1_X_SEMINORM || 00180 error_estimator.error_norm.type(var) == H1_Y_SEMINORM || 00181 error_estimator.error_norm.type(var) == H1_Z_SEMINORM || 00182 error_estimator.error_norm.type(var) == H2_SEMINORM) && 00183 (error_estimator.error_norm.type(var-1) == L2 || 00184 error_estimator.error_norm.type(var-1) == H1_SEMINORM || 00185 error_estimator.error_norm.type(var-1) == H1_X_SEMINORM || 00186 error_estimator.error_norm.type(var-1) == H1_Y_SEMINORM || 00187 error_estimator.error_norm.type(var-1) == H1_Z_SEMINORM || 00188 error_estimator.error_norm.type(var-1) == H2_SEMINORM)) || 00189 ((error_estimator.error_norm.type(var) == L_INF || 00190 error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00191 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) && 00192 (error_estimator.error_norm.type(var-1) == L_INF || 00193 error_estimator.error_norm.type(var-1) == W1_INF_SEMINORM || 00194 error_estimator.error_norm.type(var-1) == W2_INF_SEMINORM))); 00195 #endif 00196 00197 // Possibly skip this variable 00198 if (error_estimator.error_norm.weight(var) == 0.0) continue; 00199 00200 // The type of finite element to use for this variable 00201 const FEType& fe_type = dof_map.variable_type (var); 00202 00203 const Order element_order = static_cast<Order> 00204 (fe_type.order + elem->p_level()); 00205 00206 // Finite element object for use in this patch 00207 AutoPtr<FEBase> fe (FEBase::build (dim, fe_type)); 00208 00209 // Build an appropriate Gaussian quadrature rule 00210 AutoPtr<QBase> qrule (fe_type.default_quadrature_rule(dim)); 00211 00212 // Tell the finite element about the quadrature rule. 00213 fe->attach_quadrature_rule (qrule.get()); 00214 00215 // Get Jacobian values, etc.. 00216 const std::vector<Real>& JxW = fe->get_JxW(); 00217 const std::vector<Point>& q_point = fe->get_xyz(); 00218 00219 // Get whatever phi/dphi/d2phi values we need. Avoid 00220 // getting them unless the requested norm is actually going 00221 // to use them. 00222 00223 const std::vector<std::vector<Real> > *phi = NULL; 00224 // If we're using phi to assert the correct dof_indices 00225 // vector size later, then we'll need to get_phi whether we 00226 // plan to use it or not. 00227 #ifdef NDEBUG 00228 if (error_estimator.error_norm.type(var) == L2 || 00229 error_estimator.error_norm.type(var) == L_INF) 00230 #endif 00231 phi = &(fe->get_phi()); 00232 00233 const std::vector<std::vector<RealGradient> > *dphi = NULL; 00234 if (error_estimator.error_norm.type(var) == H1_SEMINORM || 00235 error_estimator.error_norm.type(var) == H1_X_SEMINORM || 00236 error_estimator.error_norm.type(var) == H1_Y_SEMINORM || 00237 error_estimator.error_norm.type(var) == H1_Z_SEMINORM || 00238 error_estimator.error_norm.type(var) == W1_INF_SEMINORM) 00239 dphi = &(fe->get_dphi()); 00240 00241 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00242 const std::vector<std::vector<RealTensor> > *d2phi = NULL; 00243 if (error_estimator.error_norm.type(var) == H2_SEMINORM || 00244 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00245 d2phi = &(fe->get_d2phi()); 00246 #endif 00247 00248 // global DOF indices 00249 std::vector<dof_id_type> dof_indices; 00250 00251 // Compute the approprite size for the patch projection matrices 00252 // and vectors; 00253 unsigned int matsize = element_order + 1; 00254 if (dim > 1) 00255 { 00256 matsize *= (element_order + 2); 00257 matsize /= 2; 00258 } 00259 if (dim > 2) 00260 { 00261 matsize *= (element_order + 3); 00262 matsize /= 3; 00263 } 00264 00265 DenseMatrix<Number> Kp(matsize,matsize); 00266 DenseVector<Number> F, Fx, Fy, Fz, Fxy, Fxz, Fyz; 00267 DenseVector<Number> Pu_h, Pu_x_h, Pu_y_h, Pu_z_h, Pu_xy_h, Pu_xz_h, Pu_yz_h; 00268 if (error_estimator.error_norm.type(var) == L2 || 00269 error_estimator.error_norm.type(var) == L_INF) 00270 { 00271 F.resize(matsize); Pu_h.resize(matsize); 00272 } 00273 else if (error_estimator.error_norm.type(var) == H1_SEMINORM || 00274 error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00275 error_estimator.error_norm.type(var) == H2_SEMINORM || 00276 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00277 { 00278 Fx.resize(matsize); Pu_x_h.resize(matsize); // stores xx in W2 cases 00279 #if LIBMESH_DIM > 1 00280 Fy.resize(matsize); Pu_y_h.resize(matsize); // stores yy in W2 cases 00281 #endif 00282 #if LIBMESH_DIM > 2 00283 Fz.resize(matsize); Pu_z_h.resize(matsize); // stores zz in W2 cases 00284 #endif 00285 } 00286 else if (error_estimator.error_norm.type(var) == H1_X_SEMINORM) 00287 { 00288 Fx.resize(matsize); Pu_x_h.resize(matsize); // Only need to compute the x gradient for the x component seminorm 00289 } 00290 else if (error_estimator.error_norm.type(var) == H1_Y_SEMINORM) 00291 { 00292 libmesh_assert(LIBMESH_DIM > 1); 00293 Fy.resize(matsize); Pu_y_h.resize(matsize); // Only need to compute the y gradient for the y component seminorm 00294 } 00295 else if (error_estimator.error_norm.type(var) == H1_Z_SEMINORM) 00296 { 00297 libmesh_assert(LIBMESH_DIM > 2); 00298 Fz.resize(matsize); Pu_z_h.resize(matsize); // Only need to compute the z gradient for the z component seminorm 00299 } 00300 00301 #if LIBMESH_DIM > 1 00302 if (error_estimator.error_norm.type(var) == H2_SEMINORM || 00303 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00304 { 00305 Fxy.resize(matsize); Pu_xy_h.resize(matsize); 00306 #if LIBMESH_DIM > 2 00307 Fxz.resize(matsize); Pu_xz_h.resize(matsize); 00308 Fyz.resize(matsize); Pu_yz_h.resize(matsize); 00309 #endif 00310 } 00311 #endif 00312 00313 00314 //------------------------------------------------------ 00315 // Loop over each element in the patch and compute their 00316 // contribution to the patch gradient projection. 00317 Patch::const_iterator patch_it = patch.begin(); 00318 const Patch::const_iterator patch_end = patch.end(); 00319 00320 for (; patch_it != patch_end; ++patch_it) 00321 { 00322 // The pth element in the patch 00323 const Elem* e_p = *patch_it; 00324 00325 // Reinitialize the finite element data for this element 00326 fe->reinit (e_p); 00327 00328 // Get the global DOF indices for the current variable 00329 // in the current element 00330 dof_map.dof_indices (e_p, dof_indices, var); 00331 libmesh_assert (dof_indices.size() == phi->size()); 00332 00333 const unsigned int n_dofs = 00334 libmesh_cast_int<unsigned int>(dof_indices.size()); 00335 const unsigned int n_qp = qrule->n_points(); 00336 00337 // Compute the weighted projection components from this cell. 00338 // \int_{Omega_e} \psi_i \psi_j = \int_{Omega_e} w * du_h/dx_k \psi_i 00339 for (unsigned int qp=0; qp<n_qp; qp++) 00340 { 00341 // Construct the shape function values for the patch projection 00342 std::vector<Real> psi(specpoly(dim, element_order, q_point[qp], matsize)); 00343 00344 // Patch matrix contribution 00345 for (unsigned int i=0; i<Kp.m(); i++) 00346 for (unsigned int j=0; j<Kp.n(); j++) 00347 Kp(i,j) += JxW[qp]*psi[i]*psi[j]; 00348 00349 if (error_estimator.error_norm.type(var) == L2 || 00350 error_estimator.error_norm.type(var) == L_INF) 00351 { 00352 // Compute the solution on the current patch element 00353 // the quadrature point 00354 Number u_h = libMesh::zero; 00355 00356 for (unsigned int i=0; i<n_dofs; i++) 00357 u_h += (*phi)[i][qp]*system.current_solution (dof_indices[i]); 00358 00359 // Patch RHS contributions 00360 for (unsigned int i=0; i<psi.size(); i++) 00361 F(i) = JxW[qp]*u_h*psi[i]; 00362 00363 } 00364 else if (error_estimator.error_norm.type(var) == H1_SEMINORM || 00365 error_estimator.error_norm.type(var) == W1_INF_SEMINORM) 00366 { 00367 // Compute the gradient on the current patch element 00368 // at the quadrature point 00369 Gradient grad_u_h; 00370 00371 for (unsigned int i=0; i<n_dofs; i++) 00372 grad_u_h.add_scaled ((*dphi)[i][qp], 00373 system.current_solution(dof_indices[i])); 00374 00375 00376 00377 // Patch RHS contributions 00378 for (unsigned int i=0; i<psi.size(); i++) 00379 { 00380 Fx(i) += JxW[qp]*grad_u_h(0)*psi[i]; 00381 #if LIBMESH_DIM > 1 00382 Fy(i) += JxW[qp]*grad_u_h(1)*psi[i]; 00383 #endif 00384 #if LIBMESH_DIM > 2 00385 Fz(i) += JxW[qp]*grad_u_h(2)*psi[i]; 00386 #endif 00387 } 00388 } 00389 else if (error_estimator.error_norm.type(var) == H1_X_SEMINORM) 00390 { 00391 // Compute the gradient on the current patch element 00392 // at the quadrature point 00393 Gradient grad_u_h; 00394 00395 for (unsigned int i=0; i<n_dofs; i++) 00396 grad_u_h.add_scaled ((*dphi)[i][qp], 00397 system.current_solution(dof_indices[i])); 00398 00399 00400 00401 // Patch RHS contributions 00402 for (unsigned int i=0; i<psi.size(); i++) 00403 { 00404 Fx(i) += JxW[qp]*grad_u_h(0)*psi[i]; 00405 } 00406 } 00407 else if (error_estimator.error_norm.type(var) == H1_Y_SEMINORM) 00408 { 00409 // Compute the gradient on the current patch element 00410 // at the quadrature point 00411 Gradient grad_u_h; 00412 00413 for (unsigned int i=0; i<n_dofs; i++) 00414 grad_u_h.add_scaled ((*dphi)[i][qp], 00415 system.current_solution(dof_indices[i])); 00416 00417 00418 00419 // Patch RHS contributions 00420 for (unsigned int i=0; i<psi.size(); i++) 00421 { 00422 Fy(i) += JxW[qp]*grad_u_h(1)*psi[i]; 00423 } 00424 } 00425 else if (error_estimator.error_norm.type(var) == H1_Z_SEMINORM) 00426 { 00427 // Compute the gradient on the current patch element 00428 // at the quadrature point 00429 Gradient grad_u_h; 00430 00431 for (unsigned int i=0; i<n_dofs; i++) 00432 grad_u_h.add_scaled ((*dphi)[i][qp], 00433 system.current_solution(dof_indices[i])); 00434 00435 00436 00437 // Patch RHS contributions 00438 for (unsigned int i=0; i<psi.size(); i++) 00439 { 00440 Fz(i) += JxW[qp]*grad_u_h(2)*psi[i]; 00441 } 00442 } 00443 else if (error_estimator.error_norm.type(var) == H2_SEMINORM || 00444 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00445 { 00446 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00447 // Compute the hessian on the current patch element 00448 // at the quadrature point 00449 Tensor hess_u_h; 00450 00451 for (unsigned int i=0; i<n_dofs; i++) 00452 hess_u_h.add_scaled ((*d2phi)[i][qp], 00453 system.current_solution(dof_indices[i])); 00454 00455 00456 00457 // Patch RHS contributions 00458 for (unsigned int i=0; i<psi.size(); i++) 00459 { 00460 Fx(i) += JxW[qp]*hess_u_h(0,0)*psi[i]; 00461 #if LIBMESH_DIM > 1 00462 Fy(i) += JxW[qp]*hess_u_h(1,1)*psi[i]; 00463 Fxy(i) += JxW[qp]*hess_u_h(0,1)*psi[i]; 00464 #endif 00465 #if LIBMESH_DIM > 2 00466 Fz(i) += JxW[qp]*hess_u_h(2,2)*psi[i]; 00467 Fxz(i) += JxW[qp]*hess_u_h(0,2)*psi[i]; 00468 Fyz(i) += JxW[qp]*hess_u_h(1,2)*psi[i]; 00469 #endif 00470 } 00471 #else 00472 libMesh::err << "ERROR: --enable-second-derivatives is required\n" 00473 << " for _sobolev_order == 2!\n"; 00474 libmesh_error(); 00475 #endif 00476 } 00477 else 00478 libmesh_error(); 00479 } // end quadrature loop 00480 } // end patch loop 00481 00482 00483 00484 //-------------------------------------------------- 00485 // Now we have fully assembled the projection system 00486 // for this patch. Project the gradient components. 00487 // MAY NEED TO USE PARTIAL PIVOTING! 00488 if (error_estimator.error_norm.type(var) == L2 || 00489 error_estimator.error_norm.type(var) == L_INF) 00490 { 00491 Kp.lu_solve(F, Pu_h); 00492 } 00493 else if (error_estimator.error_norm.type(var) == H1_SEMINORM || 00494 error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00495 error_estimator.error_norm.type(var) == H2_SEMINORM || 00496 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00497 { 00498 Kp.lu_solve (Fx, Pu_x_h); 00499 #if LIBMESH_DIM > 1 00500 Kp.lu_solve (Fy, Pu_y_h); 00501 #endif 00502 #if LIBMESH_DIM > 2 00503 Kp.lu_solve (Fz, Pu_z_h); 00504 #endif 00505 } 00506 else if (error_estimator.error_norm.type(var) == H1_X_SEMINORM) 00507 { 00508 Kp.lu_solve (Fx, Pu_x_h); 00509 } 00510 else if (error_estimator.error_norm.type(var) == H1_Y_SEMINORM) 00511 { 00512 Kp.lu_solve (Fy, Pu_y_h); 00513 } 00514 else if (error_estimator.error_norm.type(var) == H1_Z_SEMINORM) 00515 { 00516 Kp.lu_solve (Fz, Pu_z_h); 00517 } 00518 00519 #if LIBMESH_DIM > 1 00520 if (error_estimator.error_norm.type(var) == H2_SEMINORM || 00521 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00522 { 00523 Kp.lu_solve(Fxy, Pu_xy_h); 00524 #if LIBMESH_DIM > 2 00525 Kp.lu_solve(Fxz, Pu_xz_h); 00526 Kp.lu_solve(Fyz, Pu_yz_h); 00527 #endif 00528 } 00529 #endif 00530 00531 // If we are reusing patches, reuse the current patch to loop 00532 // over all elements in the current patch, otherwise build a new 00533 // patch containing just the current element and loop over it 00534 // Note that C++ will not allow patch_re_end to be a const here 00535 Patch::const_iterator patch_re_it; 00536 Patch::const_iterator patch_re_end; 00537 00538 // Declare a new patch 00539 Patch patch_re; 00540 00541 if(this->error_estimator.patch_reuse) 00542 { 00543 // Just get the iterators from the current patch 00544 patch_re_it = patch.begin(); 00545 patch_re_end = patch.end(); 00546 } 00547 else 00548 { 00549 // Use a target patch size of just 0, this will contain 00550 // just the current element 00551 patch_re.build_around_element (elem, 0, 00552 error_estimator.patch_growth_strategy); 00553 00554 // Get the iterators from this newly constructed patch 00555 patch_re_it = patch_re.begin(); 00556 patch_re_end = patch_re.end(); 00557 } 00558 00559 // If we are reusing patches, loop over all the elements 00560 // in the current patch and develop an estimate 00561 // for all the elements by computing || w * (P u_h - u_h)|| or ||w *(P grad_u_h - grad_u_h)|| 00562 // or ||w * (P hess_u_h - hess_u_h)|| according to the requested 00563 // seminorm, otherwise just compute it for the current element 00564 00565 // Get an FEMContext for this system, this will help us in 00566 // obtaining the weights from the user code 00567 FEMContext femcontext(system); 00568 error_estimator.weight_functions[var]->init_context(femcontext); 00569 00570 // Loop over every element in the patch 00571 for (unsigned int e = 0 ; patch_re_it != patch_re_end; patch_re_it++, ++e) 00572 { 00573 // Build the Finite Element for the current element 00574 00575 // The pth element in the patch 00576 const Elem* e_p = *patch_re_it; 00577 00578 // We'll need an index into the error vector for this element 00579 const dof_id_type e_p_id = e_p->id(); 00580 00581 // Get a pointer to the element, we need it to initialize 00582 // the FEMContext 00583 Elem *e_p_cast = const_cast<Elem *>(*patch_re_it); 00584 00585 // Initialize the FEMContext 00586 femcontext.pre_fe_reinit(system, e_p_cast); 00587 00588 // We will update the new_error_per_cell vector with element_error if the 00589 // error_per_cell[e_p_id] entry is non-zero, otherwise update it 00590 // with 0. i.e. leave it unchanged 00591 00592 // No need to compute the estimate if we are reusing patches and already have one 00593 if (this->error_estimator.patch_reuse && error_per_cell[e_p_id] != 0.) 00594 continue; 00595 00596 // Reinitialize the finite element data for this element 00597 fe->reinit (e_p); 00598 00599 // Get the global DOF indices for the current variable 00600 // in the current element 00601 dof_map.dof_indices (e_p, dof_indices, var); 00602 libmesh_assert (dof_indices.size() == phi->size()); 00603 00604 // The number of dofs for this variable on this element 00605 const unsigned int n_dofs = 00606 libmesh_cast_int<unsigned int>(dof_indices.size()); 00607 00608 // Variable to hold the error on the current element 00609 Real element_error = 0; 00610 00611 const Order qorder = 00612 static_cast<Order>(fe_type.order + e_p->p_level()); 00613 00614 // A quadrature rule for this element 00615 QGrid samprule (dim, qorder); 00616 00617 if (error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00618 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00619 fe->attach_quadrature_rule (&samprule); 00620 00621 // The number of points we will sample over 00622 const unsigned int n_sp = 00623 libmesh_cast_int<unsigned int>(JxW.size()); 00624 00625 // Loop over every sample point for the current element 00626 for (unsigned int sp=0; sp<n_sp; sp++) 00627 { 00628 // Compute the solution at the current sample point 00629 00630 std::vector<Number> temperr(6,0.0); // x,y,z or xx,yy,zz,xy,xz,yz 00631 00632 if (error_estimator.error_norm.type(var) == L2 || 00633 error_estimator.error_norm.type(var) == L_INF) 00634 { 00635 // Compute the value at the current sample point 00636 Number u_h = libMesh::zero; 00637 00638 for (unsigned int i=0; i<n_dofs; i++) 00639 u_h += (*phi)[i][sp]*system.current_solution (dof_indices[i]); 00640 00641 // Compute the phi values at the current sample point 00642 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00643 for (unsigned int i=0; i<matsize; i++) 00644 { 00645 temperr[0] += psi[i]*Pu_h(i); 00646 } 00647 00648 temperr[0] -= u_h; 00649 } 00650 else if (error_estimator.error_norm.type(var) == H1_SEMINORM || 00651 error_estimator.error_norm.type(var) == W1_INF_SEMINORM) 00652 { 00653 // Compute the gradient at the current sample point 00654 Gradient grad_u_h; 00655 00656 for (unsigned int i=0; i<n_dofs; i++) 00657 grad_u_h.add_scaled ((*dphi)[i][sp], 00658 system.current_solution(dof_indices[i])); 00659 00660 // Compute the phi values at the current sample point 00661 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00662 00663 for (unsigned int i=0; i<matsize; i++) 00664 { 00665 temperr[0] += psi[i]*Pu_x_h(i); 00666 #if LIBMESH_DIM > 1 00667 temperr[1] += psi[i]*Pu_y_h(i); 00668 #endif 00669 #if LIBMESH_DIM > 2 00670 temperr[2] += psi[i]*Pu_z_h(i); 00671 #endif 00672 } 00673 temperr[0] -= grad_u_h(0); 00674 #if LIBMESH_DIM > 1 00675 temperr[1] -= grad_u_h(1); 00676 #endif 00677 #if LIBMESH_DIM > 2 00678 temperr[2] -= grad_u_h(2); 00679 #endif 00680 } 00681 else if (error_estimator.error_norm.type(var) == H1_X_SEMINORM) 00682 { 00683 // Compute the gradient at the current sample point 00684 Gradient grad_u_h; 00685 00686 for (unsigned int i=0; i<n_dofs; i++) 00687 grad_u_h.add_scaled ((*dphi)[i][sp], 00688 system.current_solution(dof_indices[i])); 00689 00690 // Compute the phi values at the current sample point 00691 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00692 for (unsigned int i=0; i<matsize; i++) 00693 { 00694 temperr[0] += psi[i]*Pu_x_h(i); 00695 } 00696 00697 temperr[0] -= grad_u_h(0); 00698 } 00699 else if (error_estimator.error_norm.type(var) == H1_Y_SEMINORM) 00700 { 00701 // Compute the gradient at the current sample point 00702 Gradient grad_u_h; 00703 00704 for (unsigned int i=0; i<n_dofs; i++) 00705 grad_u_h.add_scaled ((*dphi)[i][sp], 00706 system.current_solution(dof_indices[i])); 00707 00708 // Compute the phi values at the current sample point 00709 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00710 for (unsigned int i=0; i<matsize; i++) 00711 { 00712 temperr[1] += psi[i]*Pu_y_h(i); 00713 } 00714 00715 temperr[1] -= grad_u_h(1); 00716 } 00717 else if (error_estimator.error_norm.type(var) == H1_Z_SEMINORM) 00718 { 00719 // Compute the gradient at the current sample point 00720 Gradient grad_u_h; 00721 00722 for (unsigned int i=0; i<n_dofs; i++) 00723 grad_u_h.add_scaled ((*dphi)[i][sp], 00724 system.current_solution(dof_indices[i])); 00725 00726 // Compute the phi values at the current sample point 00727 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00728 for (unsigned int i=0; i<matsize; i++) 00729 { 00730 temperr[2] += psi[i]*Pu_z_h(i); 00731 } 00732 00733 temperr[2] -= grad_u_h(2); 00734 } 00735 else if (error_estimator.error_norm.type(var) == H2_SEMINORM || 00736 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00737 { 00738 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00739 // Compute the Hessian at the current sample point 00740 Tensor hess_u_h; 00741 00742 for (unsigned int i=0; i<n_dofs; i++) 00743 hess_u_h.add_scaled ((*d2phi)[i][sp], 00744 system.current_solution(dof_indices[i])); 00745 00746 // Compute the phi values at the current sample point 00747 std::vector<Real> psi(specpoly(dim, element_order, q_point[sp], matsize)); 00748 for (unsigned int i=0; i<matsize; i++) 00749 { 00750 temperr[0] += psi[i]*Pu_x_h(i); 00751 #if LIBMESH_DIM > 1 00752 temperr[1] += psi[i]*Pu_y_h(i); 00753 temperr[3] += psi[i]*Pu_xy_h(i); 00754 #endif 00755 #if LIBMESH_DIM > 2 00756 temperr[2] += psi[i]*Pu_z_h(i); 00757 temperr[4] += psi[i]*Pu_xz_h(i); 00758 temperr[5] += psi[i]*Pu_yz_h(i); 00759 #endif 00760 } 00761 00762 temperr[0] -= hess_u_h(0,0); 00763 #if LIBMESH_DIM > 1 00764 temperr[1] -= hess_u_h(1,1); 00765 temperr[3] -= hess_u_h(0,1); 00766 #endif 00767 #if LIBMESH_DIM > 2 00768 temperr[2] -= hess_u_h(2,2); 00769 temperr[4] -= hess_u_h(0,2); 00770 temperr[5] -= hess_u_h(1,2); 00771 #endif 00772 #else 00773 libMesh::err << "ERROR: --enable-second-derivatives is required\n" 00774 << " for _sobolev_order == 2!\n"; 00775 libmesh_error(); 00776 #endif 00777 } 00778 00779 // Get the weight from the user code 00780 Number weight = (*error_estimator.weight_functions[var])(femcontext, q_point[sp], system.time); 00781 00782 // Add up relevant terms. We can easily optimize the 00783 // LIBMESH_DIM < 3 cases a little bit with the exception 00784 // of the W2 cases 00785 00786 if (error_estimator.error_norm.type(var) == L_INF) 00787 element_error = std::max(element_error, std::abs(weight*temperr[0])); 00788 else if (error_estimator.error_norm.type(var) == W1_INF_SEMINORM) 00789 for (unsigned int i=0; i != LIBMESH_DIM; ++i) 00790 element_error = std::max(element_error, std::abs(weight*temperr[i])); 00791 else if (error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00792 for (unsigned int i=0; i != 6; ++i) 00793 element_error = std::max(element_error, std::abs(weight*temperr[i])); 00794 else if (error_estimator.error_norm.type(var) == L2) 00795 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[0]); 00796 else if (error_estimator.error_norm.type(var) == H1_SEMINORM) 00797 for (unsigned int i=0; i != LIBMESH_DIM; ++i) 00798 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[i]); 00799 else if (error_estimator.error_norm.type(var) == H1_X_SEMINORM) 00800 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[0]); 00801 else if (error_estimator.error_norm.type(var) == H1_Y_SEMINORM) 00802 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[1]); 00803 else if (error_estimator.error_norm.type(var) == H1_Z_SEMINORM) 00804 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[2]); 00805 else if (error_estimator.error_norm.type(var) == H2_SEMINORM) 00806 { 00807 for (unsigned int i=0; i != LIBMESH_DIM; ++i) 00808 element_error += JxW[sp]*TensorTools::norm_sq(weight*temperr[i]); 00809 // Off diagonal terms enter into the Hessian norm twice 00810 for (unsigned int i=3; i != 6; ++i) 00811 element_error += JxW[sp]*2*TensorTools::norm_sq(weight*temperr[i]); 00812 } 00813 00814 } // End loop over sample points 00815 00816 if (error_estimator.error_norm.type(var) == L_INF || 00817 error_estimator.error_norm.type(var) == W1_INF_SEMINORM || 00818 error_estimator.error_norm.type(var) == W2_INF_SEMINORM) 00819 new_error_per_cell[e] += error_estimator.error_norm.weight(var) * element_error; 00820 else if (error_estimator.error_norm.type(var) == L2 || 00821 error_estimator.error_norm.type(var) == H1_SEMINORM || 00822 error_estimator.error_norm.type(var) == H1_X_SEMINORM || 00823 error_estimator.error_norm.type(var) == H1_Y_SEMINORM || 00824 error_estimator.error_norm.type(var) == H1_Z_SEMINORM || 00825 error_estimator.error_norm.type(var) == H2_SEMINORM) 00826 new_error_per_cell[e] += error_estimator.error_norm.weight_sq(var) * element_error; 00827 else 00828 libmesh_error(); 00829 } // End (re) loop over patch elements 00830 00831 } // end variables loop 00832 00833 // Now that we have the contributions from each variable, 00834 // we have take square roots of the entries we 00835 // added to error_per_cell to get an error norm 00836 // If we are reusing patches, once again reuse the current patch to loop 00837 // over all elements in the current patch, otherwise build a new 00838 // patch containing just the current element and loop over it 00839 Patch::const_iterator patch_re_it; 00840 Patch::const_iterator patch_re_end; 00841 00842 // Build a new patch if necessary 00843 Patch current_elem_patch; 00844 00845 if(this->error_estimator.patch_reuse) 00846 { 00847 // Just get the iterators from the current patch 00848 patch_re_it = patch.begin(); 00849 patch_re_end = patch.end(); 00850 } 00851 else 00852 { 00853 // Use a target patch size of just 0, this will contain 00854 // just the current element. 00855 current_elem_patch.build_around_element (elem, 0, 00856 error_estimator.patch_growth_strategy); 00857 00858 // Get the iterators from this newly constructed patch 00859 patch_re_it = current_elem_patch.begin(); 00860 patch_re_end = current_elem_patch.end(); 00861 } 00862 00863 // Loop over every element in the patch we just constructed 00864 for (unsigned int i = 0 ; patch_re_it != patch_re_end; ++patch_re_it, ++i) 00865 { 00866 // The pth element in the patch 00867 const Elem* e_p = *patch_re_it; 00868 00869 // We'll need an index into the error vector 00870 const dof_id_type e_p_id = e_p->id(); 00871 00872 // Update the error_per_cell vector for this element 00873 if (error_estimator.error_norm.type(0) == L2 || 00874 error_estimator.error_norm.type(0) == H1_SEMINORM || 00875 error_estimator.error_norm.type(0) == H1_X_SEMINORM || 00876 error_estimator.error_norm.type(0) == H1_Y_SEMINORM || 00877 error_estimator.error_norm.type(0) == H1_Z_SEMINORM || 00878 error_estimator.error_norm.type(0) == H2_SEMINORM) 00879 { 00880 Threads::spin_mutex::scoped_lock acquire(Threads::spin_mtx); 00881 if (!error_per_cell[e_p_id]) 00882 error_per_cell[e_p_id] = static_cast<ErrorVectorReal> 00883 (std::sqrt(new_error_per_cell[i])); 00884 } 00885 else 00886 { 00887 libmesh_assert (error_estimator.error_norm.type(0) == L_INF || 00888 error_estimator.error_norm.type(0) == W1_INF_SEMINORM || 00889 error_estimator.error_norm.type(0) == W2_INF_SEMINORM); 00890 Threads::spin_mutex::scoped_lock acquire(Threads::spin_mtx); 00891 if (!error_per_cell[e_p_id]) 00892 error_per_cell[e_p_id] = static_cast<ErrorVectorReal> 00893 (new_error_per_cell[i]); 00894 } 00895 00896 } // End loop over every element in patch 00897 00898 00899 } // end element loop 00900 00901 } // End () operator definition 00902 00903 } // namespace libMesh 00904
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:49 UTC
Hosted By: