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_GaussSeidelPreconditionerFactory.hpp"
00048
00049 #include "Teko_BlockUpperTriInverseOp.hpp"
00050 #include "Teko_BlockLowerTriInverseOp.hpp"
00051
00052 using Teuchos::rcp;
00053 using Teuchos::RCP;
00054
00055 namespace Teko {
00056
00057 GaussSeidelPreconditionerFactory::GaussSeidelPreconditionerFactory(TriSolveType solveType,const LinearOp & invD0,const LinearOp & invD1)
00058 : invOpsStrategy_(rcp(new StaticInvDiagStrategy(invD0,invD1))), solveType_(solveType)
00059 { }
00060
00061 GaussSeidelPreconditionerFactory::GaussSeidelPreconditionerFactory(TriSolveType solveType,const RCP<const BlockInvDiagonalStrategy> & strategy)
00062 : invOpsStrategy_(strategy), solveType_(solveType)
00063 { }
00064
00065 GaussSeidelPreconditionerFactory::GaussSeidelPreconditionerFactory()
00066 : solveType_(GS_UseLowerTriangle)
00067 { }
00068
00069 LinearOp GaussSeidelPreconditionerFactory::buildPreconditionerOperator(BlockedLinearOp & blo,BlockPreconditionerState & state) const
00070 {
00071 int rows = blockRowCount(blo);
00072 int cols = blockColCount(blo);
00073
00074 TEUCHOS_ASSERT(rows==cols);
00075
00076
00077 std::vector<LinearOp> invDiag;
00078 invOpsStrategy_->getInvD(blo,state,invDiag);
00079 TEUCHOS_ASSERT(rows==(int) invDiag.size());
00080
00081 if(solveType_==GS_UseUpperTriangle) {
00082
00083 BlockedLinearOp U = getUpperTriBlocks(blo);
00084
00085 return createBlockUpperTriInverseOp(U,invDiag,"Gauss Seidel");
00086 }
00087 else if(solveType_==GS_UseLowerTriangle) {
00088
00089 BlockedLinearOp L = getLowerTriBlocks(blo);
00090
00091 return createBlockLowerTriInverseOp(L,invDiag,"Gauss Seidel");
00092 }
00093
00094 TEUCHOS_ASSERT(false);
00095 }
00096
00098 void GaussSeidelPreconditionerFactory::initializeFromParameterList(const Teuchos::ParameterList & pl)
00099 {
00100 Teko_DEBUG_SCOPE("GaussSeidelPreconditionerFactory::initializeFromParameterList",10);
00101 Teko_DEBUG_MSG_BEGIN(9);
00102 DEBUG_STREAM << "Parameter list: " << std::endl;
00103 pl.print(DEBUG_STREAM);
00104 Teko_DEBUG_MSG_END();
00105
00106 const std::string inverse_type = "Inverse Type";
00107 std::vector<RCP<InverseFactory> > inverses;
00108
00109 RCP<const InverseLibrary> invLib = getInverseLibrary();
00110
00111
00112 std::string invStr ="Amesos";
00113 if(pl.isParameter(inverse_type))
00114 invStr = pl.get<std::string>(inverse_type);
00115 if(pl.isParameter("Use Upper Triangle"))
00116 solveType_ = pl.get<bool>("Use Upper Triangle") ? GS_UseUpperTriangle : GS_UseLowerTriangle;
00117
00118 Teko_DEBUG_MSG("GSPrecFact: Building default inverse \"" << invStr << "\"",5);
00119 RCP<InverseFactory> defaultInverse = invLib->getInverseFactory(invStr);
00120
00121
00122 Teuchos::ParameterList::ConstIterator itr;
00123 for(itr=pl.begin();itr!=pl.end();++itr) {
00124 std::string fieldName = itr->first;
00125 Teko_DEBUG_MSG("GSPrecFact: checking fieldName = \"" << fieldName << "\"",9);
00126
00127
00128 if(fieldName.compare(0,inverse_type.length(),inverse_type)==0 && fieldName!=inverse_type) {
00129 int position = -1;
00130 std::string inverse,type;
00131
00132
00133 std::stringstream ss(fieldName);
00134 ss >> inverse >> type >> position;
00135
00136 if(position<=0) {
00137 Teko_DEBUG_MSG("Gauss-Seidel \"Inverse Type\" must be a (strictly) positive integer",1);
00138 }
00139
00140
00141 std::string invStr = pl.get<std::string>(fieldName);
00142 Teko_DEBUG_MSG("GSPrecFact: Building inverse " << position << " \"" << invStr << "\"",5);
00143 if(position>(int) inverses.size()) {
00144 inverses.resize(position,defaultInverse);
00145 inverses[position-1] = invLib->getInverseFactory(invStr);
00146 }
00147 else
00148 inverses[position-1] = invLib->getInverseFactory(invStr);
00149 }
00150 }
00151
00152
00153 if(inverses.size()==0)
00154 inverses.push_back(defaultInverse);
00155
00156
00157 invOpsStrategy_ = rcp(new InvFactoryDiagStrategy(inverses,defaultInverse));
00158 }
00159
00160 }