This file is indexed.

/usr/include/trilinos/Ifpack_Partitioner.h is in libtrilinos-ifpack-dev 12.4.2-2.

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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
/*
//@HEADER
// ***********************************************************************
//
//       Ifpack: Object-Oriented Algebraic Preconditioner Package
//                 Copyright (2002) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// 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 Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
//@HEADER
*/

#ifndef IFPACK_PARTITIONER_H
#define IFPACK_PARTITIONER_H

#include "Ifpack_ConfigDefs.h"
#include "Teuchos_ParameterList.hpp"
class Epetra_Comm;
class Ifpack_Graph;
class Epetra_Map;
class Epetra_BlockMap;
class Epetra_Import;

//! Ifpack_Partitioner: A class to decompose local Ifpack_Graph's.

/*!

  Class Ifpack_Partitioner enables the decomposition of a local
  Ifpack_Graph's. It is supposed that the graph refers to
  a localized matrix (that is, a matrix that has been filtered
  through Ifpack_LocalFilter).

  The overloaded operator (int i) can be used to extract the local partition
  ID of local row i.

  The partitions created by Ifpack_Partitioner derived clased
  are non-overlapping in graph sense. This means that each row
  (or, more approriately, vertex)
  of \c G is assigned to exactly one partition.

  Partitioner can be extended using the functionalities of class
  Ifpack_OverlappingPartitioner (itself derived from Ifpack_Partitioner.
  This class extends the non-overlapping partitions by the required
  amount of overlap, considering local nodes only (that is, this
  overlap do \e not modify the overlap among the processes).

  Ifpack_Partitioner is a pure virtual class. Concrete implementations
  are:
  - Ifpack_LinearPartitioner, which allows the decomposition of the
    rows of the graph in simple consecutive chunks;
  - Ifpack_METISPartitioner, which calls METIS to decompose the graph
    (this requires the configuration option --enable-ifpack-metis);
  - Ifpack_GreedyPartitioner, a simple greedy algorith;
  - Ifpack_EquationPartitioner, which creates \c NumPDEEqns parts
    (where \c NumPDEEqns is the number of equations in the linear
    system). It is supposed that all the equations referring to the
    same grid node are ordered consecutively. Besides, the
    number of equations per node must be constant in the domain.

  Generically, a constructor requires an Ifpack_Graph object.
  Ifpack_Graph is a pure virtual class. Concrete implentations are:
  - Ifpack_Graph_Epetra_CrsGraph, a light-weight class to wrap
    Epetra_CrsGraph objects as Ifpack_Graph objects;
  - Ifpack_Graph_Epetra_RowMatrix, a light-weight class to
    wrap Epetra_RowMatrix objects as Ifpack_Graph objects.

  <P>An example of use is an Ifpack_Partitioner derived class is as follows:
  \code
#include "Ifpack_Partitioner.h"
#include "Ifpack_LinearPartitioner.h"
#include "Ifpack_Graph.h"
#include "Ifpack_Graph_Epetra_CrsGraph.h"
...
Epetra_CrsMatrix* A;         // A is filled
// create the wrapper from Epetra_CrsGraph
Ifpack_Graph* Graph = new Ifpack_Graph_Epetra_CrsGraph(A);

// we aim to create non-overlapping partitions only
Ifpack_Partitioner Partitioner(Graph);

Ifpack_Partitioner* Partitioner;
Partitioner = new Ifpack_Graph_Epetra_CrsGraph(&A);

// we want 16 local parts
List.set("partitioner: local parts", 16);
// and an overlap of 0 among the local parts (default option)
List.set("partitioner: overlap", 0);

// decompose the graph
Partitioner.Create(List);

// now Graph can be deleted, as Partitioner contains all the
// necessary information to use the partitions
delete Graph;

// we can get the number of parts actually created...
int NumParts = Partitioner.NumParts();

// ... and the number of rows in each of them
for (int i = 0 ; i < NumParts ; ++i) {
  cout << "rows in " << i << "=" << Partitioner.RowsInPart(i);
}

// .. and, for non-overlapping partitions only, the partition ID
// for each local row simply using:
for (int i = 0 ; i < A->NumMyRows() ; ++i)
  cout << "Partition[" << i <<"] = " << Partitioner(i) << endl;

\endcode

When overlapping partitiones are created, the user can get the
row ID contained in each partition as follows:
\code
for (int i = 0 ; i < NumParts ; ++i) {
  for (int j = 0 ; j < Partitioner.RowsInPart(i) ; ++j) {
    cout << "Partition " << i << ", contains local row "
         << Partitioner(i,j) << endl;
  }
}
\endcode

Ifpack_Partitioner is used to create the subblocks in
Ifpack_BlockJacobi, Ifpack_BlockGaussSeidel, and
Ifpack_BlockSymGaussSeidel.

\author Marzio Sala, SNL 9214.

\date Last modified on Nov-04.

*/
class Ifpack_Partitioner {

public:

  //! Destructor.
  virtual ~Ifpack_Partitioner() {};

  //! Returns the number of computed local partitions.
  virtual int NumLocalParts() const = 0;

  //! Returns the overlapping level.
  virtual int OverlappingLevel() const = 0;

  //! Returns the local non-overlapping partition ID of the specified row.
  /*! Returns the non-overlapping partition ID of the specified row.
   \param
   MyRow - (In) local row numbe

   \return
   Local ID of non-overlapping partition for \t MyRow.
   */
  virtual int operator() (int MyRow) const = 0;

  //! Returns the local overlapping partition ID of the j-th node in partition i.
  virtual int operator() (int i, int j) const = 0;

  //! Returns the number of rows contained in specified partition.
  virtual int NumRowsInPart(const int Part) const = 0;

  //! Copies into List the rows in the (overlapping) partition Part.
  virtual int RowsInPart(const int Part, int* List) const = 0;

  //! Returns a pointer to the integer vector containing the non-overlapping partition ID of each local row.
  virtual const int* NonOverlappingPartition() const = 0;

  //! Sets all the parameters for the partitioner.
  virtual int SetParameters(Teuchos::ParameterList& List) = 0;

  //! Computes the partitions. Returns 0 if successful.
  virtual int Compute() = 0;

  //! Returns true if partitions have been computed successfully.
  virtual bool IsComputed() = 0;

  //! Prints basic information about the partitioning object.
  virtual std::ostream& Print(std::ostream& os) const = 0;

}; // class Ifpack_Partitioner

inline std::ostream& operator<<(std::ostream& os, const Ifpack_Partitioner& obj)
{
  return(obj.Print(os));
}

#endif // IFPACK_PARTITIONER_H