Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_BlockInvDiagonalStrategy.cpp
1 /*
2 // @HEADER
3 //
4 // ***********************************************************************
5 //
6 // Teko: A package for block and physics based preconditioning
7 // Copyright 2010 Sandia Corporation
8 //
9 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
10 // the U.S. Government retains certain rights in this software.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
40 //
41 // ***********************************************************************
42 //
43 // @HEADER
44 
45 */
46 
47 #include "Teko_BlockInvDiagonalStrategy.hpp"
48 
49 namespace Teko {
50 
51 InvFactoryDiagStrategy::InvFactoryDiagStrategy(const Teuchos::RCP<InverseFactory> & factory)
52 {
53  // only one factory to use!
54  invDiagFact_.resize(1,factory);
55  defaultInvFact_ = factory;
56 }
57 
58 InvFactoryDiagStrategy::InvFactoryDiagStrategy(const std::vector<Teuchos::RCP<InverseFactory> > & factories,
59  const Teuchos::RCP<InverseFactory> & defaultFact)
60 {
61  invDiagFact_ = factories;
62 
63  if(defaultFact==Teuchos::null)
64  defaultInvFact_ = invDiagFact_[0];
65  else
66  defaultInvFact_ = defaultFact;
67 }
68 
72 void InvFactoryDiagStrategy::getInvD(const BlockedLinearOp & A,BlockPreconditionerState & state,std::vector<LinearOp> & invDiag) const
73 {
74  Teko_DEBUG_SCOPE("InvFactoryDiagStrategy::getInvD",10);
75 
76  // loop over diagonals, build an inverse operator for each
77  int diagCnt = A->productRange()->numBlocks();
78  int invCnt = invDiagFact_.size();
79 
80  Teko_DEBUG_MSG("# diags = " << diagCnt << ", # inverses = " << invCnt,6);
81 
82  const std::string opPrefix = "JacobiDiagOp";
83  if(diagCnt<=invCnt) {
84  for(int i=0;i<diagCnt;i++)
85  invDiag.push_back(buildInverse(*invDiagFact_[i],getBlock(i,i,A),state,opPrefix,i));
86  }
87  else {
88  for(int i=0;i<invCnt;i++)
89  invDiag.push_back(buildInverse(*invDiagFact_[i],getBlock(i,i,A),state,opPrefix,i));
90 
91  for(int i=invCnt;i<diagCnt;i++)
92  invDiag.push_back(buildInverse(*defaultInvFact_,getBlock(i,i,A),state,opPrefix,i));
93  }
94 }
95 
96 LinearOp InvFactoryDiagStrategy::buildInverse(const InverseFactory & invFact,const LinearOp & matrix,
98  const std::string & opPrefix,int i) const
99 {
100  std::stringstream ss;
101  ss << opPrefix << "_" << i;
102 
103  ModifiableLinearOp & invOp = state.getModifiableOp(ss.str());
104  if(invOp==Teuchos::null)
105  invOp = Teko::buildInverse(invFact,matrix);
106  else
107  rebuildInverse(invFact,matrix,invOp);
108 
109  return invOp;
110 }
111 
112 } // end namespace Teko
void rebuildInverse(const InverseFactory &factory, const LinearOp &A, InverseLinearOp &invA)
Abstract class for building an inverse operator.
virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state, std::vector< LinearOp > &invDiag) const
LinearOp buildInverse(const InverseFactory &invFact, const LinearOp &matrix, BlockPreconditionerState &state, const std::string &opPrefix, int i) const
Conveinence function for building inverse operators.
An implementation of a state object for block preconditioners.
InverseLinearOp buildInverse(const InverseFactory &factory, const LinearOp &A)
Build an inverse operator using a factory and a linear operator.
virtual Teko::ModifiableLinearOp & getModifiableOp(const std::string &name)
Add a named operator to the state object.