$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 // 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, 00034 const Parallel::Communicator &comm) const 00035 { 00036 // This function must be run on all processors at once 00037 // parallel_object_only(); 00038 00039 // Each processor has now computed the error contribuions 00040 // for its local elements. We may need to sum the vector to 00041 // recover the error for each element. 00042 00043 comm.sum(error_per_cell); 00044 } 00045 00046 00047 00048 void ErrorEstimator::estimate_errors(const EquationSystems& equation_systems, 00049 ErrorVector& error_per_cell, 00050 const std::map<const System*, SystemNorm>& error_norms, 00051 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00052 bool estimate_parent_error) 00053 { 00054 SystemNorm old_error_norm = this->error_norm; 00055 00056 // Sum the error values from each system 00057 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00058 { 00059 ErrorVector system_error_per_cell; 00060 const System &sys = equation_systems.get_system(s); 00061 if (error_norms.find(&sys) == error_norms.end()) 00062 this->error_norm = old_error_norm; 00063 else 00064 this->error_norm = error_norms.find(&sys)->second; 00065 00066 const NumericVector<Number>* solution_vector = NULL; 00067 if (solution_vectors && 00068 solution_vectors->find(&sys) != solution_vectors->end()) 00069 solution_vector = solution_vectors->find(&sys)->second; 00070 00071 this->estimate_error(sys, system_error_per_cell, 00072 solution_vector, estimate_parent_error); 00073 00074 if (s) 00075 { 00076 libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size()); 00077 for (unsigned int i=0; i != error_per_cell.size(); ++i) 00078 error_per_cell[i] += system_error_per_cell[i]; 00079 } 00080 else 00081 error_per_cell = system_error_per_cell; 00082 } 00083 00084 // Restore our old state before returning 00085 this->error_norm = old_error_norm; 00086 } 00087 00088 00089 00094 void ErrorEstimator::estimate_errors(const EquationSystems& equation_systems, 00095 ErrorMap& errors_per_cell, 00096 const std::map<const System*, const NumericVector<Number>* >* solution_vectors, 00097 bool estimate_parent_error) 00098 { 00099 SystemNorm old_error_norm = this->error_norm; 00100 00101 // Find the requested error values from each system 00102 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00103 { 00104 const System &sys = equation_systems.get_system(s); 00105 00106 unsigned int n_vars = sys.n_vars(); 00107 00108 for (unsigned int v = 0; v != n_vars; ++v) 00109 { 00110 // Only fill in ErrorVectors the user asks for 00111 if (errors_per_cell.find(std::make_pair(&sys, v)) == 00112 errors_per_cell.end()) 00113 continue; 00114 00115 // Calculate error in only one variable 00116 std::vector<Real> weights(n_vars, 0.0); 00117 weights[v] = 1.0; 00118 this->error_norm = 00119 SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)), 00120 weights); 00121 00122 const NumericVector<Number>* solution_vector = NULL; 00123 if (solution_vectors && 00124 solution_vectors->find(&sys) != solution_vectors->end()) 00125 solution_vector = solution_vectors->find(&sys)->second; 00126 00127 this->estimate_error 00128 (sys, *errors_per_cell[std::make_pair(&sys, v)], 00129 solution_vector, estimate_parent_error); 00130 } 00131 } 00132 00133 // Restore our old state before returning 00134 this->error_norm = old_error_norm; 00135 } 00136 00137 } // namespace libMesh