10#ifndef TPETRA_APPLY_HELPERS_HPP 
   11#define TPETRA_APPLY_HELPERS_HPP 
   13#include "Teuchos_ScalarTraits.hpp" 
   14#include "Tpetra_CrsMatrix.hpp" 
   43template <
class MatrixArray, 
class MultiVectorArray>
 
   45                  const typename std::remove_pointer<typename MultiVectorArray::value_type>::type &
X,
 
   47                  typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type 
alpha = Teuchos::ScalarTraits<
typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::one(),
 
   48                  typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type 
beta  = Teuchos::ScalarTraits<
typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::zero(),
 
   49                  Teuchos::RCP<Teuchos::ParameterList> 
params                                             = Teuchos::null) {
 
   52  using size_type   = 
typename MatrixArray::size_type;
 
   53  using matrix_type = 
typename std::remove_pointer<typename MatrixArray::value_type>::type;
 
   54  using map_type    = 
typename matrix_type::map_type;
 
   55  using import_type = 
typename matrix_type::import_type;
 
   56  using export_type = 
typename matrix_type::export_type;
 
   57  using MV          = 
typename matrix_type::MV;
 
   58  using scalar_type = 
typename matrix_type::scalar_type;
 
   61  using Teuchos::rcp_const_cast;
 
   68  int numRanks = 
X.getMap()->getComm()->getSize();
 
   71  for (size_type 
i = 0; 
i < 
N - 1; 
i++) {
 
   80  if (
params.is_null() || !
params->isParameter(
"can batch")) {
 
   90      if (!
Y[
i]->isDistributed())
 
  108        for (size_type 
i = 0; 
i < 
N; 
i++) 
Y[
i]->putScalar(
ZERO);
 
  110        for (size_type 
i = 0; 
i < 
N; 
i++) 
Y[
i]->scale(
beta);
 
  129    for (size_type 
i = 0; 
i < 
N; 
i++) {
 
  147        if (!
Y[
i]->isConstantStride() || 
X_colMap.getRawPtr() == 
Y[
i]) {
 
  163    for (size_type 
i = 0; 
i < 
N; 
i++) {
 
 
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
 
Struct that holds views of the contents of a CrsMatrix.
 
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
 
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
 
A parallel distribution of indices over processes.
 
One or more distributed dense vectors.
 
Namespace Tpetra contains the class and methods constituting the Tpetra library.
 
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
 
void batchedApply(const MatrixArray &Matrices, const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &X, MultiVectorArray &Y, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type alpha=Teuchos::ScalarTraits< typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type >::one(), typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type beta=Teuchos::ScalarTraits< typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type >::zero(), Teuchos::RCP< Teuchos::ParameterList > params=Teuchos::null)
Does multiply matrix apply() calls with a single X vector.
 
@ ADD_ASSIGN
Accumulate new values into existing values (may not be supported in all classes)
 
@ INSERT
Insert new values that don't currently exist.
 
@ ZERO
Replace old values with zero.