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