This file is indexed.

/usr/include/trilinos/Ifpack2_IlukGraph.hpp is in libtrilinos-ifpack2-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
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
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
/*@HEADER
// ***********************************************************************
//
//       Ifpack2: Tempated Object-Oriented Algebraic Preconditioner Package
//                 Copyright (2009) 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
 */

/// \file Ifpack2_IlukGraph.hpp
/// \brief Declaration and definition of IlukGraph
///
/// \warning This header file is an implementation detail of RILUK.
///   Its contents may change or it may go away at any time.

#ifndef IFPACK2_ILUK_GRAPH_HPP
#define IFPACK2_ILUK_GRAPH_HPP

#include <algorithm>
#include <vector>

#include <Ifpack2_ConfigDefs.hpp>
#include <Teuchos_ParameterList.hpp>
#include <Teuchos_CommHelpers.hpp>
#include <Tpetra_CrsGraph.hpp>
#include <Tpetra_Import.hpp>
#include <Ifpack2_CreateOverlapGraph.hpp>
#include <Ifpack2_Parameters.hpp>

namespace Ifpack2 {

/// \class IlukGraph
/// \tparam GraphType A specialization of Tpetra::CrsGraph (tested) or
///   Tpetra::RowGraph (not tested).
/// \brief Construct a level filled graph for use in computing an
///   ILU(k) incomplete factorization.
///
/// \warning This class is an implementation detail of RILUK.  Its
///   interface may change or it may go away at any time.
///
/// Ifpack2::IlukGraph enables the construction of matrix graphs using
/// level-fill algorithms.  The only function required for
/// construction is a getGlobalRowView capability, i.e., the graph
/// that is passed in to the constructor must implement the
/// Tpetra::RowGraph interface defined in Tpetra_RowGraph.hpp.
///
/// \section Ifpack2_IlukGraph_cnstr Constructing IlukGraph objects
///
/// Constructing an IlukGraph is a two step process.  First, call the
/// constructor, passing in a CrsGraph object and an integer
/// indicating the desired level of fill.  Then, call the initialize()
/// method to complete the process.  This naturally matches the
/// three-stage initialization of Preconditioner objects, and in
/// particular of RILUK.
///
/// It is worth noting that an IlukGraph object creates two
/// Tpetra::CrsGraph objects containing L and U, the graphs for the
/// lower and upper triangular parts of the ILU(k) graph.  Thus, it is
/// possible to manually insert and delete graph entries in L and U
/// via the Tpetra::CrsGraph InsertIndices and RemoveIndices
/// functions.  However, in this case FillComplete must be called
/// before the graph is used for subsequent operations.
template<class GraphType>
class IlukGraph : public Teuchos::Describable {
public:
  typedef typename GraphType::local_ordinal_type local_ordinal_type;
  typedef typename GraphType::global_ordinal_type global_ordinal_type;
  typedef typename GraphType::node_type node_type;

  //! Tpetra::RowGraph specialization used by this class.
  typedef Tpetra::RowGraph<local_ordinal_type,
                           global_ordinal_type,
                           node_type> row_graph_type;
  //! Tpetra::CrsGraph specialization used by this class.
  typedef Tpetra::CrsGraph<local_ordinal_type,
                           global_ordinal_type,
                           node_type> crs_graph_type;

  /// \brief Constructor.
  ///
  /// Create a IlukGraph object using the input graph and specified
  /// level of fill.
  ///
  /// \param G [in] An existing graph.
  /// \param levelFill [in] The level of fill to compute; the k of ILU(k).
  /// \param levelOverlap [in] The level of overlap between subdomains.
  ///
  /// \note Actual construction occurs in initialize().
  IlukGraph (const Teuchos::RCP<const GraphType>& G,
             const int levelFill,
             const int levelOverlap);

  //! IlukGraph Destructor
  virtual ~IlukGraph ();

  /// \brief Set parameters.
  ///
  /// This method recognizes two parameter names: Level_fill and
  /// Level_overlap.  Both are case insensitive, and in both cases the
  /// parameter must have type int.
  void setParameters (const Teuchos::ParameterList& parameterlist);

  /// \brief Set up the graph structure of the L and U factors.
  ///
  /// This method is called "initialize" by analogy with
  /// Preconditioner, where initialize() computes the symbolic
  /// (incomplete) factorization, and compute() computes the
  /// corresponding numeric factorization.  IlukGraph is just a graph,
  /// so it can only compute a symbolic factorization (i.e., the graph
  /// structure of the factorization).  Hence, it implements
  /// initialize(), but not compute().  RILUK calls IlukGraph's
  /// initialize() method in its own initialize() method, as one would
  /// expect.
  void initialize ();

  //! The level of fill used to construct this graph.
  int getLevelFill () const { return LevelFill_; }

  //! The level of overlap used to construct this graph.
  int getLevelOverlap () const { return LevelOverlap_; }

  //! Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph.
  Teuchos::RCP<crs_graph_type> getL_Graph () const {
    return L_Graph_;
  }

  //! Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph.
  Teuchos::RCP<crs_graph_type> getU_Graph () const {
    return U_Graph_;
  }

  //! Returns the the overlapped graph.
  Teuchos::RCP<const crs_graph_type> getOverlapGraph () const {
    return OverlapGraph_;
  }

  //! Returns the global number of diagonals in the ILU(k) graph.
  size_t getNumGlobalDiagonals() const { return NumGlobalDiagonals_; }

private:
  typedef typename GraphType::map_type map_type;

  /// \brief Copy constructor (UNIMPLEMENTED; DO NOT USE).
  ///
  /// This copy constructor is declared private and unimplemented, in
  /// order to forbid its use syntactically.  If you decide that you
  /// need to implement this method, it should do deep copies of all
  /// internal graphs.  It may do a shallow copy of the input graph
  /// (which is not modified).  Also, you must implement operator= in
  /// a similar way.
  IlukGraph (const IlukGraph<GraphType>&);

  /// \brief Assignment operator (UNIMPLEMENTED; DO NOT USE).
  ///
  /// This assignment operator is declared private and unimplemented,
  /// in order to forbid its use syntactically.  If you decide that
  /// you need to implement this method, it should do deep copies of
  /// all internal graphs.  It may do a shallow copy of the input
  /// graph (which is not modified).  Also, you must implement the
  /// copy constructor in a similar way.
  IlukGraph& operator= (const IlukGraph<GraphType>&);

  //! Construct the overlap graph.
  void constructOverlapGraph();

  Teuchos::RCP<const GraphType> Graph_;
  Teuchos::RCP<const crs_graph_type> OverlapGraph_;
  int LevelFill_;
  int LevelOverlap_;
  Teuchos::RCP<crs_graph_type> L_Graph_;
  Teuchos::RCP<crs_graph_type> U_Graph_;
  size_t NumMyDiagonals_;
  size_t NumGlobalDiagonals_;
};


template<class GraphType>
IlukGraph<GraphType>::
IlukGraph (const Teuchos::RCP<const GraphType>& G,
           const int levelFill,
           const int levelOverlap)
  : Graph_ (G),
    LevelFill_ (levelFill),
    LevelOverlap_ (levelOverlap),
    NumMyDiagonals_ (0),
    NumGlobalDiagonals_ (0)
{}


template<class GraphType>
IlukGraph<GraphType>::~IlukGraph()
{}


template<class GraphType>
void IlukGraph<GraphType>::
setParameters (const Teuchos::ParameterList& parameterlist)
{
  getParameter (parameterlist, "iluk level-of-fill", LevelFill_);
  getParameter (parameterlist, "iluk level-of-overlap", LevelOverlap_);
}


template<class GraphType>
void IlukGraph<GraphType>::constructOverlapGraph () {
  // FIXME (mfh 22 Dec 2013) This won't do if we want
  // RILUK::initialize() to do the right thing (that is,
  // unconditionally recompute the "symbolic factorization").
  if (OverlapGraph_ == Teuchos::null) {
    OverlapGraph_ = createOverlapGraph<GraphType> (Graph_, LevelOverlap_);
  }
}


template<class GraphType>
void IlukGraph<GraphType>::initialize()
{
  using Teuchos::Array;
  using Teuchos::ArrayView;
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Teuchos::REDUCE_SUM;
  using Teuchos::reduceAll;

  size_t NumIn, NumL, NumU;
  bool DiagFound;

  constructOverlapGraph();

  L_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
                                      OverlapGraph_->getRowMap (), 0));
  U_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
                                      OverlapGraph_->getRowMap (), 0));

  // Get Maximum Row length
  const int MaxNumIndices = OverlapGraph_->getNodeMaxNumRowEntries ();

  Array<local_ordinal_type> L (MaxNumIndices);
  Array<local_ordinal_type> U (MaxNumIndices);

  // First we copy the user's graph into L and U, regardless of fill level

  // FIXME (mfh 23 Dec 2013) Use size_t or whatever
  // getNodeNumElements() returns, instead of ptrdiff_t.
  const int NumMyRows = OverlapGraph_->getRowMap ()->getNodeNumElements ();
  NumMyDiagonals_ = 0;

  for (int i = 0; i< NumMyRows; ++i) {
    ArrayView<const local_ordinal_type> my_indices;
    OverlapGraph_->getLocalRowView (i, my_indices);

    // Split into L and U (we don't assume that indices are ordered).

    NumL = 0;
    NumU = 0;
    DiagFound = false;
    NumIn = my_indices.size();

    for (size_t j = 0; j < NumIn; ++j) {
      const local_ordinal_type k = my_indices[j];

      if (k<NumMyRows) { // Ignore column elements that are not in the square matrix

        if (k==i) {
          DiagFound = true;
        }
        else if (k < i) {
          L[NumL] = k;
          NumL++;
        }
        else {
          U[NumU] = k;
          NumU++;
        }
      }
    }

    // Check in things for this row of L and U

    if (DiagFound) {
      ++NumMyDiagonals_;
    }
    if (NumL) {
      ArrayView<const local_ordinal_type> L_view = L.view (0, NumL);
      L_Graph_->insertLocalIndices (i, L_view);
    }
    if (NumU) {
      ArrayView<const local_ordinal_type> U_view = U.view (0, NumU);
      U_Graph_->insertLocalIndices (i, U_view);
    }
  }

  if (LevelFill_ > 0) {
    // Complete Fill steps
    RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
    RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
    RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
    RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
    RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
    params->set ("Optimize Storage",false);
    L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
    U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
    L_Graph_->resumeFill ();
    U_Graph_->resumeFill ();

    // At this point L_Graph and U_Graph are filled with the pattern of input graph,
    // sorted and have redundant indices (if any) removed.  Indices are zero based.
    // LevelFill is greater than zero, so continue...

    int MaxRC = NumMyRows;
    std::vector<std::vector<int> > Levels(MaxRC);
    std::vector<int> LinkList(MaxRC);
    std::vector<int> CurrentLevel(MaxRC);
    Array<local_ordinal_type> CurrentRow (MaxRC + 1);
    std::vector<int> LevelsRowU(MaxRC);

    for (int i = 0; i < NumMyRows; ++i) {
      int First, Next;

      // copy column indices of row into workspace and sort them

      size_t LenL = L_Graph_->getNumEntriesInLocalRow(i);
      size_t LenU = U_Graph_->getNumEntriesInLocalRow(i);
      size_t Len = LenL + LenU + 1;

      CurrentRow.resize(Len);

      L_Graph_->getLocalRowCopy(i, CurrentRow(), LenL);      // Get L Indices
      CurrentRow[LenL] = i;                                     // Put in Diagonal
      if (LenU > 0) {
        ArrayView<local_ordinal_type> URowView = CurrentRow.view (LenL+1, LenU);
        // Get U Indices
        U_Graph_->getLocalRowCopy (i, URowView, LenU);
      }

      // Construct linked list for current row

      for (size_t j=0; j<Len-1; j++) {
        LinkList[CurrentRow[j]] = CurrentRow[j+1];
        CurrentLevel[CurrentRow[j]] = 0;
      }

      LinkList[CurrentRow[Len-1]] = NumMyRows;
      CurrentLevel[CurrentRow[Len-1]] = 0;

      // Merge List with rows in U

      First = CurrentRow[0];
      Next = First;
      while (Next < i) {
        int PrevInList = Next;
        int NextInList = LinkList[Next];
        int RowU = Next;
        // Get Indices for this row of U
        ArrayView<const local_ordinal_type> IndicesU;
        U_Graph_->getLocalRowView (RowU, IndicesU);
        // FIXME (mfh 23 Dec 2013) size() returns ptrdiff_t, not int.
        int LengthRowU = IndicesU.size ();

        int ii;

        // Scan RowU

        for (ii = 0; ii < LengthRowU; /*nop*/) {
          int CurInList = IndicesU[ii];
          if (CurInList < NextInList) {
            // new fill-in
            int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
            if (NewLevel <= LevelFill_) {
              LinkList[PrevInList]  = CurInList;
              LinkList[CurInList] = NextInList;
              PrevInList = CurInList;
              CurrentLevel[CurInList] = NewLevel;
            }
            ii++;
          }
          else if (CurInList == NextInList) {
            PrevInList = NextInList;
            NextInList = LinkList[PrevInList];
            int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
            CurrentLevel[CurInList] = std::min (CurrentLevel[CurInList], NewLevel);
            ii++;
          }
          else { // (CurInList > NextInList)
            PrevInList = NextInList;
            NextInList = LinkList[PrevInList];
          }
        }
        Next = LinkList[Next];
      }

      // Put pattern into L and U

      CurrentRow.resize (0);

      Next = First;

      // Lower

      while (Next < i) {
        CurrentRow.push_back (Next);
        Next = LinkList[Next];
      }

      // FIXME (mfh 23 Dec 2013) It's not clear to me that
      // removeLocalIndices works like people expect it to work.  In
      // particular, it does not actually change the column Map.
      L_Graph_->removeLocalIndices (i); // Delete current set of Indices
      if (CurrentRow.size() > 0) {
        L_Graph_->insertLocalIndices (i, CurrentRow ());
      }

      // Diagonal

      TEUCHOS_TEST_FOR_EXCEPTION(
        Next != i, std::runtime_error,
        "Ifpack2::IlukGraph::initialize: FATAL: U has zero diagonal")

      LevelsRowU[0] = CurrentLevel[Next];
      Next = LinkList[Next];

      // Upper

      CurrentRow.resize (0);
      LenU = 0;

      while (Next < NumMyRows) {
        LevelsRowU[LenU+1] = CurrentLevel[Next];
        CurrentRow.push_back (Next);
        ++LenU;
        Next = LinkList[Next];
      }

      // FIXME (mfh 23 Dec 2013) It's not clear to me that
      // removeLocalIndices works like people expect it to work.  In
      // particular, it does not actually change the column Map.

      U_Graph_->removeLocalIndices (i); // Delete current set of Indices
      if (LenU > 0) {
        U_Graph_->insertLocalIndices (i, CurrentRow ());
      }

      // Allocate and fill Level info for this row
      Levels[i] = std::vector<int> (LenU+1);
      for (size_t jj=0; jj<LenU+1; jj++) {
        Levels[i][jj] = LevelsRowU[jj];
      }
    }
  }

  // Complete Fill steps
  RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
  RCP<const map_type> L_RangeMap  = Graph_->getRangeMap ();
  RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
  RCP<const map_type> U_RangeMap  = OverlapGraph_->getRowMap ();
  L_Graph_->fillComplete (L_DomainMap, L_RangeMap);//DoOptimizeStorage is default here...
  U_Graph_->fillComplete (U_DomainMap, U_RangeMap);//DoOptimizeStorage is default here...

  reduceAll<int, size_t> (* (L_DomainMap->getComm ()), REDUCE_SUM, 1,
                          &NumMyDiagonals_, &NumGlobalDiagonals_);
}

} // namespace Ifpack2

#endif /* IFPACK2_ILUK_GRAPH_HPP */