numeric_vector.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 00020 // C++ includes 00021 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00022 #include <cmath> // for std::abs 00023 #include <limits> 00024 00025 // Local Includes 00026 #include "libmesh/numeric_vector.h" 00027 #include "libmesh/distributed_vector.h" 00028 #include "libmesh/laspack_vector.h" 00029 #include "libmesh/petsc_vector.h" 00030 #include "libmesh/trilinos_epetra_vector.h" 00031 #include "libmesh/shell_matrix.h" 00032 #include "libmesh/tensor_tools.h" 00033 00034 namespace libMesh 00035 { 00036 00037 00038 00039 //------------------------------------------------------------------ 00040 // NumericVector methods 00041 00042 // Full specialization for Real datatypes 00043 template <typename T> 00044 AutoPtr<NumericVector<T> > 00045 NumericVector<T>::build(const SolverPackage solver_package) 00046 { 00047 // Build the appropriate vector 00048 switch (solver_package) 00049 { 00050 00051 00052 #ifdef LIBMESH_HAVE_LASPACK 00053 case LASPACK_SOLVERS: 00054 { 00055 AutoPtr<NumericVector<T> > ap(new LaspackVector<T>); 00056 return ap; 00057 } 00058 #endif 00059 00060 00061 #ifdef LIBMESH_HAVE_PETSC 00062 case PETSC_SOLVERS: 00063 { 00064 AutoPtr<NumericVector<T> > ap(new PetscVector<T>); 00065 return ap; 00066 } 00067 #endif 00068 00069 00070 #ifdef LIBMESH_HAVE_TRILINOS 00071 case TRILINOS_SOLVERS: 00072 { 00073 AutoPtr<NumericVector<T> > ap(new EpetraVector<T>); 00074 return ap; 00075 } 00076 #endif 00077 00078 00079 default: 00080 AutoPtr<NumericVector<T> > ap(new DistributedVector<T>); 00081 return ap; 00082 } 00083 00084 AutoPtr<NumericVector<T> > ap(NULL); 00085 return ap; 00086 } 00087 00088 00089 template <typename T> 00090 int NumericVector<T>::compare (const NumericVector<T> &other_vector, 00091 const Real threshold) const 00092 { 00093 libmesh_assert (this->initialized()); 00094 libmesh_assert (other_vector.initialized()); 00095 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00096 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00097 00098 int first_different_i = std::numeric_limits<int>::max(); 00099 numeric_index_type i = first_local_index(); 00100 00101 do 00102 { 00103 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold ) 00104 first_different_i = i; 00105 else 00106 i++; 00107 } 00108 while (first_different_i==std::numeric_limits<int>::max() 00109 && i<last_local_index()); 00110 00111 // Find the correct first differing index in parallel 00112 CommWorld.min(first_different_i); 00113 00114 if (first_different_i == std::numeric_limits<int>::max()) 00115 return -1; 00116 00117 return first_different_i; 00118 } 00119 00120 00121 template <typename T> 00122 int NumericVector<T>::local_relative_compare (const NumericVector<T> &other_vector, 00123 const Real threshold) const 00124 { 00125 libmesh_assert (this->initialized()); 00126 libmesh_assert (other_vector.initialized()); 00127 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00128 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00129 00130 int first_different_i = std::numeric_limits<int>::max(); 00131 numeric_index_type i = first_local_index(); 00132 00133 do 00134 { 00135 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold * 00136 std::max(std::abs((*this)(i)), std::abs(other_vector(i)))) 00137 first_different_i = i; 00138 else 00139 i++; 00140 } 00141 while (first_different_i==std::numeric_limits<int>::max() 00142 && i<last_local_index()); 00143 00144 // Find the correct first differing index in parallel 00145 CommWorld.min(first_different_i); 00146 00147 if (first_different_i == std::numeric_limits<int>::max()) 00148 return -1; 00149 00150 return first_different_i; 00151 } 00152 00153 00154 template <typename T> 00155 int NumericVector<T>::global_relative_compare (const NumericVector<T> &other_vector, 00156 const Real threshold) const 00157 { 00158 libmesh_assert (this->initialized()); 00159 libmesh_assert (other_vector.initialized()); 00160 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00161 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00162 00163 int first_different_i = std::numeric_limits<int>::max(); 00164 numeric_index_type i = first_local_index(); 00165 00166 const Real my_norm = this->linfty_norm(); 00167 const Real other_norm = other_vector.linfty_norm(); 00168 const Real abs_threshold = std::max(my_norm, other_norm) * threshold; 00169 00170 do 00171 { 00172 if ( std::abs( (*this)(i) - other_vector(i) ) > abs_threshold ) 00173 first_different_i = i; 00174 else 00175 i++; 00176 } 00177 while (first_different_i==std::numeric_limits<int>::max() 00178 && i<last_local_index()); 00179 00180 // Find the correct first differing index in parallel 00181 CommWorld.min(first_different_i); 00182 00183 if (first_different_i == std::numeric_limits<int>::max()) 00184 return -1; 00185 00186 return first_different_i; 00187 } 00188 00189 /* 00190 // Full specialization for float datatypes (DistributedVector wants this) 00191 00192 template <> 00193 int NumericVector<float>::compare (const NumericVector<float> &other_vector, 00194 const Real threshold) const 00195 { 00196 libmesh_assert (this->initialized()); 00197 libmesh_assert (other_vector.initialized()); 00198 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00199 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00200 00201 int rvalue = -1; 00202 numeric_index_type i = first_local_index(); 00203 00204 do 00205 { 00206 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold ) 00207 rvalue = i; 00208 else 00209 i++; 00210 } 00211 while (rvalue==-1 && i<last_local_index()); 00212 00213 return rvalue; 00214 } 00215 00216 // Full specialization for double datatypes 00217 template <> 00218 int NumericVector<double>::compare (const NumericVector<double> &other_vector, 00219 const Real threshold) const 00220 { 00221 libmesh_assert (this->initialized()); 00222 libmesh_assert (other_vector.initialized()); 00223 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00224 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00225 00226 int rvalue = -1; 00227 numeric_index_type i = first_local_index(); 00228 00229 do 00230 { 00231 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold ) 00232 rvalue = i; 00233 else 00234 i++; 00235 } 00236 while (rvalue==-1 && i<last_local_index()); 00237 00238 return rvalue; 00239 } 00240 00241 #ifdef LIBMESH_DEFAULT_TRIPLE_PRECISION 00242 // Full specialization for long double datatypes 00243 template <> 00244 int NumericVector<long double>::compare (const NumericVector<long double> &other_vector, 00245 const Real threshold) const 00246 { 00247 libmesh_assert (this->initialized()); 00248 libmesh_assert (other_vector.initialized()); 00249 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00250 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00251 00252 int rvalue = -1; 00253 numeric_index_type i = first_local_index(); 00254 00255 do 00256 { 00257 if ( std::abs( (*this)(i) - other_vector(i) ) > threshold ) 00258 rvalue = i; 00259 else 00260 i++; 00261 } 00262 while (rvalue==-1 && i<last_local_index()); 00263 00264 return rvalue; 00265 } 00266 #endif 00267 00268 00269 // Full specialization for Complex datatypes 00270 template <> 00271 int NumericVector<Complex>::compare (const NumericVector<Complex> &other_vector, 00272 const Real threshold) const 00273 { 00274 libmesh_assert (this->initialized()); 00275 libmesh_assert (other_vector.initialized()); 00276 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index()); 00277 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index()); 00278 00279 int rvalue = -1; 00280 numeric_index_type i = first_local_index(); 00281 00282 do 00283 { 00284 if (( std::abs( (*this)(i).real() - other_vector(i).real() ) > threshold ) || 00285 ( std::abs( (*this)(i).imag() - other_vector(i).imag() ) > threshold )) 00286 rvalue = i; 00287 else 00288 i++; 00289 } 00290 while (rvalue==-1 && i<this->last_local_index()); 00291 00292 return rvalue; 00293 } 00294 */ 00295 00296 00297 template <class T> 00298 Real NumericVector<T>::subset_l1_norm (const std::set<numeric_index_type> & indices) const 00299 { 00300 const NumericVector<T> & v = *this; 00301 00302 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00303 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00304 00305 Real norm = 0; 00306 00307 for(; it!=it_end; ++it) 00308 norm += std::abs(v(*it)); 00309 00310 CommWorld.sum(norm); 00311 00312 return norm; 00313 } 00314 00315 template <class T> 00316 Real NumericVector<T>::subset_l2_norm (const std::set<numeric_index_type> & indices) const 00317 { 00318 const NumericVector<T> & v = *this; 00319 00320 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00321 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00322 00323 Real norm = 0; 00324 00325 for(; it!=it_end; ++it) 00326 norm += TensorTools::norm_sq(v(*it)); 00327 00328 CommWorld.sum(norm); 00329 00330 return std::sqrt(norm); 00331 } 00332 00333 template <class T> 00334 Real NumericVector<T>::subset_linfty_norm (const std::set<numeric_index_type> & indices) const 00335 { 00336 const NumericVector<T> & v = *this; 00337 00338 std::set<numeric_index_type>::const_iterator it = indices.begin(); 00339 const std::set<numeric_index_type>::const_iterator it_end = indices.end(); 00340 00341 Real norm = 0; 00342 00343 for(; it!=it_end; ++it) 00344 { 00345 Real value = std::abs(v(*it)); 00346 if(value > norm) 00347 norm = value; 00348 } 00349 00350 CommWorld.max(norm); 00351 00352 return norm; 00353 } 00354 00355 00356 00357 template <typename T> 00358 void NumericVector<T>::add_vector (const NumericVector<T>& v, 00359 const ShellMatrix<T>& a) 00360 { 00361 a.vector_mult_add(*this,v); 00362 } 00363 00364 00365 00366 //------------------------------------------------------------------ 00367 // Explicit instantiations 00368 template class NumericVector<Number>; 00369 00370 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: