44 #include <Teuchos_ArrayRCP.hpp> 45 #include <Teuchos_CommHelpers.hpp> 46 #include <Teuchos_OpaqueWrapper.hpp> 49 #include <Teuchos_DefaultMpiComm.hpp> 53 #include <Tpetra_Export.hpp> 62 const Teuchos::RCP<
const Tpetra::MultiVector<double, int, SupportId>>
64 const Teuchos::RCP<
const Tpetra::Map<int, SupportId>> &domain_map,
65 const Teuchos::RCP<
const Tpetra::Map<int, SupportId>> &range_map )
66 : d_comm( polynomial->getMap()->getComm() )
67 , d_polynomial( polynomial )
68 , d_domain_map( domain_map )
69 , d_range_map( range_map )
76 const Tpetra::MultiVector<double, int, SupportId> &X,
77 Tpetra::MultiVector<double, int, SupportId> &Y, Teuchos::ETransp mode,
78 double alpha,
double beta )
const 80 DTK_REQUIRE( d_domain_map->isSameAs( *( X.getMap() ) ) );
81 DTK_REQUIRE( d_range_map->isSameAs( *( Y.getMap() ) ) );
82 DTK_REQUIRE( X.getNumVectors() == Y.getNumVectors() );
85 int poly_size = d_polynomial->getNumVectors();
86 int num_vec = X.getNumVectors();
87 int local_length = d_polynomial->getLocalLength();
88 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const double>> poly_view =
89 d_polynomial->get2dView();
90 Teuchos::ArrayRCP<Teuchos::ArrayRCP<double>> y_view = Y.get2dViewNonConst();
96 if ( Teuchos::NO_TRANS == mode )
99 Teuchos::Array<double> x_poly( poly_size * num_vec, 0.0 );
100 if ( 0 == d_comm()->getRank() )
102 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const double>> x_view =
104 for (
int n = 0; n < num_vec; ++n )
106 std::copy( &x_view[n][0], &x_view[n][0] + poly_size,
107 &x_poly[n * poly_size] );
110 Teuchos::broadcast( *d_comm, 0, x_poly() );
114 for (
int n = 0; n < num_vec; ++n )
116 stride = n * poly_size;
118 for (
int i = 0; i < local_length; ++i )
120 for (
int p = 0; p < poly_size; ++p )
123 alpha * poly_view[p][i] * x_poly[stride + p];
130 else if ( Teuchos::TRANS == mode )
133 Tpetra::MultiVector<double, int, SupportId> work( Y.getMap(),
137 Tpetra::Export<int, SupportId> exporter( X.getMap(), work.getMap() );
138 work.doExport( X, exporter, Tpetra::INSERT );
141 Teuchos::ArrayRCP<Teuchos::ArrayRCP<double>> work_view =
142 work.get2dViewNonConst();
143 Teuchos::Array<double> products( poly_size * num_vec, 0.0 );
145 for (
int n = 0; n < num_vec; ++n )
147 stride = n * poly_size;
148 for (
int p = 0; p < poly_size; ++p )
150 for (
int i = 0; i < local_length; ++i )
152 products[stride + p] += poly_view[p][i] * work_view[n][i];
158 Teuchos::Array<double> product_sums( poly_size * num_vec, 0.0 );
160 Teuchos::RCP<const Teuchos::MpiComm<int>> mpi_comm =
161 Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm<int>>( d_comm );
162 Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm>> opaque_comm =
163 mpi_comm->getRawMpiComm();
164 MPI_Comm raw_comm = ( *opaque_comm )();
165 MPI_Reduce( products.getRawPtr(), product_sums.getRawPtr(),
166 poly_size * num_vec, MPI_DOUBLE, MPI_SUM, 0, raw_comm );
168 product_sums = products;
172 if ( 0 == d_comm->getRank() )
174 for (
int n = 0; n < num_vec; ++n )
176 for (
int p = 0; p < poly_size; ++p )
178 y_view[n][p] += alpha * product_sums[n * poly_size + p];
186 DTK_INSIST( mode == Teuchos::NO_TRANS || mode == Teuchos::TRANS );
void apply(const Tpetra::MultiVector< double, int, SupportId > &X, Tpetra::MultiVector< double, int, SupportId > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, double alpha=Teuchos::ScalarTraits< double >::one(), double beta=Teuchos::ScalarTraits< double >::zero()) const override
Computes the operator-multivector application.
PolynomialMatrix(const Teuchos::RCP< const Tpetra::MultiVector< double, int, SupportId >> &polynomial, const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &domain_map, const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &range_map)
Constructor.
Assertions and Design-by-Contract for error handling.