203void InvLSCStrategy::initializeState(
const BlockedLinearOp & A,
LSCPrecondState * state)
const
205 Teko_DEBUG_SCOPE(
"InvLSCStrategy::initializeState",10);
208 const LinearOp Bt =
getBlock(0,1,A);
213 LinearOp G = isSymmetric_ ? Bt :
adjoint(D);
215 bool isStabilized = assumeStable_ ? false : (not isZeroOp(C));
223 if(massMatrix_==Teuchos::null) {
224 Teko_DEBUG_MSG(
"LSC::initializeState Build Scaling <F> type \""
225 << getDiagonalName(scaleType_) <<
"\"" ,1);
226 state->invMass_ = getInvDiagonalOp(F,scaleType_);
228 else if(state->invMass_==Teuchos::null) {
229 Teko_DEBUG_MSG(
"LSC::initializeState Build Scaling <mass> type \""
230 << getDiagonalName(scaleType_) <<
"\"" ,1);
231 state->invMass_ = getInvDiagonalOp(massMatrix_,scaleType_);
236 state->BQBt_ = explicitMultiply(B,state->invMass_,Bt,state->BQBt_);
237 Teko_DEBUG_MSG(
"Computed BQBt",10);
240 if(wScaling_!=Teuchos::null && hScaling_==Teuchos::null) {
242 RCP<const Thyra::VectorBase<double> > w = wScaling_->col(0);
243 RCP<const Thyra::VectorBase<double> > iQu
244 = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double> >(state->invMass_)->getDiag();
245 RCP<Thyra::VectorBase<double> > h = Thyra::createMember(iQu->space());
247 Thyra::put_scalar(0.0,h.ptr());
248 Thyra::ele_wise_prod(1.0,*w,*iQu,h.ptr());
249 hScaling_ = Teuchos::rcp(
new Thyra::DefaultDiagonalLinearOp<double>(h));
252 LinearOp H = hScaling_;
253 if(H==Teuchos::null && not isSymmetric_)
258 state->BHBt_ = state->BQBt_;
260 RCP<Teuchos::Time> time = Teuchos::TimeMonitor::getNewTimer(
"InvLSCStrategy::initializeState Build BHBt");
261 Teuchos::TimeMonitor timer(*time);
264 state->BHBt_ = explicitMultiply(D,H,G,state->BHBt_);
268 if(not isStabilized) {
269 state->addInverse(
"BQBtmC",state->BQBt_);
270 state->addInverse(
"BHBtmC",state->BHBt_);
273 state->aiD_ = Teuchos::null;
275 state->setInitialized(
true);
282 if(!Teko::TpetraHelpers::isTpetraLinearOp(F)){
283#ifdef TEKO_HAVE_EPETRA
284 const RCP<const Epetra_Operator> epF = Thyra::get_Epetra_Operator(*F);
285 if(epF!=Teuchos::null && rowZeroingNeeded_) {
287 const RCP<const Epetra_CrsMatrix> crsF = rcp_dynamic_cast<const Epetra_CrsMatrix>(epF);
290 if(crsF!=Teuchos::null) {
291 std::vector<int> zeroIndices;
294 Teko::Epetra::identityRowIndices(crsF->RowMap(), *crsF,zeroIndices);
301 throw std::logic_error(
"InvLSCStrategy::initializeState is trying to use "
302 "Epetra code, but TEKO is not built with Epetra!");
307 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsF = Teko::TpetraHelpers::getTpetraCrsMatrix(F, &scalar, &transp);
309 std::vector<GO> zeroIndices;
312 Teko::TpetraHelpers::identityRowIndices(*crsF->getRowMap(), *crsF,zeroIndices);
315 modF = Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crsF->getDomainMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crsF->getRangeMap()),rcp(
new Teko::TpetraHelpers::ZeroedOperator(zeroIndices,crsF)));
319 Teko_DEBUG_MSG(
"Calculating gamma",10);
320 LinearOp iQuF = multiply(state->invMass_,modF);
323 Teko::LinearOp stabMatrix;
324 state->gamma_ = std::fabs(Teko::computeSpectralRad(iQuF,5e-2,
false,eigSolveParam_))/3.0;
325 Teko_DEBUG_MSG(
"Calculated gamma",10);
326 if(userPresStabMat_!=Teuchos::null) {
327 Teko::LinearOp invDGl = Teko::getInvDiagonalOp(userPresStabMat_);
328 Teko::LinearOp gammaOp = multiply(invDGl,C);
329 state->gamma_ *= std::fabs(Teko::computeSpectralRad(gammaOp,5e-2,
false,eigSolveParam_));
330 stabMatrix = userPresStabMat_;
336 LinearOp invDiagF = getInvDiagonalOp(F);
337 Teko::ModifiableLinearOp modB_idF_Bt = state->getInverse(
"BidFBt");
338 modB_idF_Bt = explicitMultiply(B,invDiagF,Bt,modB_idF_Bt);
339 state->addInverse(
"BidFBt",modB_idF_Bt);
340 const LinearOp B_idF_Bt = modB_idF_Bt;
342 MultiVector vec_D = getDiagonal(B_idF_Bt);
343 update(-1.0,getDiagonal(C),1.0,vec_D);
344 const LinearOp invD = buildInvDiagonal(vec_D,
"inv(D)");
346 Teko_DEBUG_MSG(
"Calculating alpha",10);
347 const LinearOp BidFBtidD = multiply<double>(B_idF_Bt,invD);
348 double num = std::fabs(Teko::computeSpectralRad(BidFBtidD,5e-2,
false,eigSolveParam_));
349 Teko_DEBUG_MSG(
"Calculated alpha",10);
350 state->alpha_ = 1.0/num;
351 state->aiD_ = Thyra::scale(state->alpha_,invD);
354 Teko::ModifiableLinearOp BQBtmC = state->getInverse(
"BQBtmC");
355 BQBtmC = explicitAdd(state->BQBt_,scale(-state->gamma_,stabMatrix),BQBtmC);
356 state->addInverse(
"BQBtmC",BQBtmC);
359 Teko::ModifiableLinearOp BHBtmC = state->getInverse(
"BHBtmC");
363 BHBtmC = explicitAdd(state->BHBt_,scale(-state->gamma_,stabMatrix),BHBtmC);
365 state->addInverse(
"BHBtmC",BHBtmC);
367 Teko_DEBUG_MSG_BEGIN(5)
368 DEBUG_STREAM <<
"LSC Gamma Parameter = " << state->gamma_ << std::endl;
369 DEBUG_STREAM <<
"LSC Alpha Parameter = " << state->alpha_ << std::endl;
372 state->setInitialized(
true);
380void InvLSCStrategy::computeInverses(
const BlockedLinearOp & A,
LSCPrecondState * state)
const
382 Teko_DEBUG_SCOPE(
"InvLSCStrategy::computeInverses",10);
383 Teko_DEBUG_EXPR(Teuchos::Time invTimer(
""));
390 Teko_DEBUG_MSG(
"LSC::computeInverses Building inv(F)",1);
391 Teko_DEBUG_EXPR(invTimer.start(
true));
392 InverseLinearOp invF = state->getInverse(
"invF");
393 if(invF==Teuchos::null) {
395 state->addInverse(
"invF",invF);
399 Teko_DEBUG_EXPR(invTimer.stop());
400 Teko_DEBUG_MSG(
"LSC::computeInverses GetInvF = " << invTimer.totalElapsedTime(),1);
405 Teko_DEBUG_MSG(
"LSC::computeInverses Building inv(BQBtmC)",1);
406 Teko_DEBUG_EXPR(invTimer.start(
true));
407 const LinearOp BQBt = state->getInverse(
"BQBtmC");
408 InverseLinearOp invBQBt = state->getInverse(
"invBQBtmC");
409 if(invBQBt==Teuchos::null) {
411 state->addInverse(
"invBQBtmC",invBQBt);
415 Teko_DEBUG_EXPR(invTimer.stop());
416 Teko_DEBUG_MSG(
"LSC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(),1);
421 ModifiableLinearOp invBHBt = state->getInverse(
"invBHBtmC");
422 if(hScaling_!=Teuchos::null || not isSymmetric_) {
424 Teko_DEBUG_MSG(
"LSC::computeInverses Building inv(BHBtmC)",1);
425 Teko_DEBUG_EXPR(invTimer.start(
true));
426 const LinearOp BHBt = state->getInverse(
"BHBtmC");
427 if(invBHBt==Teuchos::null) {
429 state->addInverse(
"invBHBtmC",invBHBt);
433 Teko_DEBUG_EXPR(invTimer.stop());
434 Teko_DEBUG_MSG(
"LSC::computeInverses GetInvBHBt = " << invTimer.totalElapsedTime(),1);
436 else if(invBHBt==Teuchos::null) {
438 state->addInverse(
"invBHBtmC",invBQBt);
443void InvLSCStrategy::initializeFromParameterList(
const Teuchos::ParameterList & pl,
const InverseLibrary & invLib)
446 std::string invStr=
"", invVStr=
"", invPStr=
"";
447 bool rowZeroing =
true;
452 if(pl.isParameter(
"Inverse Type"))
453 invStr = pl.get<std::string>(
"Inverse Type");
454 if(pl.isParameter(
"Inverse Velocity Type"))
455 invVStr = pl.get<std::string>(
"Inverse Velocity Type");
456 if(pl.isParameter(
"Inverse Pressure Type"))
457 invPStr = pl.get<std::string>(
"Inverse Pressure Type");
458 if(pl.isParameter(
"Ignore Boundary Rows"))
459 rowZeroing = pl.get<
bool>(
"Ignore Boundary Rows");
460 if(pl.isParameter(
"Use LDU"))
461 useLDU = pl.get<
bool>(
"Use LDU");
462 if(pl.isParameter(
"Use Mass Scaling"))
463 useMass_ = pl.get<
bool>(
"Use Mass Scaling");
466 if(pl.isParameter(
"Use W-Scaling"))
467 useWScaling_ = pl.get<
bool>(
"Use W-Scaling");
468 if(pl.isParameter(
"Eigen Solver Iterations"))
469 eigSolveParam_ = pl.get<
int>(
"Eigen Solver Iterations");
470 if(pl.isParameter(
"Scaling Type")) {
471 scaleType_ = getDiagonalType(pl.get<std::string>(
"Scaling Type"));
472 TEUCHOS_TEST_FOR_EXCEPT(scaleType_==
NotDiag);
474 if(pl.isParameter(
"Assume Stable Discretization"))
475 assumeStable_ = pl.get<
bool>(
"Assume Stable Discretization");
477 Teko_DEBUG_MSG_BEGIN(5)
478 DEBUG_STREAM <<
"LSC Inverse Strategy Parameters: " << std::endl;
479 DEBUG_STREAM <<
" inv type = \"" << invStr <<
"\"" << std::endl;
480 DEBUG_STREAM <<
" inv v type = \"" << invVStr <<
"\"" << std::endl;
481 DEBUG_STREAM <<
" inv p type = \"" << invPStr <<
"\"" << std::endl;
482 DEBUG_STREAM <<
" bndry rows = " << rowZeroing << std::endl;
483 DEBUG_STREAM <<
" use ldu = " << useLDU << std::endl;
484 DEBUG_STREAM <<
" use mass = " << useMass_ << std::endl;
485 DEBUG_STREAM <<
" use w-scaling = " << useWScaling_ << std::endl;
486 DEBUG_STREAM <<
" assume stable = " << assumeStable_ << std::endl;
487 DEBUG_STREAM <<
" scale type = " << getDiagonalName(scaleType_) << std::endl;
488 DEBUG_STREAM <<
"LSC Inverse Strategy Parameter list: " << std::endl;
489 pl.print(DEBUG_STREAM);
493 if(invStr==
"") invStr =
"Amesos";
494 if(invVStr==
"") invVStr = invStr;
495 if(invPStr==
"") invPStr = invStr;
498 invFactoryF_ = invLib.getInverseFactory(invVStr);
499 invFactoryS_ = invFactoryF_;
501 invFactoryS_ = invLib.getInverseFactory(invPStr);
504 setUseFullLDU(useLDU);
505 setRowZeroing(rowZeroing);
508 Teuchos::RCP<Teko::RequestHandler> rh = getRequestHandler();
509 rh->preRequest<Teko::LinearOp>(Teko::RequestMesg(
"Velocity Mass Matrix"));
511 = rh->request<Teko::LinearOp>(Teko::RequestMesg(
"Velocity Mass Matrix"));