25 label(
"Stokhos Approximate Schur Complement Preconditioner"),
28 epetraCijk(epetraCijk_),
31 prec_factory(prec_factory_),
36 Cijk(epetraCijk->getParallelCijk()),
39 upper_block_Cijk(P+1),
40 lower_block_Cijk(P+1),
43 only_use_linear(
true),
48 epetraCijk->isStochasticParallel(), std::logic_error,
49 "Stokhos::ApproxSchurComplementPreconditioner does not support " <<
50 "a parallel stochastic distribution.");
52 scale_op = params_->
get(
"Scale Operator by Inverse Basis Norms",
true);
53 symmetric = params_->
get(
"Symmetric Gauss-Seidel",
false);
63 int nj =
Cijk->num_j(k);
71 int d = prod_basis->dimension();
73 for (
int p=0; p<=
P; p++) {
100 double c = value(i_it);
109 else if (p_col < p_row) {
116 for (
int p=0; p<=
P; p++) {
134 mean_prec = prec_factory->compute(sg_poly->getCoeffPtr(0));
135 label = std::string(
"Stokhos Approximate Schur Complement Preconditioner:\n")
136 + std::string(
" ***** ") + std::string(mean_prec->Label());
143 useTranspose = UseTranspose;
144 sg_op->SetUseTranspose(UseTranspose);
145 mean_prec->SetUseTranspose(UseTranspose);
154 return sg_op->Apply(Input, Result);
161 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
168 bool made_copy =
false;
169 if (Input.Values() == Result.Values()) {
175 int m = input->NumVectors();
176 if (rhs_block ==
Teuchos::null || rhs_block->NumVectors() != m)
178 Teuchos::rcp(
new EpetraExt::BlockMultiVector(*base_map, *sg_map, m));
179 if (tmp ==
Teuchos::null || tmp->NumVectors() != m*max_num_mat_vec)
182 j_ptr.resize(m*max_num_mat_vec);
183 mj_indices.resize(m*max_num_mat_vec);
186 EpetraExt::BlockMultiVector input_block(
View, *base_map, *input);
187 EpetraExt::BlockMultiVector result_block(
View, *base_map, Result);
189 result_block.PutScalar(0.0);
192 rhs_block->Update(1.0, input_block, 0.0);
198 for (
int l=P; l>=1; l--) {
200 divide_diagonal_block(block_indices[l], block_indices[l+1],
201 *rhs_block, result_block);
204 multiply_block(upper_block_Cijk[l], -1.0, result_block, *rhs_block);
208 divide_diagonal_block(0, 1, *rhs_block, result_block);
210 for (
int l=1; l<=P; l++) {
212 multiply_block(lower_block_Cijk[l], -1.0, result_block, *rhs_block);
215 divide_diagonal_block(block_indices[l], block_indices[l+1],
216 *rhs_block, result_block);
229 return sg_op->NormInf();
237 return const_cast<char*
>(label.c_str());
251 return sg_op->HasNormInf();
279 const EpetraExt::BlockMultiVector& Input,
280 EpetraExt::BlockMultiVector& Result)
const
284 int m = Input.NumVectors();
289 k_end =
cijk->find_k(sg_basis()->dimension() + 1);
294 int nj =
cijk->num_j(k_it);
299 for (
int mm=0; mm<m; mm++) {
300 j_ptr[l*m+mm] = (*(Input.GetBlock(j)))[mm];
301 mj_indices[l*m+mm] = l*m+mm;
307 (*sg_poly)[k].Apply(input_tmp, result_tmp);
314 double c = value(i_it);
317 for (
int mm=0; mm<m; mm++)
318 (*Result.GetBlock(i))(mm)->Update(alpha*c, *result_tmp(l*m+mm), 1.0);
329 const EpetraExt::BlockMultiVector& Input,
330 EpetraExt::BlockMultiVector& Result)
const
332 for (
int i=row_begin; i<row_end; i++) {
333 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
335 "Stokhos: ASC Deterministic Preconditioner Time");
337 mean_prec->ApplyInverse(*(Input.GetBlock(i)), *(Result.GetBlock(i)));
virtual const Epetra_Comm & Comm() const
Returns a reference to the Epetra_Comm communicator associated with this operator.
int max_num_mat_vec
Maximum number of matvecs in Apply.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
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 ...
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this matrix operator.
virtual ~ApproxSchurComplementPreconditioner()
Destructor.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerking basis.
T & get(ParameterList &l, const std::string &name)
void multiply_block(const Teuchos::RCP< const Stokhos::Sparse3Tensor< int, double > > &cijk, double alpha, const EpetraExt::BlockMultiVector &Input, EpetraExt::BlockMultiVector &Result) const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Cijk_type > Cijk
Pointer to triple product.
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 ...
Bi-directional iterator for traversing a sparse array.
A multidimensional index.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()=0
Get SG polynomial.
bool symmetric
Use symmetric Gauss-Seidel.
Teuchos::Array< Teuchos::RCP< Cijk_type > > upper_block_Cijk
Triple product tensor for each sub-block.
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this matrix operator. ...
bool only_use_linear
Limit Gauss-Seidel loop to linear terms.
Teuchos::Array< Teuchos::RCP< Cijk_type > > lower_block_Cijk
void divide_diagonal_block(int row_begin, int row_end, const EpetraExt::BlockMultiVector &Input, EpetraExt::BlockMultiVector &Result) const
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)
virtual double NormInf() const
Returns an approximate infinity norm of the operator matrix.
ApproxSchurComplementPreconditioner(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 > &base_map, const Teuchos::RCP< const Epetra_Map > &sg_map, const Teuchos::RCP< Stokhos::AbstractPreconditionerFactory > &prec_factory, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor.
virtual int SetUseTranspose(bool UseTranspose)
Set to true if the transpose of the operator is requested.
virtual const char * Label() const
Returns a character string describing the operator.
virtual void setupPreconditioner(const Teuchos::RCP< Stokhos::SGOperator > &sg_op, const Epetra_Vector &x)
Setup preconditioner.
bool scale_op
Flag indicating whether operator be scaled with <^2>
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Stores Epetra Cijk tensor.
std::vector< T >::iterator iterator
int P
Total polynomial order.
Teuchos::Array< int > block_indices
Starting block indices.