error_estimator.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // Local include files 00020 #include "libmesh/libmesh_common.h" 00021 #include "libmesh/error_estimator.h" 00022 #include "libmesh/error_vector.h" 00023 #include "libmesh/equation_systems.h" 00024 #include "libmesh/parallel.h" 00025 00026 00027 00028 00029 namespace libMesh 00030 { 00031 00032 // ErrorEstimator functions 00033 void ErrorEstimator::reduce_error (std::vector<ErrorVectorReal>& error_per_cell) const 00034 { 00035 // This function must be run on all processors at once 00036 parallel_only(); 00037 00038 // Each processor has now computed the error contribuions 00039 // for its local elements. We may need to sum the vector to 00040 // recover the error for each element. 00041 00042 CommWorld.sum(error_per_cell); 00043 } 00044 00045 00046 00047 void ErrorEstimator::estimate_errors(const EquationSystems& equation_systems, 00048 ErrorVector& error_per_cell, 00049 const std::map<const System*, SystemNorm>& error_norms, 00050 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00051 bool estimate_parent_error) 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 } 00086 00087 00088 00093 void ErrorEstimator::estimate_errors(const EquationSystems& equation_systems, 00094 ErrorMap& errors_per_cell, 00095 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00096 bool estimate_parent_error) 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 } 00135 00136 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: