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