57 label(
"Stokhos Approximate Schur Complement Preconditioner"),
60 epetraCijk(epetraCijk_),
63 prec_factory(prec_factory_),
68 Cijk(epetraCijk->getParallelCijk()),
71 upper_block_Cijk(P+1),
72 lower_block_Cijk(P+1),
75 only_use_linear(
true),
81 "Stokhos::ApproxSchurComplementPreconditioner does not support " <<
82 "a parallel stochastic distribution.");
84 scale_op = params_->
get(
"Scale Operator by Inverse Basis Norms",
true);
85 symmetric = params_->
get(
"Symmetric Gauss-Seidel",
false);
95 int nj =
Cijk->num_j(k);
103 int d = prod_basis->dimension();
105 for (
int p=0; p<=
P; p++) {
132 double c = value(i_it);
141 else if (p_col < p_row) {
148 for (
int p=0; p<=
P; p++) {
166 mean_prec = prec_factory->compute(sg_poly->getCoeffPtr(0));
167 label = std::string(
"Stokhos Approximate Schur Complement Preconditioner:\n")
168 + std::string(
" ***** ") + std::string(mean_prec->Label());
175 useTranspose = UseTranspose;
176 sg_op->SetUseTranspose(UseTranspose);
177 mean_prec->SetUseTranspose(UseTranspose);
186 return sg_op->Apply(Input, Result);
193 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
200 bool made_copy =
false;
208 if (rhs_block == Teuchos::null || rhs_block->NumVectors() != m)
210 Teuchos::rcp(
new EpetraExt::BlockMultiVector(*base_map, *sg_map, m));
211 if (tmp == Teuchos::null || tmp->NumVectors() != m*max_num_mat_vec)
214 j_ptr.resize(m*max_num_mat_vec);
215 mj_indices.resize(m*max_num_mat_vec);
218 EpetraExt::BlockMultiVector input_block(
View, *base_map, *input);
219 EpetraExt::BlockMultiVector result_block(
View, *base_map, Result);
221 result_block.PutScalar(0.0);
224 rhs_block->Update(1.0, input_block, 0.0);
230 for (
int l=P; l>=1; l--) {
232 divide_diagonal_block(block_indices[l], block_indices[l+1],
233 *rhs_block, result_block);
236 multiply_block(upper_block_Cijk[l], -1.0, result_block, *rhs_block);
240 divide_diagonal_block(0, 1, *rhs_block, result_block);
242 for (
int l=1; l<=P; l++) {
244 multiply_block(lower_block_Cijk[l], -1.0, result_block, *rhs_block);
247 divide_diagonal_block(block_indices[l], block_indices[l+1],
248 *rhs_block, result_block);
261 return sg_op->NormInf();
269 return const_cast<char*>(label.c_str());
283 return sg_op->HasNormInf();
311 const EpetraExt::BlockMultiVector& Input,
312 EpetraExt::BlockMultiVector& Result)
const
316 int m = Input.NumVectors();
321 k_end =
cijk->find_k(sg_basis()->dimension() + 1);
326 int nj =
cijk->num_j(k_it);
331 for (
int mm=0; mm<m; mm++) {
332 j_ptr[l*m+mm] = (*(Input.GetBlock(
j)))[mm];
333 mj_indices[l*m+mm] = l*m+mm;
339 (*sg_poly)[k].Apply(input_tmp, result_tmp);
346 double c = value(i_it);
349 for (
int mm=0; mm<m; mm++)
350 (*Result.GetBlock(i))(mm)->Update(alpha*c, *result_tmp(l*m+mm), 1.0);
361 const EpetraExt::BlockMultiVector& Input,
362 EpetraExt::BlockMultiVector& Result)
const
364 for (
int i=row_begin; i<row_end; i++) {
365 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
367 "Stokhos: ASC Deterministic Preconditioner Time");
369 mean_prec->ApplyInverse(*(Input.GetBlock(i)), *(Result.GetBlock(i)));