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_epetra_vector.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 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
20# include <deal.II/base/mpi.h>
22
24
25# include <boost/io/ios_state.hpp>
26
27# include <Epetra_Import.h>
28# include <Epetra_Map.h>
29# include <Epetra_MpiComm.h>
30
31# include <memory>
32
33
35
36namespace LinearAlgebra
37{
38 namespace EpetraWrappers
39 {
40# ifndef DOXYGEN
41 namespace internal
42 {
43 VectorReference::operator value_type() const
44 {
45 AssertIndexRange(index, vector.size());
46
47 // Trilinos allows for vectors to be referenced by the [] or ()
48 // operators but only () checks index bounds. We check these bounds by
49 // ourselves, so we can use []. Note that we can only get local values.
50
51 const TrilinosWrappers::types::int_type local_index =
52 vector.vector->Map().LID(
53 static_cast<TrilinosWrappers::types::int_type>(index));
54
55# ifndef DEAL_II_WITH_64BIT_INDICES
56 Assert(local_index >= 0,
57 ExcAccessToNonLocalElement(index,
58 vector.vector->Map().NumMyElements(),
59 vector.vector->Map().MinMyGID(),
60 vector.vector->Map().MaxMyGID()));
61# else
62 Assert(local_index >= 0,
63 ExcAccessToNonLocalElement(index,
64 vector.vector->Map().NumMyElements(),
65 vector.vector->Map().MinMyGID64(),
66 vector.vector->Map().MaxMyGID64()));
67# endif
68
69 return (*(vector.vector))[0][local_index];
70 }
71 } // namespace internal
72# endif
73
74
75 // Check that the class we declare here satisfies the
76 // vector-space-vector concept. If we catch it here,
77 // any mistake in the vector class declaration would
78 // show up in uses of this class later on as well.
79# ifdef DEAL_II_HAVE_CXX20
80 static_assert(concepts::is_vector_space_vector<Vector>);
81# endif
82
84 : vector(new Epetra_FEVector(
85 Epetra_Map(0, 0, 0, Utilities::Trilinos::comm_self())))
86 {}
87
88
89
91 : vector(new Epetra_FEVector(V.trilinos_vector()))
92 {}
93
94
95
96 Vector::Vector(const IndexSet &parallel_partitioner,
97 const MPI_Comm communicator)
98 : vector(new Epetra_FEVector(
99 parallel_partitioner.make_trilinos_map(communicator, false)))
100 {}
101
102
103
104 void
105 Vector::reinit(const IndexSet &parallel_partitioner,
106 const MPI_Comm communicator,
107 const bool omit_zeroing_entries)
108 {
109 Epetra_Map input_map =
110 parallel_partitioner.make_trilinos_map(communicator, false);
111 if (vector->Map().SameAs(input_map) == false)
112 vector = std::make_unique<Epetra_FEVector>(input_map);
113 else if (omit_zeroing_entries == false)
114 {
115 const int ierr = vector->PutScalar(0.);
116 Assert(ierr == 0, ExcTrilinosError(ierr));
117 }
118 }
119
120
121
122 void
123 Vector::reinit(const Vector &V, const bool omit_zeroing_entries)
124 {
125 reinit(V.locally_owned_elements(),
126 V.get_mpi_communicator(),
127 omit_zeroing_entries);
128 }
129
130
131
132 void
135 const ArrayView<double> &elements) const
136 {
137 AssertDimension(indices.size(), elements.size());
138 const auto &vector = trilinos_vector();
139 const auto &map = vector.Map();
140
141 for (unsigned int i = 0; i < indices.size(); ++i)
142 {
143 AssertIndexRange(indices[i], size());
144 const auto trilinos_i =
145 map.LID(static_cast<TrilinosWrappers::types::int_type>(indices[i]));
146 elements[i] = vector[0][trilinos_i];
147 }
148 }
149
150
151
152 Vector &
154 {
155 // Distinguish three cases:
156 // - First case: both vectors have the same layout.
157 // - Second case: both vectors have the same size but different layout.
158 // - Third case: the vectors have different size.
159 if (vector->Map().SameAs(V.trilinos_vector().Map()))
160 *vector = V.trilinos_vector();
161 else
162 {
163 if (size() == V.size())
164 {
165 Epetra_Import data_exchange(vector->Map(),
166 V.trilinos_vector().Map());
167
168 const int ierr =
169 vector->Import(V.trilinos_vector(), data_exchange, Insert);
170 Assert(ierr == 0, ExcTrilinosError(ierr));
171 }
172 else
173 vector = std::make_unique<Epetra_FEVector>(V.trilinos_vector());
174 }
175
176 return *this;
177 }
178
179
180
181 Vector &
182 Vector::operator=(const double s)
183 {
184 Assert(s == 0., ExcMessage("Only 0 can be assigned to a vector."));
185
186 const int ierr = vector->PutScalar(s);
187 Assert(ierr == 0, ExcTrilinosError(ierr));
188
189 return *this;
190 }
191
192
193
194 void
197 VectorOperation::values operation,
198 const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
199 &communication_pattern)
200 {
201 // If no communication pattern is given, create one. Otherwise, use the
202 // one given.
203 if (communication_pattern == nullptr)
204 {
205 // The first time import is called, a communication pattern is
206 // created. Check if the communication pattern already exists and if
207 // it can be reused.
208 if ((source_stored_elements.size() !=
209 V.get_stored_elements().size()) ||
210 (source_stored_elements != V.get_stored_elements()))
211 {
213 V.get_stored_elements(),
214 dynamic_cast<const Epetra_MpiComm &>(vector->Comm()).Comm());
215 }
216 }
217 else
218 {
220 std::dynamic_pointer_cast<const CommunicationPattern>(
221 communication_pattern);
223 epetra_comm_pattern != nullptr,
224 ExcMessage("The communication pattern is not of type "
225 "LinearAlgebra::EpetraWrappers::CommunicationPattern."));
226 }
227
228 Epetra_Import import_map(epetra_comm_pattern->get_epetra_import());
229
230 // The TargetMap and the SourceMap have their roles inverted.
231 Epetra_FEVector source_vector(import_map.TargetMap());
232 double *values = source_vector.Values();
233 std::copy(V.begin(), V.end(), values);
234
235 if (operation == VectorOperation::insert)
236 vector->Export(source_vector, import_map, Insert);
237 else if (operation == VectorOperation::add)
238 vector->Export(source_vector, import_map, Add);
239 else if (operation == VectorOperation::max)
240 vector->Export(source_vector, import_map, Epetra_Max);
241 else if (operation == VectorOperation::min)
242 vector->Export(source_vector, import_map, Epetra_Min);
243 else
245 }
246
247
248
249 Vector &
250 Vector::operator*=(const double factor)
251 {
252 AssertIsFinite(factor);
253 vector->Scale(factor);
254
255 return *this;
256 }
257
258
259
260 Vector &
261 Vector::operator/=(const double factor)
262 {
263 AssertIsFinite(factor);
264 Assert(factor != 0., ExcZero());
265 *this *= 1. / factor;
266
267 return *this;
268 }
269
270
271
272 Vector &
274 {
275 // If the maps are the same we can Update right away.
276 if (vector->Map().SameAs(V.trilinos_vector().Map()))
277 {
278 const int ierr = vector->Update(1., V.trilinos_vector(), 1.);
279 Assert(ierr == 0, ExcTrilinosError(ierr));
280 }
281 else
282 {
283 Assert(this->size() == V.size(),
284 ExcDimensionMismatch(this->size(), V.size()));
285
286 Epetra_Import data_exchange(vector->Map(), V.trilinos_vector().Map());
287 const int ierr = vector->Import(V.trilinos_vector(),
288 data_exchange,
289 Epetra_AddLocalAlso);
290 Assert(ierr == 0, ExcTrilinosError(ierr));
291 }
292
293 return *this;
294 }
295
296
297
298 Vector &
300 {
301 this->add(-1., V);
302
303 return *this;
304 }
305
306
307
308 double
309 Vector::operator*(const Vector &V) const
310 {
311 Assert(this->size() == V.size(),
312 ExcDimensionMismatch(this->size(), V.size()));
313 Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
315
316 double result(0.);
317 const int ierr = vector->Dot(V.trilinos_vector(), &result);
318 Assert(ierr == 0, ExcTrilinosError(ierr));
319
320 return result;
321 }
322
323
324
325 void
326 Vector::add(const double a)
327 {
329 const unsigned local_size(vector->MyLength());
330 for (unsigned int i = 0; i < local_size; ++i)
331 (*vector)[0][i] += a;
332 }
333
334
335
336 void
337 Vector::add(const double a, const Vector &V)
338 {
340 Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
342
343 const int ierr = vector->Update(a, V.trilinos_vector(), 1.);
344 Assert(ierr == 0, ExcTrilinosError(ierr));
345 }
346
347
348
349 void
350 Vector::add(const double a,
351 const Vector &V,
352 const double b,
353 const Vector &W)
354 {
355 Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
357 Assert(vector->Map().SameAs(W.trilinos_vector().Map()),
361
362 const int ierr =
363 vector->Update(a, V.trilinos_vector(), b, W.trilinos_vector(), 1.);
364 Assert(ierr == 0, ExcTrilinosError(ierr));
365 }
366
367
368
369 void
370 Vector::sadd(const double s, const double a, const Vector &V)
371 {
372 *this *= s;
373 Vector tmp(V);
374 tmp *= a;
375 *this += tmp;
376 }
377
378
379
380 void
381 Vector::scale(const Vector &scaling_factors)
382 {
383 Assert(vector->Map().SameAs(scaling_factors.trilinos_vector().Map()),
385
386 const int ierr =
387 vector->Multiply(1.0, scaling_factors.trilinos_vector(), *vector, 0.0);
388 Assert(ierr == 0, ExcTrilinosError(ierr));
389 }
390
391
392
393 void
394 Vector::equ(const double a, const Vector &V)
395 {
396 // If we don't have the same map, copy.
397 if (vector->Map().SameAs(V.trilinos_vector().Map()) == false)
398 this->sadd(0., a, V);
399 else
400 {
401 // Otherwise, just update
402 int ierr = vector->Update(a, V.trilinos_vector(), 0.);
403 Assert(ierr == 0, ExcTrilinosError(ierr));
404 }
405 }
406
407
408
409 bool
411 {
412 // get a representation of the vector and
413 // loop over all the elements
414 double *start_ptr = (*vector)[0];
415 const double *ptr = start_ptr, *eptr = start_ptr + vector->MyLength();
416 unsigned int flag = 0;
417 while (ptr != eptr)
418 {
419 if (*ptr != 0)
420 {
421 flag = 1;
422 break;
423 }
424 ++ptr;
425 }
426
427 // Check that the vector is zero on _all_ processors.
428 const Epetra_MpiComm *mpi_comm =
429 dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
430 Assert(mpi_comm != nullptr, ExcInternalError());
431 unsigned int num_nonzero = Utilities::MPI::sum(flag, mpi_comm->Comm());
432
433 return num_nonzero == 0;
434 }
435
436
437
438 double
440 {
441 double mean_value(0.);
442
443 int ierr = vector->MeanValue(&mean_value);
444 Assert(ierr == 0, ExcTrilinosError(ierr));
445
446 return mean_value;
447 }
448
449
450
451 double
453 {
454 double norm(0.);
455 int ierr = vector->Norm1(&norm);
456 Assert(ierr == 0, ExcTrilinosError(ierr));
457
458 return norm;
459 }
460
461
462
463 double
465 {
466 double norm(0.);
467 int ierr = vector->Norm2(&norm);
468 Assert(ierr == 0, ExcTrilinosError(ierr));
469
470 return norm;
471 }
472
473
474
475 double
477 {
478 double norm(0.);
479 int ierr = vector->NormInf(&norm);
480 Assert(ierr == 0, ExcTrilinosError(ierr));
481
482 return norm;
483 }
484
485
486
487 double
488 Vector::add_and_dot(const double a, const Vector &V, const Vector &W)
489 {
490 this->add(a, V);
491
492 return *this * W;
493 }
494
495
496
499 {
500# ifndef DEAL_II_WITH_64BIT_INDICES
501 return vector->GlobalLength();
502# else
503 return vector->GlobalLength64();
504# endif
505 }
506
507
508
511 {
512 return vector->MyLength();
513 }
514
515
516
519 {
520 const Epetra_MpiComm *epetra_comm =
521 dynamic_cast<const Epetra_MpiComm *>(&(vector->Comm()));
522 Assert(epetra_comm != nullptr, ExcInternalError());
523 return epetra_comm->GetMpiComm();
524 }
525
526
527
530 {
531 IndexSet is(size());
532
533 // easy case: local range is contiguous
534 if (vector->Map().LinearMap())
535 {
536# ifndef DEAL_II_WITH_64BIT_INDICES
537 is.add_range(vector->Map().MinMyGID(), vector->Map().MaxMyGID() + 1);
538# else
539 is.add_range(vector->Map().MinMyGID64(),
540 vector->Map().MaxMyGID64() + 1);
541# endif
542 }
543 else if (vector->Map().NumMyElements() > 0)
544 {
545 const size_type n_indices = vector->Map().NumMyElements();
546# ifndef DEAL_II_WITH_64BIT_INDICES
547 unsigned int *vector_indices =
548 reinterpret_cast<unsigned int *>(vector->Map().MyGlobalElements());
549# else
550 size_type *vector_indices =
551 reinterpret_cast<size_type *>(vector->Map().MyGlobalElements64());
552# endif
553 is.add_indices(vector_indices, vector_indices + n_indices);
554 }
555 is.compress();
556
557 return is;
558 }
559
560
561 void
563 {}
564
565
566
567 const Epetra_FEVector &
569 {
570 return *vector;
571 }
572
573
574
575 Epetra_FEVector &
577 {
578 return *vector;
579 }
580
581
582
583 void
584 Vector::print(std::ostream &out,
585 const unsigned int precision,
586 const bool scientific,
587 const bool across) const
588 {
589 AssertThrow(out.fail() == false, ExcIO());
590 boost::io::ios_flags_saver restore_flags(out);
591
592 // Get a representation of the vector and loop over all
593 // the elements
594 double *val;
595 int leading_dimension;
596 int ierr = vector->ExtractView(&val, &leading_dimension);
597
598 Assert(ierr == 0, ExcTrilinosError(ierr));
599 out.precision(precision);
600 if (scientific)
601 out.setf(std::ios::scientific, std::ios::floatfield);
602 else
603 out.setf(std::ios::fixed, std::ios::floatfield);
604
605 if (across)
606 for (int i = 0; i < vector->MyLength(); ++i)
607 out << val[i] << ' ';
608 else
609 for (int i = 0; i < vector->MyLength(); ++i)
610 out << val[i] << std::endl;
611 out << std::endl;
612
613 // restore the representation
614 // of the vector
615 AssertThrow(out.fail() == false, ExcIO());
616 }
617
618
619
620 std::size_t
622 {
623 return sizeof(*this) +
624 vector->MyLength() *
625 (sizeof(double) + sizeof(TrilinosWrappers::types::int_type));
626 }
627
628
629
630 void
632 const MPI_Comm mpi_comm)
633 {
634 source_stored_elements = source_index_set;
636 std::make_shared<CommunicationPattern>(locally_owned_elements(),
637 source_index_set,
638 mpi_comm);
639 }
640 } // namespace EpetraWrappers
641} // namespace LinearAlgebra
642
644
645#endif
std::size_t size() const
Definition array_view.h:689
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
void add_range(const size_type begin, const size_type end)
Definition index_set.h:1818
void compress() const
Definition index_set.h:1799
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition index_set.h:1846
void equ(const double a, const Vector &V)
void compress(const VectorOperation::values operation)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
std::unique_ptr< Epetra_FEVector > vector
const Epetra_FEVector & trilinos_vector() const
void import_elements(const ReadWriteVector< double > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
void sadd(const double s, const double a, const Vector &V)
virtual size_type size() const override
void scale(const Vector &scaling_factors)
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm mpi_comm)
std::shared_ptr< const CommunicationPattern > epetra_comm_pattern
void reinit(const IndexSet &parallel_partitioner, const MPI_Comm communicator, const bool omit_zeroing_entries=false)
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, const ArrayView< double > &elements) const override
double add_and_dot(const double a, const Vector &V, const Vector &W)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIsFinite(number)
static ::ExceptionBase & ExcDifferentParallelPartitioning()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
T sum(const T &t, const MPI_Comm mpi_communicator)