libMesh::AdjointRefinementEstimator Class Reference
#include <adjoint_refinement_estimator.h>

Public Types | |
| typedef std::map< std::pair < const System *, unsigned int > , ErrorVector * > | ErrorMap |
Public Member Functions | |
| AdjointRefinementEstimator () | |
| ~AdjointRefinementEstimator () | |
| QoISet & | qoi_set () |
| const QoISet & | qoi_set () const |
| virtual void | estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=NULL, bool estimate_parent_error=false) |
| Number & | get_global_QoI_error_estimate (unsigned int qoi_index) |
| virtual void | estimate_errors (const EquationSystems &equation_systems, ErrorVector &error_per_cell, const std::map< const System *, SystemNorm > &error_norms, const std::map< const System *, const NumericVector< Number > * > *solution_vectors=NULL, bool estimate_parent_error=false) |
| virtual void | estimate_errors (const EquationSystems &equation_systems, ErrorMap &errors_per_cell, const std::map< const System *, const NumericVector< Number > * > *solution_vectors=NULL, bool estimate_parent_error=false) |
Public Attributes | |
| unsigned char | number_h_refinements |
| unsigned char | number_p_refinements |
| SystemNorm | error_norm |
Protected Member Functions | |
| void | reduce_error (std::vector< float > &error_per_cell) const |
Protected Attributes | |
| std::vector< Number > | computed_global_QoI_errors |
| QoISet | _qoi_set |
Detailed Description
This class implements a ``brute force'' goal-oriented error estimator which computes an estimate of error in a quantity of interest based on the residual of the current coarse grid primal solution as weighted against an adjoint solution on a uniformly refined (in h and/or p, for an arbitrary number of levels) grid.
Definition at line 46 of file adjoint_refinement_estimator.h.
Member Typedef Documentation
typedef std::map<std::pair<const System*, unsigned int>, ErrorVector*> libMesh::ErrorEstimator::ErrorMap [inherited] |
When calculating many error vectors at once, we need a data structure to hold them all
Definition at line 107 of file error_estimator.h.
Constructor & Destructor Documentation
| libMesh::AdjointRefinementEstimator::AdjointRefinementEstimator | ( | ) | [inline] |
Constructor. Sets the most common default parameter values.
Definition at line 53 of file adjoint_refinement_estimator.h.
References libMesh::ErrorEstimator::error_norm, and libMeshEnums::INVALID_NORM.
00053 : number_h_refinements(1), 00054 number_p_refinements(0), 00055 _qoi_set(QoISet()) 00056 { 00057 // We're not actually going to use error_norm; our norms are 00058 // absolute values of QoI error. 00059 error_norm = INVALID_NORM; 00060 }
| libMesh::AdjointRefinementEstimator::~AdjointRefinementEstimator | ( | ) | [inline] |
Member Function Documentation
| void libMesh::AdjointRefinementEstimator::estimate_error | ( | const System & | system, | |
| ErrorVector & | error_per_cell, | |||
| const NumericVector< Number > * | solution_vector = NULL, |
|||
| bool | estimate_parent_error = false | |||
| ) | [virtual] |
This function does uniform refinements and an adjoint solve to get an adjoint solution on each cell, then estimates the error by finding the weighted residual of the coarse solution with the fine adjoint solution.
system.solve() and system.assembly() must be called, and so should have no side effects.
Only the provided system is solved on the refined mesh; we don't support adjoint solves on loosely coupled collections of Systems.
The estimated error is output in the vector error_per_cell
Implements libMesh::ErrorEstimator.
Definition at line 67 of file adjoint_refinement_estimator.C.
References _qoi_set, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::System::adjoint_solve(), libMesh::MeshBase::allow_renumbering(), libMesh::NumericVector< T >::clear(), libMesh::NumericVector< T >::close(), computed_global_QoI_errors, libMesh::System::current_local_solution, libMesh::DofMap::dof_indices(), libMesh::NumericVector< T >::dot(), libMesh::ErrorVectorReal, libMesh::Elem::find_point_neighbors(), libMesh::AutoPtr< Tp >::get(), libMesh::System::get_adjoint_solution(), libMesh::System::get_dof_map(), libMesh::System::get_equation_systems(), libMesh::EquationSystems::get_mesh(), libMesh::Elem::get_node(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMeshEnums::GHOSTED, libMesh::QoISet::has_index(), libMesh::DofObject::id(), libMesh::NumericVector< T >::init(), libMesh::libmesh_real(), libMesh::NumericVector< T >::localize(), libMesh::MeshBase::max_elem_id(), mesh, libMesh::System::n_dofs(), libMesh::MeshBase::n_elem(), libMesh::System::n_local_dofs(), libMesh::Elem::n_nodes(), number_h_refinements, number_p_refinements, libMesh::Elem::parent(), libMesh::MeshBase::partitioner(), libMesh::System::project_solution_on_reinit(), libMesh::System::qoi, libMesh::Real, libMesh::EquationSystems::reinit(), libMesh::AutoPtr< Tp >::release(), libMesh::AutoPtr< Tp >::reset(), libMeshEnums::SERIAL, libMesh::System::solution, libMesh::NumericVector< T >::swap(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::MeshRefinement::uniformly_p_coarsen(), libMesh::MeshRefinement::uniformly_p_refine(), libMesh::MeshRefinement::uniformly_refine(), libMesh::System::update(), libMesh::System::vectors_begin(), libMesh::System::vectors_end(), and libMesh::QoISet::weight().
00071 { 00072 // We have to break the rules here, because we can't refine a const System 00073 System& system = const_cast<System&>(_system); 00074 00075 // An EquationSystems reference will be convenient. 00076 EquationSystems& es = system.get_equation_systems(); 00077 00078 // The current mesh 00079 MeshBase& mesh = es.get_mesh(); 00080 00081 // Resize the error_per_cell vector to be 00082 // the number of elements, initialized to 0. 00083 error_per_cell.clear(); 00084 error_per_cell.resize (mesh.max_elem_id(), 0.); 00085 00086 // We'll want to back up all coarse grid vectors 00087 std::map<std::string, NumericVector<Number> *> coarse_vectors; 00088 for (System::vectors_iterator vec = system.vectors_begin(); vec != 00089 system.vectors_end(); ++vec) 00090 { 00091 // The (string) name of this vector 00092 const std::string& var_name = vec->first; 00093 00094 coarse_vectors[var_name] = vec->second->clone().release(); 00095 } 00096 // Back up the coarse solution and coarse local solution 00097 NumericVector<Number> * coarse_solution = 00098 system.solution->clone().release(); 00099 NumericVector<Number> * coarse_local_solution = 00100 system.current_local_solution->clone().release(); 00101 // And make copies of the projected solution 00102 NumericVector<Number> * projected_solution; 00103 00104 // And we'll need to temporarily change solution projection settings 00105 bool old_projection_setting; 00106 old_projection_setting = system.project_solution_on_reinit(); 00107 00108 // Make sure the solution is projected when we refine the mesh 00109 system.project_solution_on_reinit() = true; 00110 00111 // And it'll be best to avoid any repartitioning 00112 AutoPtr<Partitioner> old_partitioner = mesh.partitioner(); 00113 mesh.partitioner().reset(NULL); 00114 00115 // And we can't allow any renumbering 00116 const bool old_renumbering_setting = mesh.allow_renumbering(); 00117 mesh.allow_renumbering(false); 00118 00119 // Use a non-standard solution vector if necessary 00120 if (solution_vector && solution_vector != system.solution.get()) 00121 { 00122 NumericVector<Number> *newsol = 00123 const_cast<NumericVector<Number>*> (solution_vector); 00124 newsol->swap(*system.solution); 00125 system.update(); 00126 } 00127 00128 #ifndef NDEBUG 00129 // n_coarse_elem is only used in an assertion later so 00130 // avoid declaring it unless asserts are active. 00131 const dof_id_type n_coarse_elem = mesh.n_elem(); 00132 #endif 00133 00134 // Uniformly refine the mesh 00135 MeshRefinement mesh_refinement(mesh); 00136 00137 libmesh_assert (number_h_refinements > 0 || number_p_refinements > 0); 00138 00139 // FIXME: this may break if there is more than one System 00140 // on this mesh but estimate_error was still called instead of 00141 // estimate_errors 00142 for (unsigned int i = 0; i != number_h_refinements; ++i) 00143 { 00144 mesh_refinement.uniformly_refine(1); 00145 es.reinit(); 00146 } 00147 00148 for (unsigned int i = 0; i != number_p_refinements; ++i) 00149 { 00150 mesh_refinement.uniformly_p_refine(1); 00151 es.reinit(); 00152 } 00153 00154 // Copy the projected coarse grid solutions, which will be 00155 // overwritten by solve() 00156 projected_solution = NumericVector<Number>::build().release(); 00157 projected_solution->init(system.solution->size(), true, SERIAL); 00158 system.solution->localize(*projected_solution, 00159 system.get_dof_map().get_send_list()); 00160 00161 // Rebuild the rhs with the projected primal solution 00162 (dynamic_cast<ImplicitSystem&>(system)).assembly(true, false); 00163 NumericVector<Number> & projected_residual = (dynamic_cast<ExplicitSystem&>(system)).get_vector("RHS Vector"); 00164 projected_residual.close(); 00165 00166 // Solve the adjoint problem on the refined FE space 00167 system.adjoint_solve(); 00168 00169 // Now that we have the refined adjoint solution and the projected primal solution, 00170 // we first compute the global QoI error estimate 00171 00172 // Resize the computed_global_QoI_errors vector to hold the error estimates for each QoI 00173 computed_global_QoI_errors.resize(system.qoi.size()); 00174 00175 // Loop over all the adjoint solutions and get the QoI error 00176 // contributions from all of them 00177 for (unsigned int j=0; j != system.qoi.size(); j++) 00178 { 00179 computed_global_QoI_errors[j] = projected_residual.dot(system.get_adjoint_solution(j)); 00180 } 00181 00182 // Done with the global error estimates, now construct the element wise error indicators 00183 00184 // We ought to account for 'spill-over' effects while computing the 00185 // element error indicators This happens because the same dof is 00186 // shared by multiple elements, one way of mitigating this is to 00187 // scale the contribution from each dof by the number of elements it 00188 // belongs to We first obtain the number of elements each node 00189 // belongs to 00190 00191 bool split_shared_dofs = false; 00192 00193 if (split_shared_dofs) { 00194 00195 // A map that relates a node id to an int that will tell us how many elements it is a node of 00196 LIBMESH_BEST_UNORDERED_MAP<dof_id_type, unsigned int>shared_element_count; 00197 00198 // To fill this map, we will loop over elements, and then in each element, we will loop 00199 // over the nodes each element contains, and then query it for the number of coarse 00200 // grid elements it was a node of 00201 00202 // We will be iterating over all the active elements in the fine mesh that live on 00203 // this processor 00204 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 00205 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 00206 00207 // Keep track of which nodes we have already dealt with 00208 std::vector<dof_id_type> processed_node_ids; 00209 00210 // Start loop over elems 00211 for(; elem_it != elem_end; ++elem_it) 00212 { 00213 // Pointer to this element 00214 const Elem* elem = *elem_it; 00215 00216 // Loop over the nodes in the element 00217 for(unsigned int n=0; n != elem->n_nodes(); ++n) 00218 { 00219 // Get a pointer to the current node 00220 Node* node = elem->get_node(n); 00221 00222 // Get the id of this node 00223 dof_id_type node_id = node->id(); 00224 00225 // A processed node flag 00226 bool processed_node = false; 00227 00228 // Loop over existing processed nodes and see if 00229 // we have already done this one 00230 for(std::size_t i = 0; i != processed_node_ids.size(); i++) 00231 { 00232 if(node_id == processed_node_ids[i]) 00233 { 00234 processed_node = true; 00235 } 00236 } 00237 00238 // If we havent already processed this node, do so now 00239 if(!processed_node) 00240 { 00241 // Declare a neighbor_set to be filled by the find_point_neighbors 00242 std::set<const Elem *> fine_grid_neighbor_set; 00243 00244 // Call find_point_neighbors to fill the neighbor_set 00245 elem->find_point_neighbors(*node, fine_grid_neighbor_set); 00246 00247 // A vector to hold the coarse grid parents neighbors 00248 std::vector<dof_id_type> coarse_grid_neighbors; 00249 00250 // Iterators over the fine grid neighbors set 00251 std::set<const Elem*>::iterator fine_neighbor_it = fine_grid_neighbor_set.begin(); 00252 const std::set<const Elem*>::iterator fine_neighbor_end = fine_grid_neighbor_set.end(); 00253 00254 // Loop over all the fine neighbors of this node 00255 for(; fine_neighbor_it != fine_neighbor_end ; ++fine_neighbor_it) 00256 { 00257 // Pointer to the current fine neighbor element 00258 const Elem* fine_elem = *fine_neighbor_it; 00259 00260 // Find the element id for the corresponding coarse grid element 00261 const Elem* coarse_elem = fine_elem; 00262 for (unsigned int j = 0; j != number_h_refinements; ++j) 00263 { 00264 libmesh_assert (coarse_elem->parent()); 00265 00266 coarse_elem = coarse_elem->parent(); 00267 } 00268 00269 // Loop over the existing coarse neighbors and check if this one is 00270 // already in there 00271 const dof_id_type coarse_id = coarse_elem->id(); 00272 std::size_t j = 0; 00273 for (; j != coarse_grid_neighbors.size(); j++) 00274 { 00275 // If the set already contains this element break out of the loop 00276 if(coarse_grid_neighbors[j] == coarse_id) 00277 { 00278 break; 00279 } 00280 } 00281 00282 // If we didn't leave the loop even at the last element, 00283 // this is a new neighbour, put in the coarse_grid_neighbor_set 00284 if(j == coarse_grid_neighbors.size()) 00285 { 00286 coarse_grid_neighbors.push_back(coarse_id); 00287 } 00288 00289 } // End loop over fine neighbors 00290 00291 // Set the shared_neighbour index for this node to the 00292 // size of the coarse grid neighbor set 00293 shared_element_count[node_id] = 00294 libmesh_cast_int<unsigned int>(coarse_grid_neighbors.size()); 00295 00296 // Add this node to processed_node_ids vector 00297 processed_node_ids.push_back(node_id); 00298 00299 } // End if processed node 00300 00301 } // End loop over nodes 00302 00303 } // End loop over elems 00304 00305 } // if (split_shared_dofs) 00306 00307 // Get a DoF map, will be used to get the nodal dof_indices for each element 00308 DofMap &dof_map = system.get_dof_map(); 00309 00310 // The global DOF indices, we will use these later on when we compute the element wise indicators 00311 std::vector<dof_id_type> dof_indices; 00312 00313 // Localize the global rhs and adjoint solution vectors (which might be shared on multiple processsors) onto a 00314 // local ghosted vector, this ensures each processor has all the dof_indices to compute an error indicator for 00315 // an element it owns 00316 AutoPtr<NumericVector<Number> > localized_projected_residual = NumericVector<Number>::build(); 00317 localized_projected_residual->init(system.n_dofs(), system.n_local_dofs(), system.get_dof_map().get_send_list(), false, GHOSTED); 00318 projected_residual.localize(*localized_projected_residual, system.get_dof_map().get_send_list()); 00319 00320 // Each adjoint solution will also require ghosting; for efficiency we'll reuse the same memory 00321 AutoPtr<NumericVector<Number> > localized_adjoint_solution = NumericVector<Number>::build(); 00322 localized_adjoint_solution->init(system.n_dofs(), system.n_local_dofs(), system.get_dof_map().get_send_list(), false, GHOSTED); 00323 00324 // We will loop over each adjoint solution, localize that adjoint 00325 // solution and then loop over local elements 00326 for (unsigned int i=0; i != system.qoi.size(); i++) 00327 { 00328 // Skip this QoI if not in the QoI Set 00329 if (_qoi_set.has_index(i)) 00330 { 00331 // Get the weight for the current QoI 00332 Real error_weight = _qoi_set.weight(i); 00333 00334 (system.get_adjoint_solution(i)).localize(*localized_adjoint_solution, system.get_dof_map().get_send_list()); 00335 00336 // Loop over elements 00337 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 00338 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 00339 00340 for(; elem_it != elem_end; ++elem_it) 00341 { 00342 // Pointer to the element 00343 const Elem* elem = *elem_it; 00344 00345 // Go up number_h_refinements levels up to find the coarse parent 00346 const Elem* coarse = elem; 00347 00348 for (unsigned int j = 0; j != number_h_refinements; ++j) 00349 { 00350 libmesh_assert (coarse->parent()); 00351 00352 coarse = coarse->parent(); 00353 } 00354 00355 const dof_id_type e_id = coarse->id(); 00356 00357 // Get the local to global degree of freedom maps for this element 00358 dof_map.dof_indices (elem, dof_indices); 00359 00360 // We will have to manually do the dot products. 00361 Number local_contribution = 0.; 00362 00363 for (unsigned int j=0; j != dof_indices.size(); j++) 00364 { 00365 // The contribution to the error indicator for this element from the current QoI 00366 local_contribution += (*localized_projected_residual)(dof_indices[j]) * (*localized_adjoint_solution)(dof_indices[j]); 00367 } 00368 00369 // Multiply by the error weight for this QoI 00370 local_contribution *= error_weight; 00371 00372 // FIXME: we're throwing away information in the 00373 // --enable-complex case 00374 error_per_cell[e_id] += static_cast<ErrorVectorReal> 00375 (libmesh_real(local_contribution)); 00376 00377 } // End loop over elements 00378 00379 } // End if belong to QoI set 00380 00381 } // End loop over QoIs 00382 00383 // Don't bother projecting the solution; we'll restore from backup 00384 // after coarsening 00385 system.project_solution_on_reinit() = false; 00386 00387 // Uniformly coarsen the mesh, without projecting the solution 00388 libmesh_assert (number_h_refinements > 0 || number_p_refinements > 0); 00389 00390 for (unsigned int i = 0; i != number_h_refinements; ++i) 00391 { 00392 mesh_refinement.uniformly_coarsen(1); 00393 // FIXME - should the reinits here be necessary? - RHS 00394 es.reinit(); 00395 } 00396 00397 for (unsigned int i = 0; i != number_p_refinements; ++i) 00398 { 00399 mesh_refinement.uniformly_p_coarsen(1); 00400 es.reinit(); 00401 } 00402 00403 // We should be back where we started 00404 libmesh_assert_equal_to (n_coarse_elem, mesh.n_elem()); 00405 00406 // Restore old solutions and clean up the heap 00407 system.project_solution_on_reinit() = old_projection_setting; 00408 00409 // Restore the coarse solution vectors and delete their copies 00410 *system.solution = *coarse_solution; 00411 delete coarse_solution; 00412 *system.current_local_solution = *coarse_local_solution; 00413 delete coarse_local_solution; 00414 delete projected_solution; 00415 00416 for (System::vectors_iterator vec = system.vectors_begin(); vec != 00417 system.vectors_end(); ++vec) 00418 { 00419 // The (string) name of this vector 00420 const std::string& var_name = vec->first; 00421 00422 system.get_vector(var_name) = *coarse_vectors[var_name]; 00423 00424 coarse_vectors[var_name]->clear(); 00425 delete coarse_vectors[var_name]; 00426 } 00427 00428 // Restore old partitioner and renumbering settings 00429 mesh.partitioner() = old_partitioner; 00430 mesh.allow_renumbering(old_renumbering_setting); 00431 00432 } // end estimate_error function
| void libMesh::ErrorEstimator::estimate_errors | ( | const EquationSystems & | equation_systems, | |
| ErrorMap & | errors_per_cell, | |||
| const std::map< const System *, const NumericVector< Number > * > * | solution_vectors = NULL, |
|||
| bool | estimate_parent_error = false | |||
| ) | [virtual, inherited] |
This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors.
Currently this function ignores the error_norm.weight() values because it calculates each variable's error individually, unscaled.
The user selects which errors get computed by filling a map with error vectors: If errors_per_cell[&system][v] exists, it will be filled with the error values in variable v of system
FIXME: This is a default implementation - derived classes should reimplement it for efficiency.
Reimplemented in libMesh::UniformRefinementEstimator.
Definition at line 93 of file error_estimator.C.
References libMesh::ErrorEstimator::error_norm, libMesh::ErrorEstimator::estimate_error(), libMesh::EquationSystems::get_system(), libMesh::EquationSystems::n_systems(), libMesh::System::n_vars(), n_vars, and libMesh::SystemNorm::type().
00097 { 00098 SystemNorm old_error_norm = this->error_norm; 00099 00100 // Find the requested error values from each system 00101 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00102 { 00103 const System &sys = equation_systems.get_system(s); 00104 00105 unsigned int n_vars = sys.n_vars(); 00106 00107 for (unsigned int v = 0; v != n_vars; ++v) 00108 { 00109 // Only fill in ErrorVectors the user asks for 00110 if (errors_per_cell.find(std::make_pair(&sys, v)) == 00111 errors_per_cell.end()) 00112 continue; 00113 00114 // Calculate error in only one variable 00115 std::vector<Real> weights(n_vars, 0.0); 00116 weights[v] = 1.0; 00117 this->error_norm = 00118 SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)), 00119 weights); 00120 00121 const NumericVector<Number>* solution_vector = NULL; 00122 if (solution_vectors && 00123 solution_vectors->find(&sys) != solution_vectors->end()) 00124 solution_vector = solution_vectors->find(&sys)->second; 00125 00126 this->estimate_error 00127 (sys, *errors_per_cell[std::make_pair(&sys, v)], 00128 solution_vector, estimate_parent_error); 00129 } 00130 } 00131 00132 // Restore our old state before returning 00133 this->error_norm = old_error_norm; 00134 }
| void libMesh::ErrorEstimator::estimate_errors | ( | const EquationSystems & | equation_systems, | |
| ErrorVector & | error_per_cell, | |||
| const std::map< const System *, SystemNorm > & | error_norms, | |||
| const std::map< const System *, const NumericVector< Number > * > * | solution_vectors = NULL, |
|||
| bool | estimate_parent_error = false | |||
| ) | [virtual, inherited] |
This virtual function can be redefined in derived classes, but by default computes the sum of the error_per_cell for each system in the equation_systems.
Currently this function ignores the error_norm member variable, and uses the function argument error_norms instead.
This function is named estimate_errors instead of estimate_error because otherwise C++ can get confused.
Reimplemented in libMesh::UniformRefinementEstimator.
Definition at line 47 of file error_estimator.C.
References libMesh::ErrorEstimator::error_norm, libMesh::ErrorEstimator::estimate_error(), libMesh::EquationSystems::get_system(), and libMesh::EquationSystems::n_systems().
00052 { 00053 SystemNorm old_error_norm = this->error_norm; 00054 00055 // Sum the error values from each system 00056 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00057 { 00058 ErrorVector system_error_per_cell; 00059 const System &sys = equation_systems.get_system(s); 00060 if (error_norms.find(&sys) == error_norms.end()) 00061 this->error_norm = old_error_norm; 00062 else 00063 this->error_norm = error_norms.find(&sys)->second; 00064 00065 const NumericVector<Number>* solution_vector = NULL; 00066 if (solution_vectors && 00067 solution_vectors->find(&sys) != solution_vectors->end()) 00068 solution_vector = solution_vectors->find(&sys)->second; 00069 00070 this->estimate_error(sys, system_error_per_cell, 00071 solution_vector, estimate_parent_error); 00072 00073 if (s) 00074 { 00075 libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size()); 00076 for (unsigned int i=0; i != error_per_cell.size(); ++i) 00077 error_per_cell[i] += system_error_per_cell[i]; 00078 } 00079 else 00080 error_per_cell = system_error_per_cell; 00081 } 00082 00083 // Restore our old state before returning 00084 this->error_norm = old_error_norm; 00085 }
| Number& libMesh::AdjointRefinementEstimator::get_global_QoI_error_estimate | ( | unsigned int | qoi_index | ) | [inline] |
This is an accessor function to access the computed global QoI error estimates
Definition at line 104 of file adjoint_refinement_estimator.h.
References computed_global_QoI_errors.
00105 { 00106 return computed_global_QoI_errors[qoi_index]; 00107 }
| const QoISet& libMesh::AdjointRefinementEstimator::qoi_set | ( | ) | const [inline] |
Access to the QoISet (default: weight all QoIs equally) to use when computing errors
Definition at line 77 of file adjoint_refinement_estimator.h.
References _qoi_set.
00077 { return _qoi_set; }
| QoISet& libMesh::AdjointRefinementEstimator::qoi_set | ( | ) | [inline] |
Access to the QoISet (default: weight all QoIs equally) to use when computing errors
Definition at line 71 of file adjoint_refinement_estimator.h.
References _qoi_set.
00071 { return _qoi_set; }
| void libMesh::ErrorEstimator::reduce_error | ( | std::vector< float > & | error_per_cell | ) | const [protected, inherited] |
This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector.
Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), and libMesh::JumpErrorEstimator::estimate_error().
Member Data Documentation
QoISet libMesh::AdjointRefinementEstimator::_qoi_set [protected] |
A QoISet to handle cases with multiple QoIs available
Definition at line 138 of file adjoint_refinement_estimator.h.
Referenced by estimate_error(), and qoi_set().
std::vector<Number> libMesh::AdjointRefinementEstimator::computed_global_QoI_errors [protected] |
Definition at line 122 of file adjoint_refinement_estimator.h.
Referenced by estimate_error(), and get_global_QoI_error_estimate().
SystemNorm libMesh::ErrorEstimator::error_norm [inherited] |
When estimating the error in a single system, the error_norm is used to control the scaling and norm choice for each variable. Not all estimators will support all norm choices. The default scaling is for all variables to be weighted equally. The default norm choice depends on the error estimator.
Part of this functionality was supported via component_scale and sobolev_order in older libMesh versions, and a small part was supported via component_mask in even older versions. Hopefully the encapsulation here will allow us to avoid changing this API again.
Definition at line 139 of file error_estimator.h.
Referenced by AdjointRefinementEstimator(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::DiscontinuityMeasure::DiscontinuityMeasure(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointResidualErrorEstimator::estimate_error(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::KellyErrorEstimator(), libMesh::LaplacianErrorEstimator::LaplacianErrorEstimator(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().
How many h refinements to perform to get the fine grid
Definition at line 112 of file adjoint_refinement_estimator.h.
Referenced by estimate_error().
How many p refinements to perform to get the fine grid
Definition at line 117 of file adjoint_refinement_estimator.h.
Referenced by estimate_error().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:06 UTC
Hosted By: