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