Flow123d  JS_before_hm-896-g486f41f
vector_mpi.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file vector_mpi.cc
15  * @brief
16  */
17 
18 #include "la/vector_mpi.hh"
19 #include <vector>
20 #include <memory>
21 #include "system/system.hh"
22 #include "system/global_defs.h"
23 #include "system/index_types.hh"
24 
25 #include <petscvec.h>
26 #include <armadillo>
27 
28 
30 : communicator_(comm)
31 {}
32 
33 
34 VectorMPI::VectorMPI(unsigned int local_size, MPI_Comm comm)
35 : communicator_(comm) {
36  resize(local_size);
37 }
38 
39 VectorMPI::VectorMPI(unsigned int local_size, std::vector<LongIdx> &ghost_idx)
40 : communicator_(PETSC_COMM_WORLD) {
41  resize(local_size, ghost_idx);
42 }
43 
45 {
46  return VectorMPI(size, PETSC_COMM_SELF);
47 }
48 
49 
50 void VectorMPI::resize(unsigned int local_size) {
51  if (data_ptr_.use_count() ==0) {
52  data_ptr_ = std::make_shared< std::vector<double> >(local_size);
53  } else {
54  ASSERT_DBG( data_ptr_.use_count() == 1 ) ( data_ptr_.use_count() ).error("Object referenced by other pointer. Can not resize.");
55  chkerr(VecDestroy(&data_petsc_));
56  data_ptr_->resize(local_size);
57  }
58  if (communicator_ == PETSC_COMM_SELF)
59  chkerr(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, local_size, &((*data_ptr_)[0]), &data_petsc_));
60  else
61  chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, local_size, PETSC_DECIDE, &((*data_ptr_)[0]), &data_petsc_));
62  chkerr(VecZeroEntries( data_petsc_ ));
63 }
64 
65 
66 void VectorMPI::resize(unsigned int local_size, std::vector<LongIdx> &ghost_idx) {
67  ASSERT_DBG(communicator_ == PETSC_COMM_WORLD).error("Cannot allocate ghost values in sequential vector.");
68  if (data_ptr_.use_count() ==0) {
69  data_ptr_ = std::make_shared< std::vector<double> >(local_size + ghost_idx.size());
70  } else {
71  ASSERT_DBG( data_ptr_.use_count() == 1 ) ( data_ptr_.use_count() ).error("Object referenced by other pointer. Can not resize.");
72  chkerr(VecDestroy(&data_petsc_));
73  data_ptr_->resize(local_size + ghost_idx.size());
74  }
75  chkerr(VecCreateGhostWithArray(PETSC_COMM_WORLD, local_size, PETSC_DECIDE, ghost_idx.size(), ghost_idx.data(), data_ptr_->data(), &data_petsc_));
76  chkerr(VecZeroEntries( data_petsc_ ));
77 }
78 
79 
82  this->resize(other.data().size());
83 }
84 
85 
86 void VectorMPI::swap(VectorMPI &other) {
88  ASSERT_EQ(this->data_ptr_->size(), other.data_ptr_->size());
89  uint size = this->data_ptr_->size();
90  std::swap(this->data_ptr_, other.data_ptr_);
91  chkerr(VecDestroy(&data_petsc_));
92  chkerr(VecDestroy(&other.data_petsc_));
93  if (communicator_ == PETSC_COMM_SELF) {
94  chkerr(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size, &((*data_ptr_)[0]), &data_petsc_));
95  chkerr(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size, &((*other.data_ptr_)[0]), &other.data_petsc_));
96  } else {
97  chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, size, PETSC_DECIDE, &((*data_ptr_)[0]), &data_petsc_));
98  chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, size, PETSC_DECIDE, &((*other.data_ptr_)[0]), &other.data_petsc_));
99  }
100 }
101 
102 
104  ASSERT_EQ(this->communicator_, other.communicator_);
105  ASSERT_EQ(this->data_ptr_->size(), other.data_ptr_->size());
106  chkerr(VecCopy(other.data_petsc_, data_petsc_));
107 }
108 
109 
111 {
113  arma::vec vec(loc_indices.n_elem);
114 
115  for(unsigned int i=0; i<loc_indices.n_elem; i++){
116  unsigned int idx = loc_indices(i);
117  ASSERT_DBG(idx < data_ptr_->size()) (idx) (data_ptr_->size());
118  vec(i) = (*data_ptr_)[idx];
119  }
120  return vec;
121 }
122 
123 
124 arma::vec VectorMPI::get_subvec(const LocDofVec& loc_indices) const
125 {
127  arma::vec vec(loc_indices.n_elem);
128 
129  for(unsigned int i=0; i<loc_indices.n_elem; i++){
130  unsigned int idx = loc_indices(i);
131  ASSERT_DBG(idx < data_ptr_->size()) (idx) (data_ptr_->size());
132  vec(i) = (*data_ptr_)[idx];
133  }
134  return vec;
135 }
136 
137 
138 void VectorMPI::set_subvec(const LocDofVec& loc_indices, const arma::vec& values)
139 {
141  ASSERT_EQ_DBG(loc_indices.n_elem, values.n_elem);
142 
143  for(unsigned int i=0; i<loc_indices.n_elem; i++){
144  unsigned int idx = loc_indices(i);
145  ASSERT_DBG(idx < data_ptr_->size()) (idx) (data_ptr_->size());
146  (*data_ptr_)[idx] = values(i);
147  }
148 }
149 
150 /// Destructor.
152 {
153  if (data_ptr_.use_count() == 1)
154  if (data_petsc_) chkerr(VecDestroy(&data_petsc_));
155 }
unsigned int size() const
Return size of output data.
Definition: vector_mpi.hh:98
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
MPI_Comm communicator_
communicator
Definition: vector_mpi.hh:182
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
ArmaVec< double, N > vec
Definition: armor.hh:861
unsigned int uint
int MPI_Comm
Definition: mpi.h:141
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:148
~VectorMPI()
Destructor.
Definition: vector_mpi.cc:151
VectorMPI(MPI_Comm comm=PETSC_COMM_SELF)
Definition: vector_mpi.cc:29
void copy_from(VectorMPI &other)
Definition: vector_mpi.cc:103
Global macros to enhance readability and debugging, general constants.
VectorDataPtr data_ptr_
shared pointer to vector of data
Definition: vector_mpi.hh:178
Vec data_petsc_
stored vector of data in PETSC format
Definition: vector_mpi.hh:180
void resize(unsigned int local_size)
Definition: vector_mpi.cc:50
void swap(nlohmann::json &j1, nlohmann::json &j2) noexcept(is_nothrow_move_constructible< nlohmann::json >::value andis_nothrow_move_assignable< nlohmann::json >::value)
exchanges the values of two JSON objects
Definition: json.hpp:8688
void swap(VectorMPI &other)
Swaps the current vector data with the other vector data.
Definition: vector_mpi.cc:86
#define ASSERT_DBG(expr)
void set_subvec(const LocDofVec &loc_indices, const arma::vec &values)
Definition: vector_mpi.cc:138
void duplicate_from(VectorMPI other)
For the current vector, it creates the same parallel structure as the other vector has...
Definition: vector_mpi.cc:80
arma::vec get_subvec(const LocDofVec &loc_indices)
Definition: vector_mpi.cc:110
VectorData & data()
Definition: vector_mpi.hh:85
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.cc:44