Teko Version of the Day
Teko_JacobiPreconditionerFactory.cpp
00001 /*
00002 // @HEADER
00003 // 
00004 // ***********************************************************************
00005 // 
00006 //      Teko: A package for block and physics based preconditioning
00007 //                  Copyright 2010 Sandia Corporation 
00008 //  
00009 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00010 // the U.S. Government retains certain rights in this software.
00011 //  
00012 // Redistribution and use in source and binary forms, with or without
00013 // modification, are permitted provided that the following conditions are
00014 // met:
00015 //  
00016 // 1. Redistributions of source code must retain the above copyright
00017 // notice, this list of conditions and the following disclaimer.
00018 //  
00019 // 2. Redistributions in binary form must reproduce the above copyright
00020 // notice, this list of conditions and the following disclaimer in the
00021 // documentation and/or other materials provided with the distribution.
00022 //  
00023 // 3. Neither the name of the Corporation nor the names of the
00024 // contributors may be used to endorse or promote products derived from
00025 // this software without specific prior written permission. 
00026 //  
00027 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00028 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00029 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00030 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00031 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00032 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00033 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00034 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00035 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 
00036 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00037 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00038 //  
00039 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
00040 // 
00041 // ***********************************************************************
00042 // 
00043 // @HEADER
00044 
00045 */
00046 
00047 #include "Teko_JacobiPreconditionerFactory.hpp"
00048 
00049 using Teuchos::rcp;
00050 
00051 namespace Teko {
00052 
00053 JacobiPreconditionerFactory::JacobiPreconditionerFactory(const LinearOp & invD0,const LinearOp & invD1)
00054       : invOpsStrategy_(rcp(new StaticInvDiagStrategy(invD0,invD1)))
00055 { }
00056 
00057 JacobiPreconditionerFactory::JacobiPreconditionerFactory(const RCP<const BlockInvDiagonalStrategy> & strategy)
00058          : invOpsStrategy_(strategy)
00059 { }
00060 
00063 JacobiPreconditionerFactory::JacobiPreconditionerFactory()
00064 { }
00065 
00066 LinearOp JacobiPreconditionerFactory::buildPreconditionerOperator(BlockedLinearOp & blo,BlockPreconditionerState & state) const
00067 {
00068    int rows = blo->productRange()->numBlocks();
00069    int cols = blo->productDomain()->numBlocks();
00070  
00071    TEUCHOS_ASSERT(rows==cols);
00072 
00073    // get diagonal blocks
00074    std::vector<LinearOp> invDiag;
00075    invOpsStrategy_->getInvD(blo,state,invDiag);
00076    TEUCHOS_ASSERT(rows==(int) invDiag.size());
00077 
00078    // create a blocked linear operator
00079    BlockedLinearOp precond = createBlockedOp();
00080    std::stringstream ss;
00081    ss << "Jacobi Preconditioner ( ";
00082 
00083    // start filling the blocked operator
00084    beginBlockFill(precond,rows,rows); // this is assuming the matrix is square
00085 
00086    // build blocked diagonal matrix
00087    for(int i=0;i<rows;i++) {
00088       ss << " op" << i << " = " << invDiag[i]->description() << ", ";
00089       precond->setBlock(i,i,invDiag[i]);
00090    }
00091    ss << " )";
00092    
00093    endBlockFill(precond);
00094    // done filling the blocked operator
00095 
00096    // precond->setObjectLabel(ss.str());
00097    precond->setObjectLabel("Jacobi");
00098    
00099    return precond; 
00100 }
00101 
00103 void JacobiPreconditionerFactory::initializeFromParameterList(const Teuchos::ParameterList & pl)
00104 #if 0
00105 {
00106    RCP<const InverseLibrary> invLib = getInverseLibrary();
00107 
00108    // get string specifying inverse
00109    std::string invStr = pl.get<std::string>("Inverse Type");
00110    if(invStr=="") invStr = "Amesos";
00111 
00112    // based on parameter type build a strategy
00113    invOpsStrategy_ = rcp(new InvFactoryDiagStrategy(invLib->getInverseFactory(invStr)));
00114 }
00115 #endif 
00116 {
00117    Teko_DEBUG_SCOPE("JacobiPreconditionerFactory::initializeFromParameterList",10);
00118    Teko_DEBUG_MSG_BEGIN(9);
00119       DEBUG_STREAM << "Parameter list: " << std::endl;
00120       pl.print(DEBUG_STREAM);
00121    Teko_DEBUG_MSG_END();
00122 
00123    const std::string inverse_type = "Inverse Type";
00124    std::vector<RCP<InverseFactory> > inverses;
00125 
00126    RCP<const InverseLibrary> invLib = getInverseLibrary();
00127 
00128    // get string specifying default inverse
00129    std::string invStr ="Amesos"; 
00130    if(pl.isParameter(inverse_type))
00131       invStr = pl.get<std::string>(inverse_type);
00132 
00133    Teko_DEBUG_MSG("JacobiPrecFact: Building default inverse \"" << invStr << "\"",5);
00134    RCP<InverseFactory> defaultInverse = invLib->getInverseFactory(invStr);
00135 
00136    // now check individual solvers
00137    Teuchos::ParameterList::ConstIterator itr;
00138    for(itr=pl.begin();itr!=pl.end();++itr) {
00139       std::string fieldName = itr->first;
00140       Teko_DEBUG_MSG("JacobiPrecFact: checking fieldName = \"" << fieldName << "\"",9);
00141 
00142       // figure out what the integer is
00143       if(fieldName.compare(0,inverse_type.length(),inverse_type)==0 && fieldName!=inverse_type) {
00144          int position = -1;
00145          std::string inverse,type;
00146 
00147          // figure out position
00148          std::stringstream ss(fieldName);
00149          ss >> inverse >> type >> position;
00150 
00151          if(position<=0) {
00152             Teko_DEBUG_MSG("Jacobi \"Inverse Type\" must be a (strictly) positive integer",1);
00153          }
00154 
00155          // inserting inverse factory into vector
00156          std::string invStr = pl.get<std::string>(fieldName);
00157          Teko_DEBUG_MSG("JacobiPrecFact: Building inverse " << position << " \"" << invStr << "\"",5);
00158          if(position>(int) inverses.size()) {
00159             inverses.resize(position,defaultInverse);
00160             inverses[position-1] = invLib->getInverseFactory(invStr);
00161          }
00162          else
00163             inverses[position-1] = invLib->getInverseFactory(invStr);
00164       }
00165    }
00166 
00167    // use default inverse
00168    if(inverses.size()==0) 
00169       inverses.push_back(defaultInverse);
00170 
00171    // based on parameter type build a strategy
00172    invOpsStrategy_ = rcp(new InvFactoryDiagStrategy(inverses,defaultInverse));
00173 }
00174 
00175 } // end namspace Teko
 All Classes Files Functions Variables