00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047 #include "Teko_LSCSIMPLECStrategy.hpp"
00048
00049 #include "Thyra_DefaultDiagonalLinearOp.hpp"
00050 #include "Thyra_EpetraThyraWrappers.hpp"
00051 #include "Thyra_get_Epetra_Operator.hpp"
00052 #include "Thyra_EpetraLinearOp.hpp"
00053 #include "Thyra_VectorStdOps.hpp"
00054
00055 #include "Epetra_Vector.h"
00056 #include "Epetra_Map.h"
00057
00058 #include "EpetraExt_RowMatrixOut.h"
00059 #include "EpetraExt_MultiVectorOut.h"
00060 #include "EpetraExt_VectorOut.h"
00061
00062 #include "Teuchos_Time.hpp"
00063 #include "Teuchos_TimeMonitor.hpp"
00064
00065
00066 #include "Teko_Utilities.hpp"
00067 #include "Teko_LSCPreconditionerFactory.hpp"
00068 #include "Teko_EpetraHelpers.hpp"
00069 #include "Teko_EpetraOperatorWrapper.hpp"
00070
00071 using Teuchos::RCP;
00072 using Teuchos::rcp_dynamic_cast;
00073 using Teuchos::rcp_const_cast;
00074
00075 namespace Teko {
00076 namespace NS {
00077
00079
00081
00082
00084 LSCSIMPLECStrategy::LSCSIMPLECStrategy()
00085 : invFactoryF_(Teuchos::null), invFactoryS_(Teuchos::null)
00086 , useFullLDU_(false), scaleType_(Diagonal)
00087 { }
00088
00090
00091 void LSCSIMPLECStrategy::buildState(BlockedLinearOp & A,BlockPreconditionerState & state) const
00092 {
00093 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::buildState",10);
00094
00095 LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
00096 TEUCHOS_ASSERT(lscState!=0);
00097
00098
00099 if(not lscState->isInitialized()) {
00100 Teko_DEBUG_EXPR(Teuchos::Time timer(""));
00101
00102
00103 {
00104 Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState constructing operators",1);
00105 Teko_DEBUG_EXPR(timer.start(true));
00106
00107 initializeState(A,lscState);
00108
00109 Teko_DEBUG_EXPR(timer.stop());
00110 Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildOpsTime = " << timer.totalElapsedTime(),1);
00111 }
00112
00113
00114 {
00115 Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState calculating inverses",1);
00116 Teko_DEBUG_EXPR(timer.start(true));
00117
00118 computeInverses(A,lscState);
00119
00120 Teko_DEBUG_EXPR(timer.stop());
00121 Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildInvTime = " << timer.totalElapsedTime(),1);
00122 }
00123 }
00124 }
00125
00126
00127 LinearOp LSCSIMPLECStrategy::getInvBQBt(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00128 {
00129 return state.getInverse("invBQBtmC");
00130 }
00131
00132 LinearOp LSCSIMPLECStrategy::getInvBHBt(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00133 {
00134 return state.getInverse("invBQBtmC").getConst();
00135 }
00136
00137 LinearOp LSCSIMPLECStrategy::getInvF(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00138 {
00139 return state.getInverse("invF");
00140 }
00141
00142 LinearOp LSCSIMPLECStrategy::getInnerStabilization(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00143 {
00144 LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
00145 TEUCHOS_ASSERT(lscState!=0);
00146 TEUCHOS_ASSERT(lscState->isInitialized())
00147
00148 const LinearOp C = getBlock(1,1,A);
00149 return scale(-1.0,C);
00150 }
00151
00152
00153 LinearOp LSCSIMPLECStrategy::getInvMass(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00154 {
00155 LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
00156 TEUCHOS_ASSERT(lscState!=0);
00157 TEUCHOS_ASSERT(lscState->isInitialized())
00158
00159 return lscState->invMass_;
00160 }
00161
00162 LinearOp LSCSIMPLECStrategy::getHScaling(const BlockedLinearOp & A,BlockPreconditionerState & state) const
00163 {
00164 return getInvMass(A,state);
00165 }
00166
00168 void LSCSIMPLECStrategy::initializeState(const BlockedLinearOp & A,LSCPrecondState * state) const
00169 {
00170 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::initializeState",10);
00171
00172 const LinearOp F = getBlock(0,0,A);
00173 const LinearOp Bt = getBlock(0,1,A);
00174 const LinearOp B = getBlock(1,0,A);
00175 const LinearOp C = getBlock(1,1,A);
00176
00177 bool isStabilized = (not isZeroOp(C));
00178
00179 state->invMass_ = getInvDiagonalOp(F,scaleType_);
00180
00181
00182 state->BQBt_ = explicitMultiply(B,state->invMass_,Bt,state->BQBt_);
00183 if(isStabilized) {
00184
00185 Teko::ModifiableLinearOp BQBtmC = state->getInverse("BQBtmC");
00186 BQBtmC = explicitAdd(state->BQBt_,scale(-1.0,C),BQBtmC);
00187 state->addInverse("BQBtmC",BQBtmC);
00188 }
00189 Teko_DEBUG_MSG("Computed BQBt",10);
00190
00191 state->setInitialized(true);
00192 }
00193
00199 void LSCSIMPLECStrategy::computeInverses(const BlockedLinearOp & A,LSCPrecondState * state) const
00200 {
00201 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::computeInverses",10);
00202 Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
00203
00204 const LinearOp F = getBlock(0,0,A);
00205
00207
00208
00209 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(F)",1);
00210 Teko_DEBUG_EXPR(invTimer.start(true));
00211 InverseLinearOp invF = state->getInverse("invF");
00212 if(invF==Teuchos::null) {
00213 invF = buildInverse(*invFactoryF_,F);
00214 state->addInverse("invF",invF);
00215 } else {
00216 rebuildInverse(*invFactoryF_,F,invF);
00217 }
00218 Teko_DEBUG_EXPR(invTimer.stop());
00219 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvF = " << invTimer.totalElapsedTime(),1);
00220
00222
00223
00224 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(BQBtmC)",1);
00225 Teko_DEBUG_EXPR(invTimer.start(true));
00226 const LinearOp BQBt = state->getInverse("BQBtmC");
00227 InverseLinearOp invBQBt = state->getInverse("invBQBtmC");
00228 if(invBQBt==Teuchos::null) {
00229 invBQBt = buildInverse(*invFactoryS_,BQBt);
00230 state->addInverse("invBQBtmC",invBQBt);
00231 } else {
00232 rebuildInverse(*invFactoryS_,BQBt,invBQBt);
00233 }
00234 Teko_DEBUG_EXPR(invTimer.stop());
00235 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(),1);
00236 }
00237
00239 void LSCSIMPLECStrategy::initializeFromParameterList(const Teuchos::ParameterList & pl,const InverseLibrary & invLib)
00240 {
00241
00242 std::string invStr="", invVStr="", invPStr="";
00243 bool useLDU = false;
00244 scaleType_ = Diagonal;
00245
00246
00247 if(pl.isParameter("Inverse Type"))
00248 invStr = pl.get<std::string>("Inverse Type");
00249 if(pl.isParameter("Inverse Velocity Type"))
00250 invVStr = pl.get<std::string>("Inverse Velocity Type");
00251 if(pl.isParameter("Inverse Pressure Type"))
00252 invPStr = pl.get<std::string>("Inverse Pressure Type");
00253 if(pl.isParameter("Use LDU"))
00254 useLDU = pl.get<bool>("Use LDU");
00255 if(pl.isParameter("Scaling Type")) {
00256 scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
00257 TEUCHOS_TEST_FOR_EXCEPT(scaleType_==NotDiag);
00258 }
00259
00260 Teko_DEBUG_MSG_BEGIN(0)
00261 DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
00262 DEBUG_STREAM << " inv type = \"" << invStr << "\"" << std::endl;
00263 DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
00264 DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
00265 DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
00266 DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
00267 DEBUG_STREAM << "LSC Inverse Strategy Parameter list: " << std::endl;
00268 pl.print(DEBUG_STREAM);
00269 Teko_DEBUG_MSG_END()
00270
00271
00272 if(invStr=="") invStr = "Amesos";
00273 if(invVStr=="") invVStr = invStr;
00274 if(invPStr=="") invPStr = invStr;
00275
00276
00277 invFactoryF_ = invLib.getInverseFactory(invVStr);
00278 invFactoryS_ = invFactoryF_;
00279 if(invVStr!=invPStr)
00280 invFactoryS_ = invLib.getInverseFactory(invPStr);
00281
00282
00283 setUseFullLDU(useLDU);
00284 }
00285
00286 }
00287 }