Reference documentation for deal.II version 8.5.1
trilinos_epetra_vector.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_epetra_vector.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 #ifdef DEAL_II_WITH_MPI
21 
22 #include <deal.II/base/index_set.h>
23 
24 #include <boost/io/ios_state.hpp>
25 
26 #include "Epetra_Import.h"
27 #include "Epetra_Map.h"
28 #include "Epetra_MpiComm.h"
29 
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 namespace LinearAlgebra
34 {
35  namespace EpetraWrappers
36  {
38  :
39  vector(new Epetra_FEVector(Epetra_Map(0,0,0,Utilities::Trilinos::comm_self())))
40  {}
41 
42 
43 
45  :
46  vector(new Epetra_FEVector(V.trilinos_vector()))
47  {}
48 
49 
50 
51  Vector::Vector(const IndexSet &parallel_partitioner,
52  const MPI_Comm &communicator)
53  :
54  vector(new Epetra_FEVector(parallel_partitioner.make_trilinos_map(communicator,false)))
55  {}
56 
57 
58 
59  void Vector::reinit(const IndexSet &parallel_partitioner,
60  const MPI_Comm &communicator,
61  const bool omit_zeroing_entries)
62  {
63  Epetra_Map input_map = parallel_partitioner.make_trilinos_map(communicator,false);
64  if (vector->Map().SameAs(input_map)==false)
65  vector.reset(new Epetra_FEVector(input_map));
66  else if (omit_zeroing_entries==false)
67  {
68  const int ierr = vector->PutScalar(0.);
69  Assert(ierr==0, ExcTrilinosError(ierr));
70  (void) ierr;
71  }
72  }
73 
74 
75 
77  {
78  // Distinguish three cases:
79  // - First case: both vectors have the same layout.
80  // - Second case: both vectors have the same size but different layout.
81  // - Third case: the vectors have different size.
82  if (vector->Map().SameAs(V.trilinos_vector().Map()))
83  *vector = V.trilinos_vector();
84  else
85  {
86  if (size()==V.size())
87  {
88  Epetra_Import data_exchange(vector->Map(), V.trilinos_vector().Map());
89 
90  const int ierr = vector->Import(V.trilinos_vector(), data_exchange, Insert);
91  Assert(ierr==0, ExcTrilinosError(ierr));
92  (void) ierr;
93  }
94  else
95  vector.reset(new Epetra_FEVector(V.trilinos_vector()));
96  }
97 
98  return *this;
99  }
100 
101 
102 
104  VectorOperation::values operation,
105  std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern)
106  {
107  // If no communication pattern is given, create one. Otherwsie, use the
108  // one given.
109  if (communication_pattern == NULL)
110  {
111  // The first time import is called, a communication pattern is created.
112  // Check if the communication pattern already exists and if it can be
113  // reused.
117  {
119  dynamic_cast<const Epetra_MpiComm &>(vector->Comm()).Comm());
120  }
121  }
122  else
123  {
125  std_cxx11::dynamic_pointer_cast<const CommunicationPattern> (communication_pattern);
127  ExcMessage(std::string("The communication pattern is not of type ") +
128  "LinearAlgebra::EpetraWrappers::CommunicationPattern."));
129  }
130 
131  Epetra_Import import(epetra_comm_pattern->get_epetra_import());
132 
133  // The TargetMap and the SourceMap have their roles inverted.
134  Epetra_FEVector source_vector(import.TargetMap());
135  double *values = source_vector.Values();
136  std::copy(V.begin(), V.end(), values);
137 
138  if (operation==VectorOperation::insert)
139  vector->Export(source_vector, import, Insert);
140  else
141  vector->Export(source_vector, import, Add);
142  }
143 
144 
145 
146  Vector &Vector::operator*= (const double factor)
147  {
148  AssertIsFinite(factor);
149  vector->Scale(factor);
150 
151  return *this;
152  }
153 
154 
155 
156  Vector &Vector::operator/= (const double factor)
157  {
158  AssertIsFinite(factor);
159  Assert(factor!=0., ExcZero());
160  *this *= 1./factor;
161 
162  return *this;
163  }
164 
165 
166 
168  {
169  // Check that casting will work.
170  Assert(dynamic_cast<const Vector *>(&V)!=NULL, ExcVectorTypeNotCompatible());
171 
172  // Downcast V. If fails, throws an exception.
173  const Vector &down_V = dynamic_cast<const Vector &>(V);
174  // If the maps are the same we can Update right away.
175  if (vector->Map().SameAs(down_V.trilinos_vector().Map()))
176  {
177  const int ierr = vector->Update(1., down_V.trilinos_vector(), 1.);
178  Assert(ierr==0, ExcTrilinosError(ierr));
179  (void) ierr;
180  }
181  else
182  {
183  Assert(this->size()==down_V.size(),
184  ExcDimensionMismatch(this->size(), down_V.size()));
185 
186 #if DEAL_II_TRILINOS_VERSION_GTE(11,11,0)
187  Epetra_Import data_exchange (vector->Map(), down_V.trilinos_vector().Map());
188  const int ierr = vector->Import(down_V.trilinos_vector(),
189  data_exchange, Epetra_AddLocalAlso);
190  Assert(ierr==0, ExcTrilinosError(ierr));
191  (void) ierr;
192 #else
193  // In versions older than 11.11 the Import function is broken for adding
194  // Hence, we provide a workaround in this case
195 
196  Epetra_MultiVector dummy(vector->Map(), 1, false);
197  Epetra_Import data_exchange(dummy.Map(), down_V.trilinos_vector().Map());
198 
199  int ierr = dummy.Import(down_V.trilinos_vector(), data_exchange, Insert);
200  Assert(ierr==0, ExcTrilinosError(ierr));
201 
202  ierr = vector->Update(1.0, dummy, 1.0);
203  Assert(ierr==0, ExcTrilinosError(ierr));
204  (void) ierr;
205 #endif
206  }
207 
208  return *this;
209  }
210 
211 
212 
214  {
215  this->add(-1.,V);
216 
217  return *this;
218  }
219 
220 
221 
223  {
224  // Check that casting will work.
225  Assert(dynamic_cast<const Vector *>(&V)!=NULL,
227 
228  // Downcast V. If fails, throws an exception.
229  const Vector &down_V = dynamic_cast<const Vector &>(V);
230  Assert(this->size()==down_V.size(),
231  ExcDimensionMismatch(this->size(), down_V.size()));
232  Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
234 
235  double result(0.);
236  const int ierr = vector->Dot(down_V.trilinos_vector(), &result);
237  Assert(ierr==0, ExcTrilinosError(ierr));
238  (void) ierr;
239 
240  return result;
241  }
242 
243 
244 
245  void Vector::add(const double a)
246  {
247  AssertIsFinite(a);
248  const unsigned local_size(vector->MyLength());
249  for (unsigned int i=0; i<local_size; ++i)
250  (*vector)[0][i] += a;
251  }
252 
253 
254 
255  void Vector::add(const double a, const VectorSpaceVector<double> &V)
256  {
257  // Check that casting will work.
258  Assert(dynamic_cast<const Vector *>(&V)!=NULL,
260 
261  // Downcast V. If fails, throws an exception.
262  const Vector &down_V = dynamic_cast<const Vector &>(V);
263  AssertIsFinite(a);
264  Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
266 
267  const int ierr = vector->Update(a, down_V.trilinos_vector(), 1.);
268  Assert(ierr==0, ExcTrilinosError(ierr));
269  (void) ierr;
270  }
271 
272 
273 
274  void Vector::add(const double a, const VectorSpaceVector<double> &V,
275  const double b, const VectorSpaceVector<double> &W)
276  {
277  // Check that casting will work.
278  Assert(dynamic_cast<const Vector *>(&V)!=NULL,
280  // Check that casting will work.
281  Assert(dynamic_cast<const Vector *>(&W)!=NULL,
283 
284  // Downcast V. If fails, throws an exception.
285  const Vector &down_V = dynamic_cast<const Vector &>(V);
286  // Downcast W. If fails, throws an exception.
287  const Vector &down_W = dynamic_cast<const Vector &>(W);
288  Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
290  Assert(vector->Map().SameAs(down_W.trilinos_vector().Map()),
292  AssertIsFinite(a);
293  AssertIsFinite(b);
294 
295  const int ierr = vector->Update(a, down_V.trilinos_vector(), b,
296  down_W.trilinos_vector(), 1.);
297  Assert(ierr==0, ExcTrilinosError(ierr));
298  (void) ierr;
299  }
300 
301 
302 
303  void Vector::sadd(const double s, const double a,
304  const VectorSpaceVector<double> &V)
305  {
306  // Check that casting will work.
307  Assert(dynamic_cast<const Vector *>(&V)!=NULL,
309 
310  *this *= s;
311  // Downcast V. It fails, throws an exception.
312  const Vector &down_V = dynamic_cast<const Vector &>(V);
313  Vector tmp(down_V);
314  tmp *= a;
315  *this += tmp;
316  }
317 
318 
319 
320  void Vector::scale(const VectorSpaceVector<double> &scaling_factors)
321  {
322  // Check that casting will work.
323  Assert(dynamic_cast<const Vector *>(&scaling_factors)!=NULL,
325 
326  // Downcast scaling_factors. If fails, throws an exception.
327  const Vector &down_scaling_factors =
328  dynamic_cast<const Vector &>(scaling_factors);
329  Assert(vector->Map().SameAs(down_scaling_factors.trilinos_vector().Map()),
331 
332  const int ierr = vector->Multiply(1.0, down_scaling_factors.trilinos_vector(),
333  *vector, 0.0);
334  Assert(ierr==0, ExcTrilinosError(ierr));
335  (void) ierr;
336  }
337 
338 
339 
340  void Vector::equ(const double a, const VectorSpaceVector<double> &V)
341  {
342  // Check that casting will work.
343  Assert(dynamic_cast<const Vector *>(&V)!=NULL,
345 
346  // Downcast V. If fails, throws an exception.
347  const Vector &down_V = dynamic_cast<const Vector &>(V);
348  // If we don't have the same map, copy.
349  if (vector->Map().SameAs(down_V.trilinos_vector().Map())==false)
350  this->sadd(0., a, V);
351  else
352  {
353  // Otherwise, just update
354  int ierr = vector->Update(a, down_V.trilinos_vector(), 0.);
355  Assert(ierr==0, ExcTrilinosError(ierr));
356  (void) ierr;
357  }
358  }
359 
360 
361 
362  double Vector::mean_value() const
363  {
364  double mean_value(0.);
365 
366  int ierr = vector->MeanValue(&mean_value);
367  Assert(ierr==0, ExcTrilinosError(ierr));
368  (void) ierr;
369 
370  return mean_value;
371  }
372 
373 
374 
375  double Vector::l1_norm() const
376  {
377  double norm(0.);
378  int ierr = vector->Norm1(&norm);
379  Assert(ierr==0, ExcTrilinosError(ierr));
380  (void) ierr;
381 
382  return norm;
383  }
384 
385 
386 
387  double Vector::l2_norm() const
388  {
389  double norm(0.);
390  int ierr = vector->Norm2(&norm);
391  Assert(ierr==0, ExcTrilinosError(ierr));
392  (void) ierr;
393 
394  return norm;
395  }
396 
397 
398 
399  double Vector::linfty_norm() const
400  {
401  double norm(0.);
402  int ierr = vector->NormInf(&norm);
403  Assert(ierr==0, ExcTrilinosError(ierr));
404  (void) ierr;
405 
406  return norm;
407  }
408 
409 
410 
411  double Vector::add_and_dot(const double a,
412  const VectorSpaceVector<double> &V,
413  const VectorSpaceVector<double> &W)
414  {
415  this->add(a, V);
416 
417  return *this * W;
418  }
419 
420 
421 
422  Vector::size_type Vector::size() const
423  {
424 #ifndef DEAL_II_WITH_64BIT_INDICES
425  return vector->GlobalLength();
426 #else
427  return vector->GlobalLength64();
428 #endif
429  }
430 
431 
432 
434  {
435  const Epetra_MpiComm *epetra_comm
436  = dynamic_cast<const Epetra_MpiComm *>(&(vector->Comm()));
437  Assert (epetra_comm != 0, ExcInternalError());
438  return epetra_comm->GetMpiComm();
439  }
440 
441 
442 
444  {
445  IndexSet is (size());
446 
447  // easy case: local range is contiguous
448  if (vector->Map().LinearMap())
449  {
450 #ifndef DEAL_II_WITH_64BIT_INDICES
451  is.add_range(vector->Map().MinMyGID(), vector->Map().MaxMyGID()+1);
452 #else
453  is.add_range(vector->Map().MinMyGID64(), vector->Map().MaxMyGID64()+1);
454 #endif
455  }
456  else if (vector->Map().NumMyElements() > 0)
457  {
458  const size_type n_indices = vector->Map().NumMyElements();
459 #ifndef DEAL_II_WITH_64BIT_INDICES
460  unsigned int *vector_indices = (unsigned int *)vector->Map().MyGlobalElements();
461 #else
462  size_type *vector_indices = (size_type *)vector->Map().MyGlobalElements64();
463 #endif
464  is.add_indices(vector_indices, vector_indices+n_indices);
465  }
466  is.compress();
467 
468  return is;
469  }
470 
471 
472 
473  const Epetra_FEVector &Vector::trilinos_vector() const
474  {
475  return *vector;
476  }
477 
478 
479 
480  Epetra_FEVector &Vector::trilinos_vector()
481  {
482  return *vector;
483  }
484 
485 
486 
487  void Vector::print(std::ostream &out,
488  const unsigned int precision,
489  const bool scientific,
490  const bool across) const
491  {
492  AssertThrow(out, ExcIO());
493  boost::io::ios_flags_saver restore_flags(out);
494 
495  // Get a representation of the vector and loop over all
496  // the elements
497  double *val;
498  int leading_dimension;
499  int ierr = vector->ExtractView(&val, &leading_dimension);
500 
501  Assert(ierr==0, ExcTrilinosError(ierr));
502  (void) ierr;
503  out.precision (precision);
504  if (scientific)
505  out.setf(std::ios::scientific, std::ios::floatfield);
506  else
507  out.setf(std::ios::fixed, std::ios::floatfield);
508 
509  if (across)
510  for (int i=0; i<vector->MyLength(); ++i)
511  out << val[i] << ' ';
512  else
513  for (int i=0; i<vector->MyLength(); ++i)
514  out << val[i] << std::endl;
515  out << std::endl;
516 
517  // restore the representation
518  // of the vector
519  AssertThrow(out, ExcIO());
520  }
521 
522 
523 
524  std::size_t Vector::memory_consumption() const
525  {
526  return sizeof(*this)
527  + vector->MyLength()*(sizeof(double)+
528  sizeof(TrilinosWrappers::types::int_type));
529  }
530 
531 
532 
533  void Vector::create_epetra_comm_pattern(const IndexSet &source_index_set,
534  const MPI_Comm &mpi_comm)
535  {
536  source_stored_elements = source_index_set;
538  source_index_set, mpi_comm));
539  }
540  }
541 }
542 
543 DEAL_II_NAMESPACE_CLOSE
544 
545 #endif
546 
547 #endif
static ::ExceptionBase & ExcDifferentParallelPartitioning()
static ::ExceptionBase & ExcIO()
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1463
virtual Vector & operator+=(const VectorSpaceVector< double > &V)
virtual std::size_t memory_consumption() const
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual void scale(const VectorSpaceVector< double > &scaling_factors)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:535
size_type size() const
Definition: index_set.h:1419
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm)
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
const Epetra_FEVector & trilinos_vector() const
virtual double operator*(const VectorSpaceVector< double > &V) const
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcVectorTypeNotCompatible()
void reinit(const IndexSet &parallel_partitioner, const MPI_Comm &communicator, const bool omit_zeroing_entries=false)
virtual ::IndexSet locally_owned_elements() const
void add_range(const size_type begin, const size_type end)
Definition: index_set.cc:88
virtual double add_and_dot(const double a, const VectorSpaceVector< double > &V, const VectorSpaceVector< double > &W)
Definition: mpi.h:48
void compress() const
Definition: index_set.h:1428
std_cxx11::shared_ptr< Epetra_FEVector > vector
virtual void equ(const double a, const VectorSpaceVector< double > &V)
virtual void sadd(const double s, const double a, const VectorSpaceVector< double > &V)
virtual Vector & operator/=(const double factor)
virtual Vector & operator-=(const VectorSpaceVector< double > &V)
static ::ExceptionBase & ExcZero()
#define AssertIsFinite(number)
Definition: exceptions.h:1197
virtual Vector & operator*=(const double factor)
const IndexSet & get_stored_elements() const
std_cxx11::shared_ptr< const CommunicationPattern > epetra_comm_pattern
virtual void import(const ReadWriteVector< double > &V, VectorOperation::values operation, std_cxx11::shared_ptr< const CommunicationPatternBase > communication_pattern=std_cxx11::shared_ptr< const CommunicationPatternBase >())
static ::ExceptionBase & ExcInternalError()