49 #include "NOX_Epetra_Interface_Jacobian.H"
54 NOX::Epetra::LinearSystemMPBD::
64 block_solver(block_solver_),
65 jacInterfacePtr(iJac),
72 num_mp_blocks = block_ops->size();
74 std::string prec_strategy = linearSolverParams.
get(
"Preconditioner Strategy",
76 if (prec_strategy ==
"Standard")
77 precStrategy = STANDARD;
78 else if (prec_strategy ==
"Mean")
80 else if (prec_strategy ==
"On the fly")
81 precStrategy = ON_THE_FLY;
84 "Invalid preconditioner strategy " << prec_strategy);
86 if (precStrategy == STANDARD)
87 precs.resize(num_mp_blocks);
90 NOX::Epetra::LinearSystemMPBD::
95 bool NOX::Epetra::LinearSystemMPBD::
96 applyJacobian(
const NOX::Epetra::Vector& input,
97 NOX::Epetra::Vector& result)
const
99 mp_op->SetUseTranspose(
false);
100 int status = mp_op->Apply(input.getEpetraVector(), result.getEpetraVector());
102 return (status == 0);
105 bool NOX::Epetra::LinearSystemMPBD::
106 applyJacobianTranspose(
const NOX::Epetra::Vector& input,
107 NOX::Epetra::Vector& result)
const
109 mp_op->SetUseTranspose(
true);
110 int status = mp_op->Apply(input.getEpetraVector(), result.getEpetraVector());
111 mp_op->SetUseTranspose(
false);
113 return (status == 0);
116 bool NOX::Epetra::LinearSystemMPBD::
118 const NOX::Epetra::Vector &input,
119 NOX::Epetra::Vector &result)
124 EpetraExt::BlockVector input_block(View, *base_map,
125 input.getEpetraVector());
126 EpetraExt::BlockVector result_block(View, *base_map,
127 result.getEpetraVector());
128 result_block.PutScalar(0.0);
132 params.
sublist(
"Deterministic Solver Parameters");
135 bool final_status =
true;
137 for (
int i=0; i<num_mp_blocks; i++) {
138 NOX::Epetra::Vector nox_input(input_block.GetBlock(i),
139 NOX::Epetra::Vector::CreateView);
140 NOX::Epetra::Vector nox_result(result_block.GetBlock(i),
141 NOX::Epetra::Vector::CreateView);
143 block_solver->setJacobianOperatorForSolve(block_ops->getCoeffPtr(i));
145 if (precStrategy == STANDARD)
146 block_solver->setPrecOperatorForSolve(precs[i]);
147 else if (precStrategy == ON_THE_FLY) {
148 block_solver->createPreconditioner(*(prec_x->GetBlock(i)),
149 block_solver_params,
false);
152 status = block_solver->applyJacobianInverse(block_solver_params, nox_input,
154 final_status = final_status && status;
160 bool NOX::Epetra::LinearSystemMPBD::
161 applyRightPreconditioning(
bool useTranspose,
163 const NOX::Epetra::Vector& input,
164 NOX::Epetra::Vector& result)
const
175 void NOX::Epetra::LinearSystemMPBD::
182 bool NOX::Epetra::LinearSystemMPBD::
183 computeJacobian(
const NOX::Epetra::Vector& x)
185 bool success = jacInterfacePtr->computeJacobian(x.getEpetraVector(),
187 block_ops = mp_op->getMPOps();
192 bool NOX::Epetra::LinearSystemMPBD::
193 createPreconditioner(
const NOX::Epetra::Vector& x,
195 bool recomputeGraph)
const
197 EpetraExt::BlockVector mp_x_block(View, *base_map, x.getEpetraVector());
199 p.
sublist(
"Deterministic Solver Parameters");
200 bool total_success =
true;
201 if (precStrategy == STANDARD) {
202 for (
int i=0; i<num_mp_blocks; i++) {
203 block_solver->setJacobianOperatorForSolve(block_ops->getCoeffPtr(i));
205 block_solver->setPrecOperatorForSolve(precs[i]);
207 block_solver->createPreconditioner(*(mp_x_block.GetBlock(i)),
208 solverParams, recomputeGraph);
209 precs[i] = block_solver->getGeneratedPrecOperator();
210 total_success = total_success && success;
214 else if (precStrategy ==
MEAN) {
215 block_solver->setJacobianOperatorForSolve(block_ops->getCoeffPtr(0));
217 block_solver->createPreconditioner(*(mp_x_block.GetBlock(0)),
218 solverParams, recomputeGraph);
219 total_success = total_success && success;
222 else if (precStrategy == ON_THE_FLY) {
224 prec_x =
Teuchos::rcp(
new EpetraExt::BlockVector(mp_x_block));
226 *prec_x = mp_x_block;
229 return total_success;
232 bool NOX::Epetra::LinearSystemMPBD::
233 destroyPreconditioner()
const
235 return block_solver->destroyPreconditioner();
238 bool NOX::Epetra::LinearSystemMPBD::
239 recomputePreconditioner(
const NOX::Epetra::Vector& x,
242 EpetraExt::BlockVector mp_x_block(View, *base_map, x.getEpetraVector());
244 p.
sublist(
"Deterministic Solver Parameters");
245 bool total_success =
true;
247 if (precStrategy == STANDARD) {
248 for (
int i=0; i<num_mp_blocks; i++) {
249 block_solver->setJacobianOperatorForSolve(block_ops->getCoeffPtr(i));
251 block_solver->setPrecOperatorForSolve(precs[i]);
253 block_solver->recomputePreconditioner(*(mp_x_block.GetBlock(i)),
255 precs[i] = block_solver->getGeneratedPrecOperator();
256 total_success = total_success && success;
260 else if (precStrategy ==
MEAN) {
261 block_solver->setJacobianOperatorForSolve(block_ops->getCoeffPtr(0));
263 block_solver->recomputePreconditioner(*(mp_x_block.GetBlock(0)),
265 total_success = total_success && success;
268 else if (precStrategy == ON_THE_FLY) {
270 prec_x =
Teuchos::rcp(
new EpetraExt::BlockVector(mp_x_block));
272 *prec_x = mp_x_block;
275 return total_success;
278 NOX::Epetra::LinearSystem::PreconditionerReusePolicyType
279 NOX::Epetra::LinearSystemMPBD::
280 getPreconditionerPolicy(
bool advanceReuseCounter)
282 return block_solver->getPreconditionerPolicy(advanceReuseCounter);
285 bool NOX::Epetra::LinearSystemMPBD::
286 isPreconditionerConstructed()
const
288 return block_solver->isPreconditionerConstructed();
291 bool NOX::Epetra::LinearSystemMPBD::
292 hasPreconditioner()
const
294 return block_solver->hasPreconditioner();
298 getJacobianOperator()
const
304 getJacobianOperator()
310 getGeneratedPrecOperator()
const
312 if (precStrategy ==
MEAN)
313 return block_solver->getGeneratedPrecOperator();
316 true, std::logic_error,
317 "Cannot call getGeneratedPrecOperator() unless prec strategy is Mean");
322 getGeneratedPrecOperator()
324 if (precStrategy ==
MEAN)
325 return block_solver->getGeneratedPrecOperator();
328 true, std::logic_error,
329 "Cannot call getGeneratedPrecOperator() unless prec strategy is Mean");
333 void NOX::Epetra::LinearSystemMPBD::
334 setJacobianOperatorForSolve(
const
344 void NOX::Epetra::LinearSystemMPBD::
347 if (precStrategy ==
MEAN)
348 block_solver->setPrecOperatorForSolve(solvePrecOp);
351 true, std::logic_error,
352 "Cannot call setPrecOperatorForSolve() unless prec strategy is Mean");
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
T & get(ParameterList &l, const std::string &name)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< ParameterList > sublist(const RCP< ParameterList > ¶mList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
An Epetra operator representing the block stochastic Galerkin operator.
virtual Teuchos::RCP< Stokhos::ProductEpetraOperator > getMPOps()
Get multi-point ops.