This file is indexed.

/usr/include/trilinos/Teko_BlockInvDiagonalStrategy.hpp is in libtrilinos-teko-dev 12.12.1-5.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
/*
// @HEADER
// 
// ***********************************************************************
// 
//      Teko: A package for block and physics based preconditioning
//                  Copyright 2010 Sandia Corporation 
//  
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//  
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//  
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//  
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//  
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission. 
//  
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//  
// Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
// 
// ***********************************************************************
// 
// @HEADER

*/

#ifndef __Teko_BlockInvDiagonalStrategy_hpp__
#define __Teko_BlockInvDiagonalStrategy_hpp__

#include <vector>

// Teuchos includes
#include "Teuchos_RCP.hpp"

// Thyra includes
#include "Thyra_LinearOpBase.hpp"

// Teko includes
#include "Teko_Utilities.hpp"
#include "Teko_InverseFactory.hpp"
#include "Teko_BlockPreconditionerFactory.hpp"

namespace Teko {

/** this should be paired with a BlockJacobiPreconditionerFactory
  * or BlockGSPreconditionerFactory.  The idea is that this object
  * provides an approximate inverse operator for each of the diagonal
  * blocks.  Then, the [whatever]PreconditionerFactory can easily
  * construct an approximate inverse. The system under consideration
  * is
  * 
  *    A = [ D0  U01 U02 ...]
  *        [ L10  D1 U12 ...]
  *        [ L20 L21  D2 ...]
  *        [  .   .   .  ...]
  *
  * where inverses of D0,D1,D2...DN are needed.
  */
class BlockInvDiagonalStrategy {
public:
   virtual ~BlockInvDiagonalStrategy() {}

   //! returns an (approximate) inverse of the diagonal blocks of A
   virtual void getInvD(const BlockedLinearOp & A,BlockPreconditionerState & state,
                        std::vector<LinearOp> & invDiag) const = 0;

};

/** This is a simple strategy for a [whatever]PreconditionerFactory
  * it simply returns statically set RCP pointers to the passed in
  * inv(D0) and inv(D1) operators. Not this will _not_ permit 
  * efficient implementations when the preconditioner has to be rebuilt
  * or reused often.
  */
class StaticInvDiagStrategy : public BlockInvDiagonalStrategy {
public:
   StaticInvDiagStrategy(const LinearOp & invD0,
                          const LinearOp & invD1)
   { invDiag_.push_back(invD0); invDiag_.push_back(invD1); }

   StaticInvDiagStrategy(const std::vector<LinearOp> & invD)
      : invDiag_(invD)
   { }

   virtual ~StaticInvDiagStrategy() {}

   /** returns an (approximate) inverse of the diagonal blocks of A
     * where A is closely related to the original source for invD0 and invD1
     */
   virtual void getInvD(const BlockedLinearOp & A, BlockPreconditionerState & state,
                        std::vector<LinearOp> & invDiag) const
   { invDiag.clear(); invDiag = invDiag_; }

protected:
   // stored inverse operators
   std::vector<Teuchos::RCP<const Thyra::LinearOpBase<double> > > invDiag_;
};

/** A simple class that takes a vector of the <code>InverseFactory</code> objects
  * and pairs each with the diagonal element of the block matrix. This provides
  * the operators needed to use Gauss-Seidel or Jacobi.
  */
class InvFactoryDiagStrategy : public BlockInvDiagonalStrategy {
public:
   /** Constructor accepting a single inverse factory that will be used
     * to invert all diagonal blocks.
     *
     * \param[in] factory Factory to be used to invert each diagonal block.
     */
   InvFactoryDiagStrategy(const Teuchos::RCP<InverseFactory> & factory);

   /** Constructor that lets the inverse of each block be set individually.
     *
     * \param[in] factories A vector of <code>InverseFactory</code> objects
     *                      which should be the same length as the number of
     *                      diagonal blocks.
     * \param[in] defaultFact The default factory to use if none is specified.
     */
   InvFactoryDiagStrategy(const std::vector<Teuchos::RCP<InverseFactory> > & factories,
                          const Teuchos::RCP<InverseFactory> & defaultFact=Teuchos::null);

   virtual ~InvFactoryDiagStrategy() {}

   /** returns an (approximate) inverse of the diagonal blocks of A
     * where A is closely related to the original source for invD0 and invD1
     * 
     * \param[in]     A       Operator to extract the block diagonals from.
     * \param[in]     state   State object for this operator.
     * \param[in,out] invDiag Vector eventually containing the inverse operators
     */
   virtual void getInvD(const BlockedLinearOp & A, BlockPreconditionerState & state,
                        std::vector<LinearOp> & invDiag) const;

   //! Get factories for testing purposes.
   const std::vector<Teuchos::RCP<InverseFactory> > & getFactories() const
   { return invDiagFact_; }

protected:
   // stored inverse operators
   std::vector<Teuchos::RCP<InverseFactory> > invDiagFact_;
   Teuchos::RCP<InverseFactory> defaultInvFact_;

   //! Conveinence function for building inverse operators
   LinearOp buildInverse(const InverseFactory & invFact,const LinearOp & matrix,
                         BlockPreconditionerState & state,
                         const std::string & opPrefix,int i) const;

private:
   InvFactoryDiagStrategy();
   InvFactoryDiagStrategy(const InvFactoryDiagStrategy &);
};

} // end namespace Teko

#endif