42 #ifndef AMESOS2_SOLVER_UQ_PCE_HPP
43 #define AMESOS2_SOLVER_UQ_PCE_HPP
45 #include "Amesos2_Solver.hpp"
46 #include "Amesos2_Factory.hpp"
52 #if defined(TPETRA_HAVE_KOKKOS_REFACTOR)
53 template <
class S,
class LO,
class GO,
class D>
79 class Device,
template<
class,
class>
class ConcreteSolver>
80 class PCESolverAdapter :
81 public Solver< Tpetra::CrsMatrix<Sacado::UQ::PCE<Storage>,
84 Kokkos::Compat::KokkosDeviceWrapperNode<Device> >,
85 Tpetra::MultiVector<Sacado::UQ::PCE<Storage>,
88 Kokkos::Compat::KokkosDeviceWrapperNode<Device> >
94 typedef Kokkos::Compat::KokkosDeviceWrapperNode<Device>
Node;
95 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Matrix;
96 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Vector;
99 typedef Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Map;
100 typedef Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> FlatGraph;
101 typedef Tpetra::CrsMatrix<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> FlatMatrix;
102 typedef Tpetra::MultiVector<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> FlatVector;
103 typedef ConcreteSolver<FlatMatrix,FlatVector> FlatConcreteSolver;
104 typedef Solver<FlatMatrix,FlatVector> FlatSolver;
106 typedef Solver<Matrix,Vector> solver_type;
107 typedef typename solver_type::type type;
108 typedef typename Scalar::cijk_type cijk_type;
115 A(A_), X(X_),
B(B_) {
116 cijk = get_pce_cijk(
A, X,
B);
118 Stokhos::create_flat_pce_graph(*(
A->getCrsGraph()),
131 create_solver_with_supported_type<ConcreteSolver,FlatMatrix,FlatVector>::apply(flat_A, flat_X, flat_B);
144 virtual type& preOrdering(
void ) {
145 flat_solver->preOrdering();
155 virtual type& symbolicFactorization(
void ) {
156 flat_solver->symbolicFactorization();
173 virtual type& numericFactorization(
void ) {
174 flat_solver->numericFactorization();
191 virtual void solve(
void ) {
192 flat_solver->solve();
232 virtual void solve(Vector* XX,
const Vector* BB)
const {
257 virtual type& setParameters(
259 flat_solver->setParameters(parameterList);
269 getValidParameters(
void )
const {
270 return flat_solver->getValidParameters();
300 EPhase keep_phase = CLEAN ) {
304 cijk = get_pce_cijk(
A);
305 if (keep_phase <= CLEAN) {
310 Stokhos::create_flat_pce_graph(*(
A->getCrsGraph()),
317 if (keep_phase <= SYMBFACT)
320 flat_solver->setA(flat_A, keep_phase);
342 virtual void setA(
const Matrix* a, EPhase keep_phase = CLEAN ) {
348 virtual bool matrixShapeOK(
void ) {
349 return flat_solver->matrixShapeOK();
360 flat_solver->setX(flat_X);
365 virtual void setX( Vector* x ) {
374 flat_solver->setX(flat_X);
385 virtual Vector* getXRaw(
void ) {
397 flat_solver->setB(flat_B);
402 virtual void setB(
const Vector* b ) {
411 flat_solver->setB(flat_B);
422 virtual const Vector* getBRaw(
void ) {
429 return flat_solver->getComm();
434 virtual Status& getStatus()
const {
435 return flat_solver->getStatus();
440 virtual std::string name(
void )
const {
441 return flat_solver->name();
451 virtual std::string description(
void )
const {
453 return flat_solver->description();
461 flat_solver->describe(out, verbLevel);
474 flat_solver->printTiming(out, verbLevel);
487 flat_solver->getTiming(timingParameterList);
510 template <
template <
class,
class>
class ConcreteSolver,
511 class ST,
class LO,
class GO,
class D >
512 struct create_solver_with_supported_type<
514 Tpetra::CrsMatrix<Sacado::UQ::PCE<ST>,LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<D> >,
515 Tpetra::MultiVector<Sacado::UQ::PCE<ST>,LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<D> > > {
517 typedef Kokkos::Compat::KokkosDeviceWrapperNode<D>
NO;
518 typedef Tpetra::CrsMatrix<SC,LO,GO,NO> Matrix;
519 typedef Tpetra::MultiVector<SC,LO,GO,NO> Vector;
527 typename MatrixTraits<Matrix>::scalar_t,
528 typename MultiVecAdapter<Vector>::scalar_t
530 > same_scalar_assertion;
531 (void)same_scalar_assertion;
534 return Teuchos::rcp(
new PCESolverAdapter<ST,LO,GO,D,ConcreteSolver>(A, X, B) );
541 template <
template <
class,
class>
class ConcreteSolver,
543 struct solver_supports_scalar<ConcreteSolver, Sacado::UQ::PCE<Storage> > {
546 typedef typename solver_traits<ConcreteSolver>::supported_scalars supported_scalars;
547 static const bool value =
548 Meta::if_then_else<Meta::is_same<supported_scalars, Meta::nil_t>::value,
550 Meta::type_list_contains<supported_scalars,
551 BaseScalar> >::type::value;
558 #endif // AMESOS2_SOLVER_UQ_PCE_HPP
Stokhos::StandardStorage< int, double > Storage
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< View< T, P...> >::value, unsigned >::type dimension_scalar(const View< T, P...> &view)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const Tpetra::MultiVector< typename Storage::value_type, LocalOrdinal, GlobalOrdinal, Node > > create_flat_vector_view(const Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LocalOrdinal, GlobalOrdinal, Node > &vec_const, const Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &flat_map)
KokkosClassic::DefaultNode::DefaultNodeType Node
Teuchos::RCP< Tpetra::CrsMatrix< typename Storage::value_type, LocalOrdinal, GlobalOrdinal, Node > > create_flat_matrix(const Tpetra::CrsMatrix< Sacado::MP::Vector< Storage >, LocalOrdinal, GlobalOrdinal, Node > &mat, const Teuchos::RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &flat_graph, const LocalOrdinal block_size)
static const EVerbosityLevel verbLevel_default
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< view_type >::value, typename CijkType< view_type >::type >::type cijk(const view_type &view)