distributed_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 #include "libmesh/libmesh_common.h" 00021 00022 // C++ includes 00023 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00024 #include <cmath> // for std::abs 00025 #include <limits> // std::numeric_limits<T>::min() 00026 00027 // Local Includes 00028 #include "libmesh/distributed_vector.h" 00029 #include "libmesh/dense_vector.h" 00030 #include "libmesh/dense_subvector.h" 00031 #include "libmesh/parallel.h" 00032 #include "libmesh/tensor_tools.h" 00033 00034 namespace libMesh 00035 { 00036 00037 00038 00039 //-------------------------------------------------------------------------- 00040 // DistributedVector methods 00041 template <typename T> 00042 T DistributedVector<T>::sum () const 00043 { 00044 // This function must be run on all processors at once 00045 parallel_only(); 00046 00047 libmesh_assert (this->initialized()); 00048 libmesh_assert_equal_to (_values.size(), _local_size); 00049 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00050 00051 T local_sum = 0.; 00052 00053 for (numeric_index_type i=0; i<local_size(); i++) 00054 local_sum += _values[i]; 00055 00056 CommWorld.sum(local_sum); 00057 00058 return local_sum; 00059 } 00060 00061 00062 00063 template <typename T> 00064 Real DistributedVector<T>::l1_norm () const 00065 { 00066 // This function must be run on all processors at once 00067 parallel_only(); 00068 00069 libmesh_assert (this->initialized()); 00070 libmesh_assert_equal_to (_values.size(), _local_size); 00071 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00072 00073 double local_l1 = 0.; 00074 00075 for (numeric_index_type i=0; i<local_size(); i++) 00076 local_l1 += std::abs(_values[i]); 00077 00078 CommWorld.sum(local_l1); 00079 00080 return local_l1; 00081 } 00082 00083 00084 00085 template <typename T> 00086 Real DistributedVector<T>::l2_norm () const 00087 { 00088 // This function must be run on all processors at once 00089 parallel_only(); 00090 00091 libmesh_assert (this->initialized()); 00092 libmesh_assert_equal_to (_values.size(), _local_size); 00093 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00094 00095 double local_l2 = 0.; 00096 00097 for (numeric_index_type i=0; i<local_size(); i++) 00098 local_l2 += TensorTools::norm_sq(_values[i]); 00099 00100 CommWorld.sum(local_l2); 00101 00102 return std::sqrt(local_l2); 00103 } 00104 00105 00106 00107 template <typename T> 00108 Real DistributedVector<T>::linfty_norm () const 00109 { 00110 // This function must be run on all processors at once 00111 parallel_only(); 00112 00113 libmesh_assert (this->initialized()); 00114 libmesh_assert_equal_to (_values.size(), _local_size); 00115 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00116 00117 Real local_linfty = 0.; 00118 00119 for (numeric_index_type i=0; i<local_size(); i++) 00120 local_linfty = std::max(local_linfty, 00121 static_cast<Real>(std::abs(_values[i])) 00122 ); // Note we static_cast so that both 00123 // types are the same, as required 00124 // by std::max 00125 00126 CommWorld.max(local_linfty); 00127 00128 return local_linfty; 00129 } 00130 00131 00132 00133 template <typename T> 00134 NumericVector<T>& DistributedVector<T>::operator += (const NumericVector<T>& v) 00135 { 00136 libmesh_assert (this->closed()); 00137 libmesh_assert (this->initialized()); 00138 libmesh_assert_equal_to (_values.size(), _local_size); 00139 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00140 00141 add(1., v); 00142 00143 return *this; 00144 } 00145 00146 00147 00148 template <typename T> 00149 NumericVector<T>& DistributedVector<T>::operator -= (const NumericVector<T>& v) 00150 { 00151 libmesh_assert (this->closed()); 00152 libmesh_assert (this->initialized()); 00153 libmesh_assert_equal_to (_values.size(), _local_size); 00154 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00155 00156 add(-1., v); 00157 00158 return *this; 00159 } 00160 00161 00162 00163 template <typename T> 00164 void DistributedVector<T>::add_vector (const std::vector<T>& v, 00165 const std::vector<numeric_index_type>& dof_indices) 00166 { 00167 libmesh_assert (!v.empty()); 00168 libmesh_assert_equal_to (v.size(), dof_indices.size()); 00169 libmesh_assert (this->initialized()); 00170 libmesh_assert_equal_to (_values.size(), _local_size); 00171 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00172 00173 for (std::size_t i=0; i<v.size(); i++) 00174 add (dof_indices[i], v[i]); 00175 } 00176 00177 00178 00179 template <typename T> 00180 void DistributedVector<T>::add_vector (const NumericVector<T>& V, 00181 const std::vector<numeric_index_type>& dof_indices) 00182 { 00183 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00184 libmesh_assert (this->initialized()); 00185 libmesh_assert_equal_to (_values.size(), _local_size); 00186 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00187 00188 for (std::size_t i=0; i<V.size(); i++) 00189 add (dof_indices[i], V(i)); 00190 } 00191 00192 00193 00194 template <typename T> 00195 void DistributedVector<T>::add_vector (const DenseVector<T>& V, 00196 const std::vector<numeric_index_type>& dof_indices) 00197 { 00198 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00199 libmesh_assert (this->initialized()); 00200 libmesh_assert_equal_to (_values.size(), _local_size); 00201 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00202 00203 for (unsigned int i=0; i<V.size(); i++) 00204 add (dof_indices[i], V(i)); 00205 } 00206 00207 00208 00209 00210 template <typename T> 00211 void DistributedVector<T>::reciprocal() 00212 { 00213 for (numeric_index_type i=0; i<local_size(); i++) 00214 { 00215 // Don't divide by zero 00216 libmesh_assert_not_equal_to (_values[i], T(0)); 00217 00218 _values[i] = 1. / _values[i]; 00219 } 00220 } 00221 00222 00223 00224 00225 template <typename T> 00226 void DistributedVector<T>::add (const T v) 00227 { 00228 libmesh_assert (this->initialized()); 00229 libmesh_assert_equal_to (_values.size(), _local_size); 00230 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00231 00232 for (numeric_index_type i=0; i<local_size(); i++) 00233 _values[i] += v; 00234 } 00235 00236 00237 00238 template <typename T> 00239 void DistributedVector<T>::add (const NumericVector<T>& v) 00240 { 00241 libmesh_assert (this->initialized()); 00242 libmesh_assert_equal_to (_values.size(), _local_size); 00243 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00244 00245 add (1., v); 00246 } 00247 00248 00249 00250 template <typename T> 00251 void DistributedVector<T>::add (const T a, const NumericVector<T>& v) 00252 { 00253 libmesh_assert (this->initialized()); 00254 libmesh_assert_equal_to (_values.size(), _local_size); 00255 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00256 00257 add(a, v); 00258 } 00259 00260 00261 00262 template <typename T> 00263 void DistributedVector<T>::insert (const std::vector<T>& v, 00264 const std::vector<numeric_index_type>& dof_indices) 00265 { 00266 libmesh_assert (!v.empty()); 00267 libmesh_assert_equal_to (v.size(), dof_indices.size()); 00268 libmesh_assert (this->initialized()); 00269 libmesh_assert_equal_to (_values.size(), _local_size); 00270 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00271 00272 for (std::size_t i=0; i<v.size(); i++) 00273 this->set (dof_indices[i], v[i]); 00274 } 00275 00276 00277 00278 template <typename T> 00279 void DistributedVector<T>::insert (const NumericVector<T>& V, 00280 const std::vector<numeric_index_type>& dof_indices) 00281 { 00282 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00283 libmesh_assert (this->initialized()); 00284 libmesh_assert_equal_to (_values.size(), _local_size); 00285 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00286 00287 for (std::size_t i=0; i<V.size(); i++) 00288 this->set (dof_indices[i], V(i)); 00289 } 00290 00291 00292 00293 template <typename T> 00294 void DistributedVector<T>::insert (const DenseVector<T>& V, 00295 const std::vector<numeric_index_type>& dof_indices) 00296 { 00297 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00298 libmesh_assert (this->initialized()); 00299 libmesh_assert_equal_to (_values.size(), _local_size); 00300 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00301 00302 for (unsigned int i=0; i<V.size(); i++) 00303 this->set (dof_indices[i], V(i)); 00304 } 00305 00306 00307 00308 template <typename T> 00309 void DistributedVector<T>::insert (const DenseSubVector<T>& V, 00310 const std::vector<numeric_index_type>& dof_indices) 00311 { 00312 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00313 libmesh_assert (this->initialized()); 00314 libmesh_assert_equal_to (_values.size(), _local_size); 00315 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00316 00317 for (unsigned int i=0; i<V.size(); i++) 00318 this->set (dof_indices[i], V(i)); 00319 } 00320 00321 00322 00323 template <typename T> 00324 void DistributedVector<T>::scale (const T factor) 00325 { 00326 libmesh_assert (this->initialized()); 00327 libmesh_assert_equal_to (_values.size(), _local_size); 00328 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00329 00330 for (std::size_t i=0; i<local_size(); i++) 00331 _values[i] *= factor; 00332 } 00333 00334 template <typename T> 00335 void DistributedVector<T>::abs() 00336 { 00337 libmesh_assert (this->initialized()); 00338 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00339 00340 for (std::size_t i=0; i<local_size(); i++) 00341 this->set(i,std::abs(_values[i])); 00342 } 00343 00344 00345 00346 00347 00348 template <typename T> 00349 T DistributedVector<T>::dot (const NumericVector<T>& V) const 00350 { 00351 // This function must be run on all processors at once 00352 parallel_only(); 00353 00354 // Make sure the NumericVector passed in is really a DistributedVector 00355 const DistributedVector<T>* v = libmesh_cast_ptr<const DistributedVector<T>*>(&V); 00356 00357 // Make sure that the two vectors are distributed in the same way. 00358 libmesh_assert_equal_to ( this->first_local_index(), v->first_local_index() ); 00359 libmesh_assert_equal_to ( this->last_local_index(), v->last_local_index() ); 00360 00361 // The result of dotting together the local parts of the vector. 00362 T local_dot = 0; 00363 00364 for (std::size_t i=0; i<this->local_size(); i++) 00365 local_dot += this->_values[i] * v->_values[i]; 00366 00367 // The local dot products are now summed via MPI 00368 CommWorld.sum(local_dot); 00369 00370 return local_dot; 00371 } 00372 00373 00374 00375 template <typename T> 00376 NumericVector<T>& 00377 DistributedVector<T>::operator = (const T s) 00378 { 00379 libmesh_assert (this->initialized()); 00380 libmesh_assert_equal_to (_values.size(), _local_size); 00381 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00382 00383 for (std::size_t i=0; i<local_size(); i++) 00384 _values[i] = s; 00385 00386 return *this; 00387 } 00388 00389 00390 00391 template <typename T> 00392 NumericVector<T>& 00393 DistributedVector<T>::operator = (const NumericVector<T>& v_in) 00394 { 00395 // Make sure the NumericVector passed in is really a DistributedVector 00396 const DistributedVector<T>* v = libmesh_cast_ptr<const DistributedVector<T>*>(&v_in); 00397 00398 *this = *v; 00399 00400 return *this; 00401 } 00402 00403 00404 00405 template <typename T> 00406 DistributedVector<T>& 00407 DistributedVector<T>::operator = (const DistributedVector<T>& v) 00408 { 00409 this->_is_initialized = v._is_initialized; 00410 this->_is_closed = v._is_closed; 00411 00412 _global_size = v._global_size; 00413 _local_size = v._local_size; 00414 _first_local_index = v._first_local_index; 00415 _last_local_index = v._last_local_index; 00416 00417 if (v.local_size() == this->local_size()) 00418 { 00419 _values = v._values; 00420 } 00421 else 00422 { 00423 libmesh_error(); 00424 } 00425 00426 return *this; 00427 } 00428 00429 00430 00431 template <typename T> 00432 NumericVector<T>& 00433 DistributedVector<T>::operator = (const std::vector<T>& v) 00434 { 00435 libmesh_assert (this->initialized()); 00436 libmesh_assert_equal_to (_values.size(), _local_size); 00437 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00438 00439 if (v.size() == local_size()) 00440 _values = v; 00441 00442 else if (v.size() == size()) 00443 for (std::size_t i=first_local_index(); i<last_local_index(); i++) 00444 _values[i-first_local_index()] = v[i]; 00445 00446 else 00447 { 00448 libmesh_error(); 00449 } 00450 00451 00452 return *this; 00453 } 00454 00455 00456 00457 template <typename T> 00458 void DistributedVector<T>::localize (NumericVector<T>& v_local_in) const 00459 00460 { 00461 libmesh_assert (this->initialized()); 00462 libmesh_assert_equal_to (_values.size(), _local_size); 00463 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00464 00465 DistributedVector<T>* v_local = libmesh_cast_ptr<DistributedVector<T>*>(&v_local_in); 00466 00467 v_local->_first_local_index = 0; 00468 00469 v_local->_global_size = 00470 v_local->_local_size = 00471 v_local->_last_local_index = size(); 00472 00473 v_local->_is_initialized = 00474 v_local->_is_closed = true; 00475 00476 // Call localize on the vector's values. This will help 00477 // prevent code duplication 00478 localize (v_local->_values); 00479 00480 #ifndef LIBMESH_HAVE_MPI 00481 00482 libmesh_assert_equal_to (local_size(), size()); 00483 00484 #endif 00485 } 00486 00487 00488 00489 template <typename T> 00490 void DistributedVector<T>::localize (NumericVector<T>& v_local_in, 00491 const std::vector<numeric_index_type>&) const 00492 { 00493 libmesh_assert (this->initialized()); 00494 libmesh_assert_equal_to (_values.size(), _local_size); 00495 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00496 00497 // TODO: We don't yet support the send list; this is inefficient: 00498 localize (v_local_in); 00499 } 00500 00501 00502 00503 template <typename T> 00504 void DistributedVector<T>::localize (const numeric_index_type first_local_idx, 00505 const numeric_index_type last_local_idx, 00506 const std::vector<numeric_index_type>& send_list) 00507 { 00508 // Only good for serial vectors 00509 libmesh_assert_equal_to (this->size(), this->local_size()); 00510 libmesh_assert_greater (last_local_idx, first_local_idx); 00511 libmesh_assert_less_equal (send_list.size(), this->size()); 00512 libmesh_assert_less (last_local_idx, this->size()); 00513 00514 const numeric_index_type my_size = this->size(); 00515 const numeric_index_type my_local_size = (last_local_idx - first_local_idx + 1); 00516 00517 // Don't bother for serial cases 00518 if ((first_local_idx == 0) && 00519 (my_local_size == my_size)) 00520 return; 00521 00522 00523 // Build a parallel vector, initialize it with the local 00524 // parts of (*this) 00525 DistributedVector<T> parallel_vec; 00526 00527 parallel_vec.init (my_size, my_local_size, true, PARALLEL); 00528 00529 // Copy part of *this into the parallel_vec 00530 for (numeric_index_type i=first_local_idx; i<=last_local_idx; i++) 00531 parallel_vec._values[i-first_local_idx] = _values[i]; 00532 00533 // localize like normal 00534 parallel_vec.localize (*this, send_list); 00535 } 00536 00537 00538 00539 template <typename T> 00540 void DistributedVector<T>::localize (std::vector<T>& v_local) const 00541 { 00542 // This function must be run on all processors at once 00543 parallel_only(); 00544 00545 libmesh_assert (this->initialized()); 00546 libmesh_assert_equal_to (_values.size(), _local_size); 00547 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00548 00549 v_local = this->_values; 00550 00551 CommWorld.allgather (v_local); 00552 00553 #ifndef LIBMESH_HAVE_MPI 00554 libmesh_assert_equal_to (local_size(), size()); 00555 #endif 00556 } 00557 00558 00559 00560 template <typename T> 00561 void DistributedVector<T>::localize_to_one (std::vector<T>& v_local, 00562 const processor_id_type pid) const 00563 { 00564 // This function must be run on all processors at once 00565 parallel_only(); 00566 00567 libmesh_assert (this->initialized()); 00568 libmesh_assert_equal_to (_values.size(), _local_size); 00569 libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size); 00570 00571 v_local = this->_values; 00572 00573 CommWorld.gather (pid, v_local); 00574 00575 #ifndef LIBMESH_HAVE_MPI 00576 libmesh_assert_equal_to (local_size(), size()); 00577 #endif 00578 } 00579 00580 00581 00582 template <typename T> 00583 void DistributedVector<T>::pointwise_mult (const NumericVector<T>&, 00584 const NumericVector<T>&) 00585 //void DistributedVector<T>::pointwise_mult (const NumericVector<T>& vec1, 00586 // const NumericVector<T>& vec2) 00587 { 00588 libmesh_not_implemented(); 00589 } 00590 00591 00592 00593 //-------------------------------------------------------------- 00594 // Explicit instantiations 00595 template class DistributedVector<Number>; 00596 00597 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: