Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_MultiVecAdapter_decl.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Amesos2: Templated Direct Sparse Solver Package
4 //
5 // Copyright 2011 NTESS and the Amesos2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
21 #ifndef AMESOS2_MULTIVEC_ADAPTER_DECL_HPP
22 #define AMESOS2_MULTIVEC_ADAPTER_DECL_HPP
23 
24 #include <Teuchos_RCP.hpp>
25 #include <Teuchos_Ptr.hpp>
26 #include <Teuchos_ArrayView.hpp>
27 #include <Tpetra_Map.hpp>
28 
29 #include "Amesos2_TypeDecl.hpp"
30 #include "Amesos2_VectorTraits.hpp"
31 
32 namespace Amesos2 {
33 
34 
141  template <class MV>
142  struct MultiVecAdapter {};
143 
144 
152  template <class MV>
153  Teuchos::RCP<MultiVecAdapter<MV> >
154  createMultiVecAdapter(Teuchos::RCP<MV> mv){
155  using Teuchos::rcp;
156 
157  if(mv.is_null()) return Teuchos::null;
158  return( rcp(new MultiVecAdapter<MV>(mv)) );
159  }
160 
161  template <class MV>
162  Teuchos::RCP<const MultiVecAdapter<MV> >
163  createConstMultiVecAdapter(Teuchos::RCP<const MV> mv){
164  using Teuchos::rcp;
165  using Teuchos::rcp_const_cast;
166 
167  if(mv.is_null()) return Teuchos::null;
168  return( rcp(new MultiVecAdapter<MV>(Teuchos::rcp_const_cast<MV,const MV>(mv))).getConst() );
169  }
170 
171 
173  // Utilities for getting and putting data from MultiVecs //
175 
176  namespace Util {
177 
184  template <typename MV, typename V>
186 
187  typedef typename VectorTraits<V>::ptr_scalar_type ptr_return_type ;
188 
189  static ptr_return_type * get_pointer_to_vector ( const Teuchos::Ptr< MV> &mv ) ;
190 
191  static ptr_return_type * get_pointer_to_vector ( Teuchos::Ptr< MV> &mv ) ;
192 
193  static ptr_return_type * get_pointer_to_vector ( const Teuchos::Ptr< const MV > &mv ) ;
194 
195  static ptr_return_type * get_pointer_to_vector ( Teuchos::Ptr< const MV > &mv ) ;
196  };
197 
198  /*
199  * If the multivector scalar type and the desired scalar tpye are
200  * the same, then we can do a simple straight copy.
201  */
202  template <typename MV>
203  struct same_type_get_copy {
204  static void apply(const Teuchos::Ptr<const MV> mv,
205  const Teuchos::ArrayView<typename MV::scalar_t>& v,
206  const size_t ldx,
207  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
208  EDistribution distribution );
209  };
210 
211  /*
212  * In the case where the scalar type of the multi-vector and the
213  * corresponding S type are different, then we need to first get a
214  * copy of the scalar values, then convert each one into the S
215  * type before inserting into the vals array.
216  */
217  template <typename MV, typename S>
218  struct diff_type_get_copy {
219  static void apply(const Teuchos::Ptr<const MV> mv,
220  const Teuchos::ArrayView<S>& v,
221  const size_t ldx,
222  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
223  EDistribution distribution );
224  };
225 
232  template <class MV, typename S>
234  static void
235  do_get (const Teuchos::Ptr<const MV>& mv,
236  const Teuchos::ArrayView<S>& vals,
237  const size_t ldx,
238  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
239  EDistribution distribution = ROOTED);
240 
241  static void
242  do_get (const Teuchos::Ptr<const MV>& mv,
243  const Teuchos::ArrayView<S>& vals,
244  const size_t ldx,
245  EDistribution distribution,
246  typename MV::global_ordinal_t indexBase = 0);
247 
248  static void
249  do_get (const Teuchos::Ptr<const MV>& mv,
250  const Teuchos::ArrayView<S>& vals,
251  const size_t ldx);
252  };
253 
254  /*
255  do_get
256 
257  Return type (bool):
258  true: The input kokkos_vals view is now pointing directly to the adapter's data (same memory and type).
259  If this is x for an Ax=b solve, you don't need 'do_put x' after the solve since you modified the adapter directly.
260  false: The input kokkos_vals view is now resized to match the adapter's size.
261  kokkos_vals will only have the adapter values deep_copied if bInitialize is true (see below).
262  If this is x for an Ax=b solve, you must call 'do_put x' after the solve to deep copy back to the adapter.
263 
264  Inputs
265  bInitialize (bool): tells the adapter whether kokkos_vals needs to have the values of the adapter.
266  true: We require kokkos_vals to have the same size and values of the adapter.
267  For b in Ax=b solves, set bInitialize true because you need the size and values of the adapter.
268  false: We require kokkos_vals to have the same size as the adapter but we don't need the values.
269  For x in Ax=b solves, set bInitialize false because you just need the size, not the values.
270 
271  Note: When this method returns true, meaning direct assignment of the view occurred,
272  bInitialize is not used because you already have the values whether you need them or not.
273 
274  kokkos_vals (View<scalar_t**>): The view which will contain the x or b data.
275  Do not allocate the size of kokkos_vals, let the do_get method do it for you.
276  This is because kokkos_vals may be set to point directly to the adapter memory
277  and then any pre-allocation of size will have been wasted.
278  */
279  template <class MV, typename KV>
280  struct get_1d_copy_helper_kokkos_view {
281  static bool
282  do_get (bool bInitialize,
283  const Teuchos::Ptr<const MV>& mv,
284  KV& kokkos_vals,
285  const size_t ldx,
286  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
287  EDistribution distribution = ROOTED);
288 
289  static bool
290  do_get (bool bInitialize,
291  const Teuchos::Ptr<const MV>& mv,
292  KV& kokkos_vals,
293  const size_t ldx,
294  EDistribution distribution,
295  typename MV::global_ordinal_t indexBase = 0);
296 
297  static bool
298  do_get (bool bInitialize,
299  const Teuchos::Ptr<const MV>& mv,
300  KV& kokkos_vals,
301  const size_t ldx);
302  };
303 
304  /*
305  * If the multivector scalar type and the desired scalar tpye are
306  * the same, then we can do a simple straight copy.
307  */
308  template <typename MV>
309  struct same_type_data_put {
310  static void apply(const Teuchos::Ptr<MV>& mv,
311  const Teuchos::ArrayView<typename MV::scalar_t>& data,
312  const size_t ldx,
313  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
314  EDistribution distribution );
315  };
316 
317  /*
318  * In the case where the scalar type of the multi-vector and the
319  * corresponding S type are different, then we need to first get a
320  * copy of the scalar values, then convert each one into the S
321  * type before inserting into the vals array.
322  */
323  template <typename MV, typename S>
324  struct diff_type_data_put {
325  static void apply(const Teuchos::Ptr<MV>& mv,
326  const Teuchos::ArrayView<S>& data,
327  const size_t ldx,
328  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
329  EDistribution distribution );
330  };
331 
338  template <class MV, typename S>
340  static void do_put(const Teuchos::Ptr<MV>& mv,
341  const Teuchos::ArrayView<S>& data,
342  const size_t ldx,
343  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
344  EDistribution distribution = ROOTED);
345 
346  static void do_put(const Teuchos::Ptr<MV>& mv,
347  const Teuchos::ArrayView<S>& data,
348  const size_t ldx,
349  EDistribution distribution, typename MV::global_ordinal_t indexBase = 0);
350 
351  static void do_put(const Teuchos::Ptr<MV>& mv,
352  const Teuchos::ArrayView<S>& data,
353  const size_t ldx);
354  };
355 
356  template <class MV, typename KV>
357  struct put_1d_data_helper_kokkos_view {
358  static void do_put(const Teuchos::Ptr<MV>& mv,
359  KV& kokkos_data,
360  const size_t ldx,
361  Teuchos::Ptr<const Tpetra::Map<typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t> > distribution_map,
362  EDistribution distribution = ROOTED);
363 
364  static void do_put(const Teuchos::Ptr<MV>& mv,
365  KV& kokkos_data,
366  const size_t ldx,
367  EDistribution distribution, typename MV::global_ordinal_t indexBase = 0);
368 
369  static void do_put(const Teuchos::Ptr<MV>& mv,
370  KV& kokkos_data,
371  const size_t ldx);
372  };
373  }
374 } // end namespace Amesos2
375 
378 #ifdef HAVE_AMESOS2_EPETRA
380 #endif
381 
382 #endif // AMESOS2_MULTIVEC_ADAPTER_DECL_HPP
Teuchos::RCP< MultiVecAdapter< MV > > createMultiVecAdapter(Teuchos::RCP< MV > mv)
Factory creation method for MultiVecAdapters.
Definition: Amesos2_MultiVecAdapter_decl.hpp:154
Amesos2::MultiVecAdapter specialization for the Kokkos::View class.
Helper class for getting 1-D copies of multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:233
Helper struct for getting pointers to the MV data - only used when number of vectors = 1 and single M...
Definition: Amesos2_MultiVecAdapter_decl.hpp:185
static void do_get(const Teuchos::Ptr< const MV > &mv, const Teuchos::ArrayView< S > &vals, const size_t ldx, Teuchos::Ptr< const Tpetra::Map< typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t > > distribution_map, EDistribution distribution=ROOTED)
Helper class for getting 1-D copies of multivectors.
Definition: Amesos2_MultiVecAdapter_def.hpp:112
Amesos2::MultiVecAdapter specialization for the Tpetra::MultiVector class.
Amesos2::MultiVecAdapter specialization for the Epetra_MultiVector class.
Enum and other types declarations for Amesos2.
static void do_put(const Teuchos::Ptr< MV > &mv, const Teuchos::ArrayView< S > &data, const size_t ldx, Teuchos::Ptr< const Tpetra::Map< typename MV::local_ordinal_t, typename MV::global_ordinal_t, typename MV::node_t > > distribution_map, EDistribution distribution=ROOTED)
Helper class for putting 1-D data arrays into multivectors.
Definition: Amesos2_MultiVecAdapter_def.hpp:287
EDistribution
Definition: Amesos2_TypeDecl.hpp:89
Helper class for putting 1-D data arrays into multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:339
A templated MultiVector class adapter for Amesos2.
Definition: Amesos2_MultiVecAdapter_decl.hpp:142