Reference documentation for deal.II version 8.4.2
trilinos_vector_base.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2015 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/base/memory_consumption.h>
17 #include <deal.II/lac/trilinos_vector_base.h>
18 
19 #ifdef DEAL_II_WITH_TRILINOS
20 
21 # include <cmath>
22 
23 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
24 # include <Epetra_Import.h>
25 # include <Epetra_Export.h>
27 
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 namespace TrilinosWrappers
32 {
33  namespace
34  {
35 #ifndef DEAL_II_WITH_64BIT_INDICES
36  // define a helper function that queries the global vector length of an
37  // Epetra_FEVector object by calling either the 32- or 64-bit
38  // function necessary.
39  int global_length(const Epetra_FEVector &vector)
40  {
41  return vector.GlobalLength();
42  }
43 #else
44  // define a helper function that queries the global vector length of an
45  // Epetra_FEVector object by calling either the 32- or 64-bit
46  // function necessary.
47  long long int global_length(const Epetra_FEVector &vector)
48  {
49  return vector.GlobalLength64();
50  }
51 #endif
52  }
53 
54  namespace internal
55  {
56  VectorReference::operator TrilinosScalar () const
57  {
58  Assert (index < vector.size(),
59  ExcIndexRange (index, 0, vector.size()));
60 
61  // Trilinos allows for vectors to be referenced by the [] or ()
62  // operators but only () checks index bounds. We check these bounds by
63  // ourselves, so we can use []. Note that we can only get local values.
64 
65  const TrilinosWrappers::types::int_type local_index =
66  vector.vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
67  Assert (local_index >= 0,
68  VectorBase::ExcAccessToNonLocalElement (index, vector.local_size(),
69  vector.vector->Map().MinMyGID(),
70  vector.vector->Map().MaxMyGID()));
71 
72 
73  return (*(vector.vector))[0][local_index];
74  }
75  }
76 
77 
78 
80  :
81  last_action (Zero),
82  compressed (true),
83  has_ghosts (false),
84 #ifdef DEAL_II_WITH_MPI
85  vector(new Epetra_FEVector(
86  Epetra_Map(0,0,Epetra_MpiComm(MPI_COMM_SELF))))
87 #else
88  vector(new Epetra_FEVector(
89  Epetra_Map(0,0,Epetra_SerialComm())))
90 #endif
91  {}
92 
93 
94 
96  :
97  Subscriptor(),
98  last_action (Zero),
99  compressed (true),
101  vector(new Epetra_FEVector(*v.vector))
102  {}
103 
104 
105 
107  {}
108 
109 
110 
111  void
113  {
114  // When we clear the vector, reset the pointer and generate an empty
115  // vector.
116 #ifdef DEAL_II_WITH_MPI
117  Epetra_Map map (0, 0, Epetra_MpiComm(MPI_COMM_SELF));
118 #else
119  Epetra_Map map (0, 0, Epetra_SerialComm());
120 #endif
121 
122  has_ghosts = false;
123  vector.reset (new Epetra_FEVector(map));
124  last_action = Zero;
125  }
126 
127 
128 
129  VectorBase &
131  {
132  Assert (vector.get() != 0,
133  ExcMessage("Vector is not constructed properly."));
134 
135  if (local_range() != v.local_range())
136  {
137  last_action = Zero;
138  vector.reset (new Epetra_FEVector(*v.vector));
140  }
141  else
142  {
143  Assert (vector->Map().SameAs(v.vector->Map()) == true,
144  ExcMessage ("The Epetra maps in the assignment operator ="
145  " do not match, even though the local_range "
146  " seems to be the same. Check vector setup!"));
147  int ierr;
148  ierr = vector->GlobalAssemble(last_action);
149  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
150 
151  ierr = vector->Update(1.0, *v.vector, 0.0);
152  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
153 
154  last_action = Zero;
155  }
156 
157  return *this;
158  }
159 
160 
161 
162  template <typename number>
163  VectorBase &
164  VectorBase::operator = (const ::Vector<number> &v)
165  {
166  Assert (size() == v.size(),
167  ExcDimensionMismatch(size(), v.size()));
168 
169  // this is probably not very efficient
170  // but works. in particular, we could do
171  // better if we know that
172  // number==TrilinosScalar because then we
173  // could elide the copying of elements
174  //
175  // let's hope this isn't a
176  // particularly frequent operation
177  std::pair<size_type, size_type>
178  local_range = this->local_range ();
179  for (size_type i=local_range.first; i<local_range.second; ++i)
180  (*vector)[0][i-local_range.first] = v(i);
181 
182  return *this;
183  }
184 
185 
186 
187  void
188  VectorBase::compress (::VectorOperation::values given_last_action)
189  {
190  //Select which mode to send to Trilinos. Note that we use last_action if
191  //available and ignore what the user tells us to detect wrongly mixed
192  //operations. Typically given_last_action is only used on machines that do
193  //not execute an operation (because they have no own cells for example).
194  Epetra_CombineMode mode = last_action;
195  if (last_action == Zero)
196  {
197  if (given_last_action==::VectorOperation::add)
198  mode = Add;
199  else if (given_last_action==::VectorOperation::insert)
200  mode = Insert;
201  }
202 
203 #ifdef DEBUG
204 # ifdef DEAL_II_WITH_MPI
205  // check that every process has decided to use the same mode. This will
206  // otherwise result in undefined behaviour in the call to
207  // GlobalAssemble().
208  double double_mode = mode;
210  = Utilities::MPI::min_max_avg (double_mode,
211  dynamic_cast<const Epetra_MpiComm *>
212  (&vector_partitioner().Comm())->GetMpiComm());
213  Assert(result.max-result.min<1e-5,
214  ExcMessage ("Not all processors agree whether the last operation on "
215  "this vector was an addition or a set operation. This will "
216  "prevent the compress() operation from succeeding."));
217 
218 # endif
219 #endif
220 
221  // Now pass over the information about what we did last to the vector.
222  int ierr = 0;
223  if (nonlocal_vector.get() == 0 || mode != Add)
224  ierr = vector->GlobalAssemble(mode);
225  else
226  {
227  Epetra_Export exporter(nonlocal_vector->Map(), vector->Map());
228  ierr = vector->Export(*nonlocal_vector, exporter, mode);
229  nonlocal_vector->PutScalar(0.);
230  }
231  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
232  last_action = Zero;
233 
234  compressed = true;
235  }
236 
237 
238 
239  TrilinosScalar
240  VectorBase::el (const size_type index) const
241  {
242  // Extract local indices in the vector.
243  TrilinosWrappers::types::int_type trilinos_i =
244  vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
245 
246  // If the element is not present on the current processor, we can't
247  // continue. Just print out 0 as opposed to the () method below.
248  if (trilinos_i == -1)
249  return 0.;
250  else
251  return (*vector)[0][trilinos_i];
252  }
253 
254 
255 
256  TrilinosScalar
257  VectorBase::operator () (const size_type index) const
258  {
259  // Extract local indices in the vector.
260  TrilinosWrappers::types::int_type trilinos_i =
261  vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
262  TrilinosScalar value = 0.;
263 
264  // If the element is not present on the current processor, we can't
265  // continue. This is the main difference to the el() function.
266  if (trilinos_i == -1)
267  {
268  Assert (false, ExcAccessToNonLocalElement(index, local_size(),
269  vector->Map().MinMyGID(),
270  vector->Map().MaxMyGID()));
271  }
272  else
273  value = (*vector)[0][trilinos_i];
274 
275  return value;
276  }
277 
278 
279 
280  void
282  const bool allow_different_maps)
283  {
284  if (allow_different_maps == false)
285  *this += v;
286  else
287  {
288  Assert (!has_ghost_elements(), ExcGhostsPresent());
289  AssertThrow (size() == v.size(),
290  ExcDimensionMismatch (size(), v.size()));
291 
292 #if DEAL_II_TRILINOS_VERSION_GTE(11,11,0)
293  Epetra_Import data_exchange (vector->Map(), v.vector->Map());
294  int ierr = vector->Import(*v.vector, data_exchange, Epetra_AddLocalAlso);
295  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
296  last_action = Add;
297 #else
298  // In versions older than 11.11 the Import function is broken for adding
299  // Hence, we provide a workaround in this case
300 
301  Epetra_MultiVector dummy(vector->Map(), 1, false);
302  Epetra_Import data_exchange (dummy.Map(), v.vector->Map());
303 
304  int ierr = dummy.Import(*v.vector, data_exchange, Insert);
305  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
306 
307  ierr = vector->Update (1.0, dummy, 1.0);
308  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
309 #endif
310  }
311  }
312 
313 
314 
315  bool
317  {
318  Assert (size() == v.size(),
319  ExcDimensionMismatch(size(), v.size()));
320  if (local_size() != v.local_size())
321  return false;
322 
323  size_type i;
324  for (i=0; i<local_size(); i++)
325  if ((*(v.vector))[0][i]!=(*vector)[0][i]) return false;
326 
327  return true;
328  }
329 
330 
331 
332  bool
334  {
335  Assert (size() == v.size(),
336  ExcDimensionMismatch(size(), v.size()));
337 
338  return (!(*this==v));
339  }
340 
341 
342 
343  bool
345  {
346  // get a representation of the vector and
347  // loop over all the elements
348  TrilinosScalar *start_ptr = (*vector)[0];
349  const TrilinosScalar *ptr = start_ptr,
350  *eptr = start_ptr + local_size();
351  unsigned int flag = 0;
352  while (ptr != eptr)
353  {
354  if (*ptr != 0)
355  {
356  flag = 1;
357  break;
358  }
359  ++ptr;
360  }
361 
362 #ifdef DEAL_II_WITH_MPI
363  // in parallel, check that the vector
364  // is zero on _all_ processors.
365  const Epetra_MpiComm *mpi_comm
366  = dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
367  unsigned int num_nonzero = Utilities::MPI::sum(flag, mpi_comm->Comm());
368  return num_nonzero == 0;
369 #else
370  return flag == 0;
371 #endif
372 
373  }
374 
375 
376 
377  bool
379  {
380 #ifdef DEAL_II_WITH_MPI
381  // if this vector is a parallel one, then
382  // we need to communicate to determine
383  // the answer to the current
384  // function. this still has to be
385  // implemented
386  AssertThrow(local_size() == size(), ExcNotImplemented());
387 #endif
388  // get a representation of the vector and
389  // loop over all the elements
390  TrilinosScalar *start_ptr;
391  int leading_dimension;
392  int ierr = vector->ExtractView (&start_ptr, &leading_dimension);
393  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
394 
395  // TODO: This
396  // won't work in parallel like
397  // this. Find out a better way to
398  // this in that case.
399  const TrilinosScalar *ptr = start_ptr,
400  *eptr = start_ptr + size();
401  bool flag = true;
402  while (ptr != eptr)
403  {
404  if (*ptr < 0.0)
405  {
406  flag = false;
407  break;
408  }
409  ++ptr;
410  }
411 
412  return flag;
413  }
414 
415 
416 
417  void
418  VectorBase::equ (const TrilinosScalar a,
419  const VectorBase &v,
420  const TrilinosScalar b,
421  const VectorBase &w)
422  {
423  // if we have ghost values, do not allow
424  // writing to this vector at all.
425  Assert (!has_ghost_elements(), ExcGhostsPresent());
426  Assert (v.local_size() == w.local_size(),
427  ExcDimensionMismatch (v.local_size(), w.local_size()));
428 
429  AssertIsFinite(a);
430  AssertIsFinite(b);
431 
432  // If we don't have the same map, copy.
433  if (vector->Map().SameAs(v.vector->Map())==false)
434  {
435  sadd(0., a, v, b, w);
436  }
437  else
438  {
439  // Otherwise, just update. verify
440  // that *this does not only have
441  // the same map as v (the
442  // if-condition above) but also as
443  // w
444  Assert (vector->Map().SameAs(w.vector->Map()),
445  ExcDifferentParallelPartitioning());
446  int ierr = vector->Update(a, *v.vector, b, *w.vector, 0.0);
447  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
448 
449  last_action = Zero;
450  }
451  }
452 
453 
454 
455  // TODO: up to now only local
456  // data printed out! Find a
457  // way to neatly output
458  // distributed data...
459  void
460  VectorBase::print (const char *format) const
461  {
462  Assert (global_length(*vector)!=0, ExcEmptyObject());
463  (void)global_length;
464 
465  for (size_type j=0; j<size(); ++j)
466  {
467  double t = (*vector)[0][j];
468 
469  if (format != 0)
470  std::printf (format, t);
471  else
472  std::printf (" %5.2f", double(t));
473  }
474  std::printf ("\n");
475  }
476 
477 
478 
479  void
480  VectorBase::print (std::ostream &out,
481  const unsigned int precision,
482  const bool scientific,
483  const bool across) const
484  {
485  AssertThrow (out, ExcIO());
486 
487  // get a representation of the
488  // vector and loop over all
489  // the elements TODO: up to
490  // now only local data printed
491  // out! Find a way to neatly
492  // output distributed data...
493  TrilinosScalar *val;
494  int leading_dimension;
495  int ierr = vector->ExtractView (&val, &leading_dimension);
496 
497  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
498  out.precision (precision);
499  if (scientific)
500  out.setf (std::ios::scientific, std::ios::floatfield);
501  else
502  out.setf (std::ios::fixed, std::ios::floatfield);
503 
504  if (across)
505  for (size_type i=0; i<size(); ++i)
506  out << static_cast<double>(val[i]) << ' ';
507  else
508  for (size_type i=0; i<size(); ++i)
509  out << static_cast<double>(val[i]) << std::endl;
510  out << std::endl;
511 
512  // restore the representation
513  // of the vector
514  AssertThrow (out, ExcIO());
515  }
516 
517 
518 
519  void
521  {
522  std::swap(last_action, v.last_action);
523  std::swap(compressed, v.compressed);
524  std::swap(vector, v.vector);
525  }
526 
527 
528 
529  std::size_t
531  {
532  //TODO[TH]: No accurate memory
533  //consumption for Trilinos vectors
534  //yet. This is a rough approximation with
535  //one index and the value per local
536  //entry.
537  return sizeof(*this)
538  + this->local_size()*( sizeof(double)+
539  sizeof(TrilinosWrappers::types::int_type) );
540  }
541 
542 } /* end of namespace TrilinosWrappers */
543 
544 
545 namespace TrilinosWrappers
546 {
547 #include "trilinos_vector_base.inst"
548 }
549 
550 DEAL_II_NAMESPACE_CLOSE
551 
552 #endif // DEAL_II_WITH_TRILINOS
TrilinosScalar el(const size_type index) const DEAL_II_DEPRECATED
reference operator()(const size_type index)
void sadd(const TrilinosScalar s, const VectorBase &V)
::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_MultiVector > nonlocal_vector
std::pair< size_type, size_type > local_range() const
void print(const char *format=0) const DEAL_II_DEPRECATED
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
T sum(const T &t, const MPI_Comm &mpi_communicator)
Definition: mpi.h:611
std_cxx11::shared_ptr< Epetra_FEVector > vector
#define Assert(cond, exc)
Definition: exceptions.h:294
const Epetra_Map & vector_partitioner() const
std::size_t memory_consumption() const
void compress(::VectorOperation::values operation)
VectorBase & operator=(const TrilinosScalar s)
bool operator!=(const VectorBase &v) const
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:239
bool operator==(const VectorBase &v) const
#define AssertIsFinite(number)
Definition: exceptions.h:1096
size_type local_size() const
void equ(const TrilinosScalar a, const VectorBase &V)