Reference documentation for deal.II version 9.3.3
\(\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_epetra_vector.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2015 - 2021 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# ifdef DEAL_II_WITH_MPI
21
23
25
26# include <boost/io/ios_state.hpp>
27
28# include <Epetra_Import.h>
29# include <Epetra_Map.h>
30# include <Epetra_MpiComm.h>
31
32# include <memory>
33
34
36
37namespace LinearAlgebra
38{
39 namespace EpetraWrappers
40 {
42 : Subscriptor()
43 , vector(new Epetra_FEVector(
44 Epetra_Map(0, 0, 0, Utilities::Trilinos::comm_self())))
45 {}
46
47
48
50 : Subscriptor()
51 , vector(new Epetra_FEVector(V.trilinos_vector()))
52 {}
53
54
55
56 Vector::Vector(const IndexSet &parallel_partitioner,
57 const MPI_Comm &communicator)
58 : Subscriptor()
59 , vector(new Epetra_FEVector(
60 parallel_partitioner.make_trilinos_map(communicator, false)))
61 {}
62
63
64
65 void
66 Vector::reinit(const IndexSet &parallel_partitioner,
67 const MPI_Comm &communicator,
68 const bool omit_zeroing_entries)
69 {
70 Epetra_Map input_map =
71 parallel_partitioner.make_trilinos_map(communicator, false);
72 if (vector->Map().SameAs(input_map) == false)
73 vector = std::make_unique<Epetra_FEVector>(input_map);
74 else if (omit_zeroing_entries == false)
75 {
76 const int ierr = vector->PutScalar(0.);
77 Assert(ierr == 0, ExcTrilinosError(ierr));
78 (void)ierr;
79 }
80 }
81
82
83
84 void
86 const bool omit_zeroing_entries)
87 {
88 // Check that casting will work.
89 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
91
92 // Downcast V. If fails, throws an exception.
93 const Vector &down_V = dynamic_cast<const Vector &>(V);
94
96 down_V.get_mpi_communicator(),
97 omit_zeroing_entries);
98 }
99
100
101
102 Vector &
104 {
105 // Distinguish three cases:
106 // - First case: both vectors have the same layout.
107 // - Second case: both vectors have the same size but different layout.
108 // - Third case: the vectors have different size.
109 if (vector->Map().SameAs(V.trilinos_vector().Map()))
110 *vector = V.trilinos_vector();
111 else
112 {
113 if (size() == V.size())
114 {
115 Epetra_Import data_exchange(vector->Map(),
116 V.trilinos_vector().Map());
117
118 const int ierr =
119 vector->Import(V.trilinos_vector(), data_exchange, Insert);
120 Assert(ierr == 0, ExcTrilinosError(ierr));
121 (void)ierr;
122 }
123 else
124 vector = std::make_unique<Epetra_FEVector>(V.trilinos_vector());
125 }
126
127 return *this;
128 }
129
130
131
132 Vector &
133 Vector::operator=(const double s)
134 {
135 Assert(s == 0., ExcMessage("Only 0 can be assigned to a vector."));
136
137 const int ierr = vector->PutScalar(s);
138 Assert(ierr == 0, ExcTrilinosError(ierr));
139 (void)ierr;
140
141 return *this;
142 }
143
144
145
146 void
149 VectorOperation::values operation,
150 std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
151 communication_pattern)
152 {
153 // If no communication pattern is given, create one. Otherwise, use the
154 // one given.
155 if (communication_pattern == nullptr)
156 {
157 // The first time import is called, a communication pattern is
158 // created. Check if the communication pattern already exists and if
159 // it can be reused.
161 V.get_stored_elements().size()) ||
162 (source_stored_elements != V.get_stored_elements()))
163 {
165 V.get_stored_elements(),
166 dynamic_cast<const Epetra_MpiComm &>(vector->Comm()).Comm());
167 }
168 }
169 else
170 {
172 std::dynamic_pointer_cast<const CommunicationPattern>(
173 communication_pattern);
175 epetra_comm_pattern != nullptr,
177 std::string("The communication pattern is not of type ") +
178 "LinearAlgebra::EpetraWrappers::CommunicationPattern."));
179 }
180
181 Epetra_Import import(epetra_comm_pattern->get_epetra_import());
182
183 // The TargetMap and the SourceMap have their roles inverted.
184 Epetra_FEVector source_vector(import.TargetMap());
185 double * values = source_vector.Values();
186 std::copy(V.begin(), V.end(), values);
187
188 if (operation == VectorOperation::insert)
189 vector->Export(source_vector, import, Insert);
190 else if (operation == VectorOperation::add)
191 vector->Export(source_vector, import, Add);
192 else if (operation == VectorOperation::max)
193 vector->Export(source_vector, import, Epetra_Max);
194 else if (operation == VectorOperation::min)
195 vector->Export(source_vector, import, Epetra_Min);
196 else
198 }
199
200
201
202 Vector &
203 Vector::operator*=(const double factor)
204 {
205 AssertIsFinite(factor);
206 vector->Scale(factor);
207
208 return *this;
209 }
210
211
212
213 Vector &
214 Vector::operator/=(const double factor)
215 {
216 AssertIsFinite(factor);
217 Assert(factor != 0., ExcZero());
218 *this *= 1. / factor;
219
220 return *this;
221 }
222
223
224
225 Vector &
227 {
228 // Check that casting will work.
229 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
231
232 // Downcast V. If fails, throws an exception.
233 const Vector &down_V = dynamic_cast<const Vector &>(V);
234 // If the maps are the same we can Update right away.
235 if (vector->Map().SameAs(down_V.trilinos_vector().Map()))
236 {
237 const int ierr = vector->Update(1., down_V.trilinos_vector(), 1.);
238 Assert(ierr == 0, ExcTrilinosError(ierr));
239 (void)ierr;
240 }
241 else
242 {
243 Assert(this->size() == down_V.size(),
244 ExcDimensionMismatch(this->size(), down_V.size()));
245
246# if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
247 Epetra_Import data_exchange(vector->Map(),
248 down_V.trilinos_vector().Map());
249 const int ierr = vector->Import(down_V.trilinos_vector(),
250 data_exchange,
251 Epetra_AddLocalAlso);
252 Assert(ierr == 0, ExcTrilinosError(ierr));
253 (void)ierr;
254# else
255 // In versions older than 11.11 the Import function is broken for
256 // adding Hence, we provide a workaround in this case
257
258 Epetra_MultiVector dummy(vector->Map(), 1, false);
259 Epetra_Import data_exchange(dummy.Map(),
260 down_V.trilinos_vector().Map());
261
262 int ierr =
263 dummy.Import(down_V.trilinos_vector(), data_exchange, Insert);
264 Assert(ierr == 0, ExcTrilinosError(ierr));
265
266 ierr = vector->Update(1.0, dummy, 1.0);
267 Assert(ierr == 0, ExcTrilinosError(ierr));
268 (void)ierr;
269# endif
270 }
271
272 return *this;
273 }
274
275
276
277 Vector &
279 {
280 this->add(-1., V);
281
282 return *this;
283 }
284
285
286
288 {
289 // Check that casting will work.
290 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
292
293 // Downcast V. If fails, throws an exception.
294 const Vector &down_V = dynamic_cast<const Vector &>(V);
295 Assert(this->size() == down_V.size(),
296 ExcDimensionMismatch(this->size(), down_V.size()));
297 Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
299
300 double result(0.);
301 const int ierr = vector->Dot(down_V.trilinos_vector(), &result);
302 Assert(ierr == 0, ExcTrilinosError(ierr));
303 (void)ierr;
304
305 return result;
306 }
307
308
309
310 void
311 Vector::add(const double a)
312 {
314 const unsigned local_size(vector->MyLength());
315 for (unsigned int i = 0; i < local_size; ++i)
316 (*vector)[0][i] += a;
317 }
318
319
320
321 void
323 {
324 // Check that casting will work.
325 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
327
328 // Downcast V. If fails, throws an exception.
329 const Vector &down_V = dynamic_cast<const Vector &>(V);
331 Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
333
334 const int ierr = vector->Update(a, down_V.trilinos_vector(), 1.);
335 Assert(ierr == 0, ExcTrilinosError(ierr));
336 (void)ierr;
337 }
338
339
340
341 void
342 Vector::add(const double a,
344 const double b,
346 {
347 // Check that casting will work.
348 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
350 // Check that casting will work.
351 Assert(dynamic_cast<const Vector *>(&W) != nullptr,
353
354 // Downcast V. If fails, throws an exception.
355 const Vector &down_V = dynamic_cast<const Vector &>(V);
356 // Downcast W. If fails, throws an exception.
357 const Vector &down_W = dynamic_cast<const Vector &>(W);
358 Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
360 Assert(vector->Map().SameAs(down_W.trilinos_vector().Map()),
364
365 const int ierr = vector->Update(
366 a, down_V.trilinos_vector(), b, down_W.trilinos_vector(), 1.);
367 Assert(ierr == 0, ExcTrilinosError(ierr));
368 (void)ierr;
369 }
370
371
372
373 void
374 Vector::sadd(const double s,
375 const double a,
377 {
378 // Check that casting will work.
379 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
381
382 *this *= s;
383 // Downcast V. It fails, throws an exception.
384 const Vector &down_V = dynamic_cast<const Vector &>(V);
385 Vector tmp(down_V);
386 tmp *= a;
387 *this += tmp;
388 }
389
390
391
392 void
394 {
395 // Check that casting will work.
396 Assert(dynamic_cast<const Vector *>(&scaling_factors) != nullptr,
398
399 // Downcast scaling_factors. If fails, throws an exception.
400 const Vector &down_scaling_factors =
401 dynamic_cast<const Vector &>(scaling_factors);
402 Assert(vector->Map().SameAs(down_scaling_factors.trilinos_vector().Map()),
404
405 const int ierr = vector->Multiply(1.0,
406 down_scaling_factors.trilinos_vector(),
407 *vector,
408 0.0);
409 Assert(ierr == 0, ExcTrilinosError(ierr));
410 (void)ierr;
411 }
412
413
414
415 void
417 {
418 // Check that casting will work.
419 Assert(dynamic_cast<const Vector *>(&V) != nullptr,
421
422 // Downcast V. If fails, throws an exception.
423 const Vector &down_V = dynamic_cast<const Vector &>(V);
424 // If we don't have the same map, copy.
425 if (vector->Map().SameAs(down_V.trilinos_vector().Map()) == false)
426 this->sadd(0., a, V);
427 else
428 {
429 // Otherwise, just update
430 int ierr = vector->Update(a, down_V.trilinos_vector(), 0.);
431 Assert(ierr == 0, ExcTrilinosError(ierr));
432 (void)ierr;
433 }
434 }
435
436
437
438 bool
440 {
441 // get a representation of the vector and
442 // loop over all the elements
443 double * start_ptr = (*vector)[0];
444 const double *ptr = start_ptr, *eptr = start_ptr + vector->MyLength();
445 unsigned int flag = 0;
446 while (ptr != eptr)
447 {
448 if (*ptr != 0)
449 {
450 flag = 1;
451 break;
452 }
453 ++ptr;
454 }
455
456 // Check that the vector is zero on _all_ processors.
457 const Epetra_MpiComm *mpi_comm =
458 dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
459 Assert(mpi_comm != nullptr, ExcInternalError());
460 unsigned int num_nonzero = Utilities::MPI::sum(flag, mpi_comm->Comm());
461
462 return num_nonzero == 0;
463 }
464
465
466
467 double
469 {
470 double mean_value(0.);
471
472 int ierr = vector->MeanValue(&mean_value);
473 Assert(ierr == 0, ExcTrilinosError(ierr));
474 (void)ierr;
475
476 return mean_value;
477 }
478
479
480
481 double
483 {
484 double norm(0.);
485 int ierr = vector->Norm1(&norm);
486 Assert(ierr == 0, ExcTrilinosError(ierr));
487 (void)ierr;
488
489 return norm;
490 }
491
492
493
494 double
496 {
497 double norm(0.);
498 int ierr = vector->Norm2(&norm);
499 Assert(ierr == 0, ExcTrilinosError(ierr));
500 (void)ierr;
501
502 return norm;
503 }
504
505
506
507 double
509 {
510 double norm(0.);
511 int ierr = vector->NormInf(&norm);
512 Assert(ierr == 0, ExcTrilinosError(ierr));
513 (void)ierr;
514
515 return norm;
516 }
517
518
519
520 double
521 Vector::add_and_dot(const double a,
524 {
525 this->add(a, V);
526
527 return *this * W;
528 }
529
530
531
534 {
535# ifndef DEAL_II_WITH_64BIT_INDICES
536 return vector->GlobalLength();
537# else
538 return vector->GlobalLength64();
539# endif
540 }
541
542
543
546 {
547 return vector->MyLength();
548 }
549
550
551
554 {
555 const Epetra_MpiComm *epetra_comm =
556 dynamic_cast<const Epetra_MpiComm *>(&(vector->Comm()));
557 Assert(epetra_comm != nullptr, ExcInternalError());
558 return epetra_comm->GetMpiComm();
559 }
560
561
562
565 {
566 IndexSet is(size());
567
568 // easy case: local range is contiguous
569 if (vector->Map().LinearMap())
570 {
571# ifndef DEAL_II_WITH_64BIT_INDICES
572 is.add_range(vector->Map().MinMyGID(), vector->Map().MaxMyGID() + 1);
573# else
574 is.add_range(vector->Map().MinMyGID64(),
575 vector->Map().MaxMyGID64() + 1);
576# endif
577 }
578 else if (vector->Map().NumMyElements() > 0)
579 {
580 const size_type n_indices = vector->Map().NumMyElements();
581# ifndef DEAL_II_WITH_64BIT_INDICES
582 unsigned int *vector_indices =
583 reinterpret_cast<unsigned int *>(vector->Map().MyGlobalElements());
584# else
585 size_type *vector_indices =
586 reinterpret_cast<size_type *>(vector->Map().MyGlobalElements64());
587# endif
588 is.add_indices(vector_indices, vector_indices + n_indices);
589 }
590 is.compress();
591
592 return is;
593 }
594
595
596
597 const Epetra_FEVector &
599 {
600 return *vector;
601 }
602
603
604
605 Epetra_FEVector &
607 {
608 return *vector;
609 }
610
611
612
613 void
614 Vector::print(std::ostream & out,
615 const unsigned int precision,
616 const bool scientific,
617 const bool across) const
618 {
619 AssertThrow(out, ExcIO());
620 boost::io::ios_flags_saver restore_flags(out);
621
622 // Get a representation of the vector and loop over all
623 // the elements
624 double *val;
625 int leading_dimension;
626 int ierr = vector->ExtractView(&val, &leading_dimension);
627
628 Assert(ierr == 0, ExcTrilinosError(ierr));
629 (void)ierr;
630 out.precision(precision);
631 if (scientific)
632 out.setf(std::ios::scientific, std::ios::floatfield);
633 else
634 out.setf(std::ios::fixed, std::ios::floatfield);
635
636 if (across)
637 for (int i = 0; i < vector->MyLength(); ++i)
638 out << val[i] << ' ';
639 else
640 for (int i = 0; i < vector->MyLength(); ++i)
641 out << val[i] << std::endl;
642 out << std::endl;
643
644 // restore the representation
645 // of the vector
646 AssertThrow(out, ExcIO());
647 }
648
649
650
651 std::size_t
653 {
654 return sizeof(*this) +
655 vector->MyLength() *
656 (sizeof(double) + sizeof(TrilinosWrappers::types::int_type));
657 }
658
659
660
661 void
663 const MPI_Comm &mpi_comm)
664 {
665 source_stored_elements = source_index_set;
667 std::make_shared<CommunicationPattern>(locally_owned_elements(),
668 source_index_set,
669 mpi_comm);
670 }
671 } // namespace EpetraWrappers
672} // namespace LinearAlgebra
673
675
676# endif
677
678#endif
size_type size() const
Definition: index_set.h:1634
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:596
void add_range(const size_type begin, const size_type end)
Definition: index_set.h:1673
void compress() const
Definition: index_set.h:1642
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1703
virtual double mean_value() const override
virtual Vector & operator/=(const double factor) override
virtual double l2_norm() const override
std::unique_ptr< Epetra_FEVector > vector
const Epetra_FEVector & trilinos_vector() const
virtual double l1_norm() const override
virtual ::IndexSet locally_owned_elements() const override
virtual double add_and_dot(const double a, const VectorSpaceVector< double > &V, const VectorSpaceVector< double > &W) override
virtual void scale(const VectorSpaceVector< double > &scaling_factors) override
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm)
virtual double operator*(const VectorSpaceVector< double > &V) const override
void reinit(const IndexSet &parallel_partitioner, const MPI_Comm &communicator, const bool omit_zeroing_entries=false)
virtual size_type size() const override
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override
virtual Vector & operator-=(const VectorSpaceVector< double > &V) override
virtual double linfty_norm() const override
virtual void import(const ReadWriteVector< double > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) override
virtual void add(const double a) override
virtual void equ(const double a, const VectorSpaceVector< double > &V) override
virtual std::size_t memory_consumption() const override
std::shared_ptr< const CommunicationPattern > epetra_comm_pattern
virtual Vector & operator*=(const double factor) override
virtual Vector & operator+=(const VectorSpaceVector< double > &V) override
virtual void sadd(const double s, const double a, const VectorSpaceVector< double > &V) override
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcVectorTypeNotCompatible()
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertIsFinite(number)
Definition: exceptions.h:1721
static ::ExceptionBase & ExcDifferentParallelPartitioning()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
static const char V
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition: divergence.h:472
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
T sum(const T &t, const MPI_Comm &mpi_communicator)
const Epetra_Comm & comm_self()
Definition: utilities.cc:1113
void copy(const T *begin, const T *end, U *dest)