Reference documentation for deal.II version GIT 112f7bbc69 2023-05-28 21:25:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
trilinos_vector.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2022 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/base/mpi.h>
22 
27 
28 # include <boost/io/ios_state.hpp>
29 
30 # include <Epetra_Export.h>
31 # include <Epetra_Import.h>
32 # include <Epetra_Vector.h>
33 
34 # include <cmath>
35 # include <memory>
36 
37 
39 
40 namespace TrilinosWrappers
41 {
42 # ifndef DOXYGEN
43  namespace internal
44  {
45  VectorReference::operator TrilinosScalar() const
46  {
47  AssertIndexRange(index, vector.size());
48 
49  // Trilinos allows for vectors to be referenced by the [] or ()
50  // operators but only () checks index bounds. We check these bounds by
51  // ourselves, so we can use []. Note that we can only get local values.
52 
53  const TrilinosWrappers::types::int_type local_index =
54  vector.vector->Map().LID(
55  static_cast<TrilinosWrappers::types::int_type>(index));
56  Assert(local_index >= 0,
58  index,
59  vector.local_size(),
60  vector.vector->Map().MinMyGID(),
61  vector.vector->Map().MaxMyGID()));
62 
63 
64  return (*(vector.vector))[0][local_index];
65  }
66  } // namespace internal
67 # endif
68 
69  namespace MPI
70  {
72  : Subscriptor()
73  , last_action(Zero)
74  , compressed(true)
75  , has_ghosts(false)
76  , vector(new Epetra_FEVector(
77  Epetra_Map(0, 0, 0, Utilities::Trilinos::comm_self())))
78  {}
79 
80 
81 
82  Vector::Vector(const IndexSet &parallel_partitioning,
83  const MPI_Comm communicator)
84  : Vector()
85  {
86  reinit(parallel_partitioning, communicator);
87  }
88 
89 
90 
92  : Vector()
93  {
95  vector = std::make_unique<Epetra_FEVector>(*v.vector);
97  }
98 
99 
100 
101  Vector::Vector(Vector &&v) // NOLINT
102  : Vector()
103  {
104  // initialize a minimal, valid object and swap
105  static_cast<Subscriptor &>(*this) = static_cast<Subscriptor &&>(v);
106  swap(v);
107  }
108 
109 
110 
111  Vector::Vector(const IndexSet &parallel_partitioner,
112  const Vector & v,
113  const MPI_Comm communicator)
114  : Vector()
115  {
116  AssertThrow(parallel_partitioner.size() ==
117  static_cast<size_type>(
119  ExcDimensionMismatch(parallel_partitioner.size(),
121  v.vector->Map())));
122 
123  vector = std::make_unique<Epetra_FEVector>(
124  parallel_partitioner.make_trilinos_map(communicator, true));
125  reinit(v, false, true);
126  }
127 
128 
129 
130  Vector::Vector(const IndexSet &local,
131  const IndexSet &ghost,
132  const MPI_Comm communicator)
133  : Vector()
134  {
135  reinit(local, ghost, communicator, false);
136  }
137 
138 
139 
140  void
142  {
143  // When we clear the vector, reset the pointer and generate an empty
144  // vector.
145  Epetra_Map map(0, 0, Epetra_MpiComm(MPI_COMM_SELF));
146 
147  has_ghosts = false;
148  vector = std::make_unique<Epetra_FEVector>(map);
149  last_action = Zero;
150  }
151 
152 
153 
154  void
155  Vector::reinit(const IndexSet &parallel_partitioner,
156  const MPI_Comm communicator,
157  const bool /*omit_zeroing_entries*/)
158  {
159  nonlocal_vector.reset();
160 
161  const bool overlapping =
162  !parallel_partitioner.is_ascending_and_one_to_one(communicator);
163 
164  Epetra_Map map =
165  parallel_partitioner.make_trilinos_map(communicator, overlapping);
166 
167  vector = std::make_unique<Epetra_FEVector>(map);
168 
169  has_ghosts = vector->Map().UniqueGIDs() == false;
170 
171  // If the IndexSets are overlapping, we don't really know
172  // which process owns what. So we decide that no process
173  // owns anything in that case. In particular asking for
174  // the locally owned elements is not allowed.
175  if (has_ghosts)
176  {
179  }
180  else
181  owned_elements = parallel_partitioner;
182 
183 # ifdef DEBUG
184  const size_type n_elements_global =
186 
187  Assert(has_ghosts || n_elements_global == size(), ExcInternalError());
188 # endif
189 
190  last_action = Zero;
191  }
192 
193 
194 
195  void
197  const bool omit_zeroing_entries,
198  const bool allow_different_maps)
199  {
200  nonlocal_vector.reset();
201 
202  // In case we do not allow to have different maps, this call means that
203  // we have to reset the vector. So clear the vector, initialize our map
204  // with the map in v, and generate the vector.
205  if (allow_different_maps == false)
206  {
207  // check equality for MPI communicators: We can only choose the fast
208  // version in case the underlying Epetra_MpiComm object is the same,
209  // otherwise we might access an MPI_Comm object that has been
210  // deleted
211  const Epetra_MpiComm *my_comm =
212  dynamic_cast<const Epetra_MpiComm *>(&vector->Comm());
213  const Epetra_MpiComm *v_comm =
214  dynamic_cast<const Epetra_MpiComm *>(&v.vector->Comm());
215  const bool same_communicators =
216  my_comm != nullptr && v_comm != nullptr &&
217  my_comm->DataPtr() == v_comm->DataPtr();
218  if (!same_communicators ||
219  vector->Map().SameAs(v.vector->Map()) == false)
220  {
221  vector = std::make_unique<Epetra_FEVector>(v.vector->Map());
223  last_action = Zero;
225  }
226  else if (omit_zeroing_entries == false)
227  {
228  // old and new vectors have exactly the same map, i.e. size and
229  // parallel distribution
230  int ierr;
231  ierr = vector->GlobalAssemble(last_action);
232  (void)ierr;
233  Assert(ierr == 0, ExcTrilinosError(ierr));
234 
235  ierr = vector->PutScalar(0.0);
236  Assert(ierr == 0, ExcTrilinosError(ierr));
237 
238  last_action = Zero;
239  }
240  }
241 
242  // Otherwise, we have to check that the two vectors are already of the
243  // same size, create an object for the data exchange and then insert all
244  // the data. The first assertion is only a check whether the user knows
245  // what they are doing.
246  else
247  {
248  Assert(omit_zeroing_entries == false,
249  ExcMessage(
250  "It is not possible to exchange data with the "
251  "option 'omit_zeroing_entries' set, which would not write "
252  "elements."));
253 
254  AssertThrow(size() == v.size(),
255  ExcDimensionMismatch(size(), v.size()));
256 
257  Epetra_Import data_exchange(vector->Map(), v.vector->Map());
258 
259  const int ierr = vector->Import(*v.vector, data_exchange, Insert);
260  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
261 
262  last_action = Insert;
263  }
264 # ifdef DEBUG
265  const Epetra_MpiComm *comm_ptr =
266  dynamic_cast<const Epetra_MpiComm *>(&(v.vector->Comm()));
267  Assert(comm_ptr != nullptr, ExcInternalError());
268  const size_type n_elements_global =
269  Utilities::MPI::sum(owned_elements.n_elements(), comm_ptr->Comm());
270  Assert(has_ghosts || n_elements_global == size(), ExcInternalError());
271 # endif
272  }
273 
274 
275 
276  void
277  Vector::reinit(const MPI::BlockVector &v, const bool import_data)
278  {
279  nonlocal_vector.reset();
282 
283  // In case we do not allow to have different maps, this call means that
284  // we have to reset the vector. So clear the vector, initialize our map
285  // with the map in v, and generate the vector.
286  if (v.n_blocks() == 0)
287  return;
288 
289  // create a vector that holds all the elements contained in the block
290  // vector. need to manually create an Epetra_Map.
291  size_type n_elements = 0, added_elements = 0, block_offset = 0;
292  for (size_type block = 0; block < v.n_blocks(); ++block)
293  n_elements += v.block(block).local_size();
294  std::vector<TrilinosWrappers::types::int_type> global_ids(n_elements, -1);
295  for (size_type block = 0; block < v.n_blocks(); ++block)
296  {
297  TrilinosWrappers::types::int_type *glob_elements =
299  v.block(block).trilinos_partitioner());
300  for (size_type i = 0; i < v.block(block).local_size(); ++i)
301  global_ids[added_elements++] = glob_elements[i] + block_offset;
302  owned_elements.add_indices(v.block(block).owned_elements,
303  block_offset);
304  block_offset += v.block(block).size();
305  }
306 
307  Assert(n_elements == added_elements, ExcInternalError());
308  Epetra_Map new_map(v.size(),
309  n_elements,
310  global_ids.data(),
311  0,
312  v.block(0).trilinos_partitioner().Comm());
313 
314  auto actual_vec = std::make_unique<Epetra_FEVector>(new_map);
315 
316  TrilinosScalar *entries = (*actual_vec)[0];
317  for (size_type block = 0; block < v.n_blocks(); ++block)
318  {
319  v.block(block).trilinos_vector().ExtractCopy(entries, 0);
320  entries += v.block(block).local_size();
321  }
322 
323  if (import_data == true)
324  {
326  *actual_vec)) == v.size(),
328  *actual_vec),
329  v.size()));
330 
331  Epetra_Import data_exchange(vector->Map(), actual_vec->Map());
332 
333  const int ierr = vector->Import(*actual_vec, data_exchange, Insert);
334  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
335 
336  last_action = Insert;
337  }
338  else
339  vector = std::move(actual_vec);
340 # ifdef DEBUG
341  const Epetra_MpiComm *comm_ptr =
342  dynamic_cast<const Epetra_MpiComm *>(&(vector->Comm()));
343  Assert(comm_ptr != nullptr, ExcInternalError());
344  const size_type n_elements_global =
345  Utilities::MPI::sum(owned_elements.n_elements(), comm_ptr->Comm());
346 
347  Assert(has_ghosts || n_elements_global == size(), ExcInternalError());
348 # endif
349  }
350 
351 
352 
353  void
354  Vector::reinit(const IndexSet &locally_owned_entries,
355  const IndexSet &ghost_entries,
356  const MPI_Comm communicator,
357  const bool vector_writable)
358  {
359  nonlocal_vector.reset();
360  owned_elements = locally_owned_entries;
361  if (vector_writable == false)
362  {
363  IndexSet parallel_partitioner = locally_owned_entries;
364  parallel_partitioner.add_indices(ghost_entries);
365  Epetra_Map map =
366  parallel_partitioner.make_trilinos_map(communicator, true);
367  vector = std::make_unique<Epetra_FEVector>(map);
368  }
369  else
370  {
371  Epetra_Map map =
372  locally_owned_entries.make_trilinos_map(communicator, true);
373  Assert(map.IsOneToOne(),
374  ExcMessage("A writable vector must not have ghost entries in "
375  "its parallel partitioning"));
376 
377  if (vector->Map().SameAs(map) == false)
378  vector = std::make_unique<Epetra_FEVector>(map);
379  else
380  {
381  const int ierr = vector->PutScalar(0.);
382  (void)ierr;
383  Assert(ierr == 0, ExcTrilinosError(ierr));
384  }
385 
386  IndexSet nonlocal_entries(ghost_entries);
387  nonlocal_entries.subtract_set(locally_owned_entries);
388  if (Utilities::MPI::n_mpi_processes(communicator) > 1)
389  {
390  Epetra_Map nonlocal_map =
391  nonlocal_entries.make_trilinos_map(communicator, true);
393  std::make_unique<Epetra_MultiVector>(nonlocal_map, 1);
394  }
395  }
396 
397  has_ghosts = vector->Map().UniqueGIDs() == false;
398 
399  last_action = Zero;
400 
401 # ifdef DEBUG
402  const size_type n_elements_global =
404 
405  Assert(has_ghosts || n_elements_global == size(), ExcInternalError());
406 # endif
407  }
408 
409 
410 
411  void
413  const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
414  const bool vector_writable)
415  {
416  this->reinit(partitioner->locally_owned_range(),
417  partitioner->ghost_indices(),
418  partitioner->get_mpi_communicator(),
419  vector_writable);
420  }
421 
422 
423 
424  Vector &
426  {
427  Assert(vector.get() != nullptr,
428  ExcMessage("Vector is not constructed properly."));
429 
430  // check equality for MPI communicators to avoid accessing a possibly
431  // invalid MPI_Comm object
432  const Epetra_MpiComm *my_comm =
433  dynamic_cast<const Epetra_MpiComm *>(&vector->Comm());
434  const Epetra_MpiComm *v_comm =
435  dynamic_cast<const Epetra_MpiComm *>(&v.vector->Comm());
436  const bool same_communicators = my_comm != nullptr && v_comm != nullptr &&
437  my_comm->DataPtr() == v_comm->DataPtr();
438  // Need to ask MPI whether the communicators are the same. We would like
439  // to use the following checks but currently we cannot make sure the
440  // memory of my_comm is not stale from some MPI_Comm_free
441  // somewhere. This can happen when a vector lives in GrowingVectorMemory
442  // data structures. Thus, the following code is commented out.
443  //
444  // if (my_comm != nullptr &&
445  // v_comm != nullptr &&
446  // my_comm->DataPtr() != v_comm->DataPtr())
447  // {
448  // int communicators_same = 0;
449  // const int ierr = MPI_Comm_compare (my_comm->GetMpiComm(),
450  // v_comm->GetMpiComm(),
451  // &communicators_same);
452  // AssertThrowMPI(ierr);
453  // if (!(communicators_same == MPI_IDENT ||
454  // communicators_same == MPI_CONGRUENT))
455  // same_communicators = false;
456  // else
457  // same_communicators = true;
458  // }
459 
460  // distinguish three cases. First case: both vectors have the same
461  // layout (just need to copy the local data, not reset the memory and
462  // the underlying Epetra_Map). The third case means that we have to
463  // rebuild the calling vector.
464  if (same_communicators && v.vector->Map().SameAs(vector->Map()))
465  {
466  *vector = *v.vector;
467  if (v.nonlocal_vector.get() != nullptr)
469  std::make_unique<Epetra_MultiVector>(v.nonlocal_vector->Map(), 1);
470  last_action = Zero;
471  }
472  // Second case: vectors have the same global
473  // size, but different parallel layouts (and
474  // one of them a one-to-one mapping). Then we
475  // can call the import/export functionality.
476  else if (size() == v.size() &&
477  (v.vector->Map().UniqueGIDs() || vector->Map().UniqueGIDs()))
478  {
479  reinit(v, false, true);
480  }
481  // Third case: Vectors do not have the same
482  // size.
483  else
484  {
485  vector = std::make_unique<Epetra_FEVector>(*v.vector);
486  last_action = Zero;
489  }
490 
491  if (v.nonlocal_vector.get() != nullptr)
493  std::make_unique<Epetra_MultiVector>(v.nonlocal_vector->Map(), 1);
494 
495  return *this;
496  }
497 
498 
499 
500  Vector &
501  Vector::operator=(Vector &&v) noexcept
502  {
503  static_cast<Subscriptor &>(*this) = static_cast<Subscriptor &&>(v);
504  swap(v);
505  return *this;
506  }
507 
508 
509 
510  template <typename number>
511  Vector &
512  Vector::operator=(const ::Vector<number> &v)
513  {
514  Assert(size() == v.size(), ExcDimensionMismatch(size(), v.size()));
515 
516  // this is probably not very efficient but works. in particular, we could
517  // do better if we know that number==TrilinosScalar because then we could
518  // elide the copying of elements
519  //
520  // let's hope this isn't a particularly frequent operation
521  std::pair<size_type, size_type> local_range = this->local_range();
522  for (size_type i = local_range.first; i < local_range.second; ++i)
523  (*vector)[0][i - local_range.first] = v(i);
524 
525  return *this;
526  }
527 
528 
529 
530  void
532  const Vector & v)
533  {
534  Assert(m.trilinos_matrix().Filled() == true,
535  ExcMessage("Matrix is not compressed. "
536  "Cannot find exchange information!"));
537  Assert(v.vector->Map().UniqueGIDs() == true,
538  ExcMessage("The input vector has overlapping data, "
539  "which is not allowed."));
540 
541  if (vector->Map().SameAs(m.trilinos_matrix().ColMap()) == false)
542  vector =
543  std::make_unique<Epetra_FEVector>(m.trilinos_matrix().ColMap());
544 
545  Epetra_Import data_exchange(vector->Map(), v.vector->Map());
546  const int ierr = vector->Import(*v.vector, data_exchange, Insert);
547 
548  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
549 
550  last_action = Insert;
551  }
552 
553 
554  void
556  const VectorOperation::values operation)
557  {
558  Assert(
559  this->size() == rwv.size(),
560  ExcMessage(
561  "Both vectors need to have the same size for import_elements() to work!"));
562  // TODO: a generic import_elements() function should handle any kind of
563  // data layout in ReadWriteVector, but this function is of limited use as
564  // this class will (hopefully) be retired eventually.
567 
568  if (operation == VectorOperation::insert)
569  {
570  for (const auto idx : this->locally_owned_elements())
571  (*this)[idx] = rwv[idx];
572  }
573  else if (operation == VectorOperation::add)
574  {
575  for (const auto idx : this->locally_owned_elements())
576  (*this)[idx] += rwv[idx];
577  }
578  else
579  AssertThrow(false, ExcNotImplemented());
580 
581  this->compress(operation);
582  }
583 
584 
585  void
587  {
588  // Select which mode to send to Trilinos. Note that we use last_action if
589  // available and ignore what the user tells us to detect wrongly mixed
590  // operations. Typically given_last_action is only used on machines that
591  // do not execute an operation (because they have no own cells for
592  // example).
593  Epetra_CombineMode mode = last_action;
594  if (last_action == Zero)
595  {
596  if (given_last_action == VectorOperation::add)
597  mode = Add;
598  else if (given_last_action == VectorOperation::insert)
599  mode = Insert;
600  else
601  Assert(
602  false,
603  ExcMessage(
604  "compress() can only be called with VectorOperation add, insert, or unknown"));
605  }
606  else
607  {
608  Assert(
609  ((last_action == Add) &&
610  (given_last_action == VectorOperation::add)) ||
611  ((last_action == Insert) &&
612  (given_last_action == VectorOperation::insert)),
613  ExcMessage(
614  "The last operation on the Vector and the given last action in the compress() call do not agree!"));
615  }
616 
617 
618 # ifdef DEBUG
619  // check that every process has decided to use the same mode. This will
620  // otherwise result in undefined behavior in the call to
621  // GlobalAssemble().
622  const double double_mode = mode;
623  const Epetra_MpiComm *comm_ptr =
624  dynamic_cast<const Epetra_MpiComm *>(&(trilinos_partitioner().Comm()));
625  Assert(comm_ptr != nullptr, ExcInternalError());
626 
627  const Utilities::MPI::MinMaxAvg result =
628  Utilities::MPI::min_max_avg(double_mode, comm_ptr->GetMpiComm());
629  Assert(result.max == result.min,
630  ExcMessage(
631  "Not all processors agree whether the last operation on "
632  "this vector was an addition or a set operation. This will "
633  "prevent the compress() operation from succeeding."));
634 
635 # endif
636 
637  // Now pass over the information about what we did last to the vector.
638  int ierr = 0;
639  if (nonlocal_vector.get() == nullptr || mode != Add)
640  ierr = vector->GlobalAssemble(mode);
641  else
642  {
643  Epetra_Export exporter(nonlocal_vector->Map(), vector->Map());
644  ierr = vector->Export(*nonlocal_vector, exporter, mode);
645  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
646  ierr = nonlocal_vector->PutScalar(0.);
647  }
648  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
649  last_action = Zero;
650 
651  compressed = true;
652  }
653 
654 
655 
657  Vector::operator()(const size_type index) const
658  {
659  // Extract local indices in the vector.
660  TrilinosWrappers::types::int_type trilinos_i = vector->Map().LID(
661  static_cast<TrilinosWrappers::types::int_type>(index));
662  TrilinosScalar value = 0.;
663 
664  // If the element is not present on the current processor, we can't
665  // continue. This is the main difference to the el() function.
666  if (trilinos_i == -1)
667  {
668  Assert(false,
670  local_size(),
671  vector->Map().MinMyGID(),
672  vector->Map().MaxMyGID()));
673  }
674  else
675  value = (*vector)[0][trilinos_i];
676 
677  return value;
678  }
679 
680 
681 
682  void
683  Vector::add(const Vector &v, const bool allow_different_maps)
684  {
685  if (allow_different_maps == false)
686  *this += v;
687  else
688  {
690  AssertThrow(size() == v.size(),
691  ExcDimensionMismatch(size(), v.size()));
692 
693 # if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
694  Epetra_Import data_exchange(vector->Map(), v.vector->Map());
695  int ierr =
696  vector->Import(*v.vector, data_exchange, Epetra_AddLocalAlso);
697  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
698  last_action = Add;
699 # else
700  // In versions older than 11.11 the Import function is broken for
701  // adding Hence, we provide a workaround in this case
702 
703  Epetra_MultiVector dummy(vector->Map(), 1, false);
704  Epetra_Import data_exchange(dummy.Map(), v.vector->Map());
705 
706  int ierr = dummy.Import(*v.vector, data_exchange, Insert);
707  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
708 
709  ierr = vector->Update(1.0, dummy, 1.0);
710  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
711 # endif
712  }
713  }
714 
715 
716 
717  bool
718  Vector::operator==(const Vector &v) const
719  {
720  Assert(size() == v.size(), ExcDimensionMismatch(size(), v.size()));
721  if (local_size() != v.local_size())
722  return false;
723 
724  size_type i;
725  for (i = 0; i < local_size(); ++i)
726  if ((*(v.vector))[0][i] != (*vector)[0][i])
727  return false;
728 
729  return true;
730  }
731 
732 
733 
734  bool
735  Vector::operator!=(const Vector &v) const
736  {
737  Assert(size() == v.size(), ExcDimensionMismatch(size(), v.size()));
738 
739  return (!(*this == v));
740  }
741 
742 
743 
744  bool
746  {
747  // get a representation of the vector and
748  // loop over all the elements
749  TrilinosScalar * start_ptr = (*vector)[0];
750  const TrilinosScalar *ptr = start_ptr, *eptr = start_ptr + local_size();
751  unsigned int flag = 0;
752  while (ptr != eptr)
753  {
754  if (*ptr != 0)
755  {
756  flag = 1;
757  break;
758  }
759  ++ptr;
760  }
761 
762  // in parallel, check that the vector
763  // is zero on _all_ processors.
764  const Epetra_MpiComm *mpi_comm =
765  dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
766  Assert(mpi_comm != nullptr, ExcInternalError());
767  unsigned int num_nonzero = Utilities::MPI::sum(flag, mpi_comm->Comm());
768  return num_nonzero == 0;
769  }
770 
771 
772 
773  bool
775  {
776  // get a representation of the vector and
777  // loop over all the elements
778  TrilinosScalar * start_ptr = (*vector)[0];
779  const TrilinosScalar *ptr = start_ptr, *eptr = start_ptr + local_size();
780  unsigned int flag = 0;
781  while (ptr != eptr)
782  {
783  if (*ptr < 0.0)
784  {
785  flag = 1;
786  break;
787  }
788  ++ptr;
789  }
790 
791  // in parallel, check that the vector
792  // is zero on _all_ processors.
793  const auto max_n_negative =
795  return max_n_negative == 0;
796  }
797 
798 
799 
800  void
801  Vector::print(std::ostream & out,
802  const unsigned int precision,
803  const bool scientific,
804  const bool across) const
805  {
806  AssertThrow(out.fail() == false, ExcIO());
807  boost::io::ios_flags_saver restore_flags(out);
808 
809 
810  out.precision(precision);
811  if (scientific)
812  out.setf(std::ios::scientific, std::ios::floatfield);
813  else
814  out.setf(std::ios::fixed, std::ios::floatfield);
815 
816  if (size() != local_size())
817  {
818  auto global_id = [&](const size_type index) {
819  return gid(vector->Map(), index);
820  };
821  out << "size:" << size() << " local_size:" << local_size() << " :"
822  << std::endl;
823  for (size_type i = 0; i < local_size(); ++i)
824  out << "[" << global_id(i) << "]: " << (*(vector))[0][i]
825  << std::endl;
826  }
827  else
828  {
829  TrilinosScalar *val;
830  int leading_dimension;
831  int ierr = vector->ExtractView(&val, &leading_dimension);
832 
833  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
834  if (across)
835  for (size_type i = 0; i < size(); ++i)
836  out << static_cast<double>(val[i]) << ' ';
837  else
838  for (size_type i = 0; i < size(); ++i)
839  out << static_cast<double>(val[i]) << std::endl;
840  out << std::endl;
841  }
842 
843  AssertThrow(out.fail() == false, ExcIO());
844  }
845 
846 
847 
848  void
850  {
854  std::swap(vector, v.vector);
857  }
858 
859 
860 
861  std::size_t
863  {
864  // TODO[TH]: No accurate memory
865  // consumption for Trilinos vectors
866  // yet. This is a rough approximation with
867  // one index and the value per local
868  // entry.
869  return sizeof(*this) +
870  this->local_size() *
871  (sizeof(double) + sizeof(TrilinosWrappers::types::int_type));
872  }
873 
874  // explicit instantiations
875 # ifndef DOXYGEN
876 # include "trilinos_vector.inst"
877 # endif
878  } // namespace MPI
879 } // namespace TrilinosWrappers
880 
882 
883 #endif // DEAL_II_WITH_TRILINOS
unsigned int n_blocks() const
std::size_t size() const
BlockType & block(const unsigned int i)
bool is_ascending_and_one_to_one(const MPI_Comm communicator) const
Definition: index_set.cc:878
size_type size() const
Definition: index_set.h:1661
size_type n_elements() const
Definition: index_set.h:1816
void set_size(const size_type size)
Definition: index_set.h:1649
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:789
void clear()
Definition: index_set.h:1637
void subtract_set(const IndexSet &other)
Definition: index_set.cc:268
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1727
const IndexSet & get_stored_elements() const
void compress(VectorOperation::values operation)
const Epetra_BlockMap & trilinos_partitioner() const
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
void import_elements(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
MPI_Comm get_mpi_communicator() const
reference operator()(const size_type index)
void import_nonlocal_data_for_fe(const ::TrilinosWrappers::SparseMatrix &matrix, const Vector &vector)
std::unique_ptr< Epetra_FEVector > vector
size_type local_size() const
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
IndexSet locally_owned_elements() const
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
bool operator!=(const Vector &v) const
bool operator==(const Vector &v) const
std::pair< size_type, size_type > local_range() const
Vector & operator=(const TrilinosScalar s)
::types::global_dof_index size_type
std::size_t memory_consumption() const
const Epetra_CrsMatrix & trilinos_matrix() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define Assert(cond, exc)
Definition: exceptions.h:1614
static ::ExceptionBase & ExcNotImplemented()
#define AssertIndexRange(index, range)
Definition: exceptions.h:1855
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcGhostsPresent()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1703
void swap(MemorySpaceData< T, MemorySpace > &u, MemorySpaceData< T, MemorySpace > &v)
void swap(BlockVector &u, BlockVector &v)
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
int gid(const Epetra_BlockMap &map, int i)
TrilinosWrappers::types::int_type * my_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition: mpi.cc:150
T max(const T &t, const MPI_Comm mpi_communicator)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm mpi_communicator)
Definition: mpi.cc:124
const Epetra_Comm & comm_self()
double TrilinosScalar
Definition: types.h:175