Reference documentation for deal.II version 9.1.0-pre
copy.cc
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 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 
16 #include <deal.II/sundials/copy.h>
17 
18 #ifdef DEAL_II_WITH_SUNDIALS
19 
20 DEAL_II_NAMESPACE_OPEN
21 namespace SUNDIALS
22 {
23  namespace internal
24  {
32  inline std::size_t
33  N_Vector_length(const N_Vector &vec)
34  {
35  const N_Vector_ID id = N_VGetVectorID(vec);
36  long int length = -1;
37  switch (id)
38  {
39  case SUNDIALS_NVEC_SERIAL:
40  {
41  length = NV_LENGTH_S(vec);
42  break;
43  }
44 # ifdef DEAL_II_WITH_MPI
45  case SUNDIALS_NVEC_PARALLEL:
46  {
47  length = NV_LOCLENGTH_P(vec);
48  break;
49  }
50 # endif
51  default:
52  Assert(false, ExcNotImplemented());
53  }
54 
55  Assert(length >= 0, ExcInternalError());
56  return static_cast<std::size_t>(length);
57  }
58 
59 # ifdef DEAL_II_WITH_MPI
60 
61 # ifdef DEAL_II_WITH_TRILINOS
62 
63 
64  void
65  copy(TrilinosWrappers::MPI::Vector &dst, const N_Vector &src)
66  {
67  const IndexSet is = dst.locally_owned_elements();
68  const size_t N = is.n_elements();
69  AssertDimension(N, N_Vector_length(src));
70  for (size_t i = 0; i < N; ++i)
71  {
72  dst[is.nth_index_in_set(i)] = NV_Ith_P(src, i);
73  }
75  }
76 
77  void
78  copy(N_Vector &dst, const TrilinosWrappers::MPI::Vector &src)
79  {
80  const IndexSet is = src.locally_owned_elements();
81  const size_t N = is.n_elements();
82  AssertDimension(N, N_Vector_length(dst));
83  for (size_t i = 0; i < N; ++i)
84  {
85  NV_Ith_P(dst, i) = src[is.nth_index_in_set(i)];
86  }
87  }
88 
89  void
90  copy(TrilinosWrappers::MPI::BlockVector &dst, const N_Vector &src)
91  {
92  const IndexSet is = dst.locally_owned_elements();
93  const size_t N = is.n_elements();
94  AssertDimension(N, N_Vector_length(src));
95  for (size_t i = 0; i < N; ++i)
96  {
97  dst[is.nth_index_in_set(i)] = NV_Ith_P(src, i);
98  }
100  }
101 
102  void
103  copy(N_Vector &dst, const TrilinosWrappers::MPI::BlockVector &src)
104  {
105  IndexSet is = src.locally_owned_elements();
106  const size_t N = is.n_elements();
107  AssertDimension(N, N_Vector_length(dst));
108  for (size_t i = 0; i < N; ++i)
109  {
110  NV_Ith_P(dst, i) = src[is.nth_index_in_set(i)];
111  }
112  }
113 
114 # endif // DEAL_II_WITH_TRILINOS
115 
116 # ifdef DEAL_II_WITH_PETSC
117 # ifndef PETSC_USE_COMPLEX
118 
119  void
120  copy(PETScWrappers::MPI::Vector &dst, const N_Vector &src)
121  {
122  const IndexSet is = dst.locally_owned_elements();
123  const size_t N = is.n_elements();
124  AssertDimension(N, N_Vector_length(src));
125  for (size_t i = 0; i < N; ++i)
126  {
127  dst[is.nth_index_in_set(i)] = NV_Ith_P(src, i);
128  }
130  }
131 
132  void
133  copy(N_Vector &dst, const PETScWrappers::MPI::Vector &src)
134  {
135  const IndexSet is = src.locally_owned_elements();
136  const size_t N = is.n_elements();
137  AssertDimension(N, N_Vector_length(dst));
138  for (size_t i = 0; i < N; ++i)
139  {
140  NV_Ith_P(dst, i) = src[is.nth_index_in_set(i)];
141  }
142  }
143 
144  void
145  copy(PETScWrappers::MPI::BlockVector &dst, const N_Vector &src)
146  {
147  const IndexSet is = dst.locally_owned_elements();
148  const size_t N = is.n_elements();
149  AssertDimension(N, N_Vector_length(src));
150  for (size_t i = 0; i < N; ++i)
151  {
152  dst[is.nth_index_in_set(i)] = NV_Ith_P(src, i);
153  }
155  }
156 
157  void
158  copy(N_Vector &dst, const PETScWrappers::MPI::BlockVector &src)
159  {
160  const IndexSet is = src.locally_owned_elements();
161  const size_t N = is.n_elements();
162  AssertDimension(N, N_Vector_length(dst));
163  for (size_t i = 0; i < N; ++i)
164  {
165  NV_Ith_P(dst, i) = src[is.nth_index_in_set(i)];
166  }
167  }
168 
169 # endif // PETSC_USE_COMPLEX
170 # endif // DEAL_II_WITH_PETSC
171 
172 # endif // mpi
173 
174  void
175  copy(BlockVector<double> &dst, const N_Vector &src)
176  {
177  const size_t N = dst.size();
178  AssertDimension(N_Vector_length(src), N);
179  for (size_t i = 0; i < N; ++i)
180  {
181  dst[i] = NV_Ith_S(src, i);
182  }
183  }
184 
185  void
186  copy(N_Vector &dst, const BlockVector<double> &src)
187  {
188  const size_t N = src.size();
189  AssertDimension(N_Vector_length(dst), N);
190  for (size_t i = 0; i < N; ++i)
191  {
192  NV_Ith_S(dst, i) = src[i];
193  }
194  }
195 
196  void
197  copy(Vector<double> &dst, const N_Vector &src)
198  {
199  const size_t N = dst.size();
200  AssertDimension(N_Vector_length(src), N);
201  for (size_t i = 0; i < N; ++i)
202  {
203  dst[i] = NV_Ith_S(src, i);
204  }
205  }
206 
207  void
208  copy(N_Vector &dst, const Vector<double> &src)
209  {
210  const size_t N = src.size();
211  AssertDimension(N_Vector_length(dst), N);
212  for (size_t i = 0; i < N; ++i)
213  {
214  NV_Ith_S(dst, i) = src[i];
215  }
216  }
217  } // namespace internal
218 } // namespace SUNDIALS
219 DEAL_II_NAMESPACE_CLOSE
220 
221 #endif
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
void compress(::VectorOperation::values operation)
size_type size() const
size_type n_elements() const
Definition: index_set.h:1743
void compress(const VectorOperation::values operation)
#define Assert(cond, exc)
Definition: exceptions.h:1227
void compress(::VectorOperation::values operation)
IndexSet locally_owned_elements() const
IndexSet locally_owned_elements() const
static::ExceptionBase & ExcNotImplemented()
size_type nth_index_in_set(const unsigned int local_index) const
Definition: index_set.h:1793
IndexSet locally_owned_elements() const
static::ExceptionBase & ExcInternalError()