Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_MatrixFreeOperator.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_MATRIX_FREE_OPERATOR_HPP
11 #define STOKHOS_MATRIX_FREE_OPERATOR_HPP
12 
13 #include "Stokhos_SGOperator.hpp"
14 #include "EpetraExt_MultiComm.h"
17 #include "Epetra_Map.h"
18 #include "Epetra_Import.h"
20 #include "Teuchos_Array.hpp"
21 
22 namespace Stokhos {
23 
29 
30  public:
31 
42 
44  virtual ~MatrixFreeOperator();
45 
47  double countApplyFlops() const;
48 
51 
53  virtual void setupOperator(
55 
59 
62  getSGPolynomial() const;
63 
65 
68 
70  virtual int SetUseTranspose(bool UseTranspose);
71 
76  virtual int Apply(const Epetra_MultiVector& Input,
77  Epetra_MultiVector& Result) const;
78 
83  virtual int ApplyInverse(const Epetra_MultiVector& X,
84  Epetra_MultiVector& Y) const;
85 
87  virtual double NormInf() const;
88 
90  virtual const char* Label () const;
91 
93  virtual bool UseTranspose() const;
94 
99  virtual bool HasNormInf() const;
100 
105  virtual const Epetra_Comm & Comm() const;
106 
111  virtual const Epetra_Map& OperatorDomainMap () const;
112 
117  virtual const Epetra_Map& OperatorRangeMap () const;
118 
120 
121  private:
122 
125 
128 
129  protected:
130 
132  std::string label;
133 
136 
139 
142 
145 
148 
151 
154 
157 
160 
163 
166 
169 
172 
175 
178 
181 
183  bool scale_op;
184 
187 
190 
193 
196 
199 
202 
205 
208 
211 
214 
217 
220 
223 
226 
229 
230  }; // class MatrixFreeOperator
231 
232 } // namespace Stokhos
233 
234 #endif // STOKHOS_MATRIX_FREE_OPERATOR_HPP
MatrixFreeOperator & operator=(const MatrixFreeOperator &)
Private to prohibit copying.
MatrixFreeOperator(const Teuchos::RCP< const EpetraExt::MultiComm > &sg_comm, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > &epetraCijk, const Teuchos::RCP< const Epetra_Map > &domain_base_map, const Teuchos::RCP< const Epetra_Map > &range_base_map, const Teuchos::RCP< const Epetra_Map > &domain_sg_map, const Teuchos::RCP< const Epetra_Map > &range_sg_map, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor.
An Epetra operator representing the block stochastic Galerkin operator.
Teuchos::RCP< const Cijk_type > Cijk
Stores triple product tensor.
Teuchos::RCP< Epetra_MultiVector > tmp
Temporary multivector used in Apply()
virtual double NormInf() const
Returns an approximate infinity norm of the operator matrix.
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this matrix operator.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of the inverse of the operator applied to a Epetra_MultiVector Input in Result as ...
virtual const char * Label() const
Returns a character string describing the operator.
Cijk_type::k_iterator k_end
Ending k iterator.
int expansion_size
Number of terms in expansion.
Stokhos::Sparse3Tensor< int, double > Cijk_type
Short-hand for Cijk.
Teuchos::RCP< Epetra_Map > global_col_map_trans
Stores operator column SG map for transpose.
Teuchos::RCP< const Epetra_BlockMap > stoch_col_map
Stores stochastic part of column map.
Teuchos::RCP< Epetra_MultiVector > input_col_trans
Temporary to store result of importing input into column map (transpose)
bool useTranspose
Flag indicating whether transpose was selected.
virtual int SetUseTranspose(bool UseTranspose)
Set to true if the transpose of the operator is requested.
Teuchos::RCP< Epetra_MultiVector > input_col
Temporary to store result of importing input into column map.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerking basis.
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Bi-directional iterator for traversing a sparse array.
Teuchos::RCP< Epetra_MultiVector > tmp_trans
Temporary multivector used in Apply() for transpose.
bool is_stoch_parallel
Whether we have parallelism over stochastic blocks.
Teuchos::RCP< Epetra_Map > global_col_map
Stores operator column SG map.
An abstract class to represent a generic stochastic Galerkin operator as an Epetra_Operator.
int max_num_mat_vec
Maximum number of matvecs in Apply.
Cijk_type::k_iterator k_begin
Starting k iterator.
Teuchos::RCP< const Epetra_Map > range_base_map
Stores range base map.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Stores Epetra Cijk tensor.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
bool include_mean
Flag indicating whether to include mean term.
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > block_ops
Stores operators.
virtual const Epetra_Comm & Comm() const
Returns a reference to the Epetra_Comm communicator associated with this operator.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Stores SG parallel communicator.
virtual int Apply(const Epetra_MultiVector &Input, Epetra_MultiVector &Result) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector Input in Result as described ...
Teuchos::Array< Teuchos::RCP< Epetra_MultiVector > > result_block
MultiVectors for each block for Apply() result.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()
Get SG polynomial.
Teuchos::RCP< Epetra_Import > col_importer
Importer from domain map to column map.
Teuchos::RCP< const Epetra_Map > domain_base_map
Stores domain base map.
bool scale_op
Flag indicating whether operator be scaled with &lt;^2&gt;
Teuchos::RCP< const Epetra_Map > domain_sg_map
Stores domain SG map.
int num_blocks
Number of Jacobian blocks (not necessarily equal to expansion_size)
bool use_block_apply
Flag indicating whether to use block Epetra_MultiVector apply.
bool only_use_linear
Flag indicating whether to only use linear terms.
double countApplyFlops() const
Return number of FLOPS for each call to Apply()
Teuchos::RCP< const Epetra_Map > range_sg_map
Stores range SG map.
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this matrix operator. ...
Teuchos::RCP< Epetra_Import > col_importer_trans
Importer from range map to column map.
virtual void setupOperator(const Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > &poly)
Setup operator.
std::string label
Label for operator.
Teuchos::Array< Teuchos::RCP< const Epetra_MultiVector > > input_block
MultiVectors for each block for Apply() input.