This file is indexed.

/usr/include/trilinos/MueLu_Zoltan2Interface_decl.hpp is in libtrilinos-muelu-dev 12.10.1-3.

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
// @HEADER
//
// ***********************************************************************
//
//        MueLu: A package for multigrid based preconditioning
//                  Copyright 2012 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
//                    Jonathan Hu       (jhu@sandia.gov)
//                    Andrey Prokopenko (aprokop@sandia.gov)
//                    Ray Tuminaro      (rstumin@sandia.gov)
//
// ***********************************************************************
//
// @HEADER
#ifndef MUELU_ZOLTAN2INTERFACE_DECL_HPP
#define MUELU_ZOLTAN2INTERFACE_DECL_HPP

#include "MueLu_ConfigDefs.hpp"

#if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)

#include <Xpetra_Matrix.hpp>
#include <Xpetra_VectorFactory.hpp>

#include "MueLu_SingleLevelFactoryBase.hpp"
#include "MueLu_Zoltan2Interface_fwd.hpp"

#include "MueLu_Level_fwd.hpp"
#include "MueLu_FactoryBase_fwd.hpp"
#include "MueLu_Utilities_fwd.hpp"

#if defined(HAVE_MUELU_ZOLTAN)
#include "MueLu_ZoltanInterface.hpp"
#endif

namespace MueLu {


  /*!
    @class Zoltan2Interface
    @brief Interface to Zoltan2 library.
    @ingroup Rebalancing

    This interface provides access to partitioning methods in Zoltan2.
    Currently, it supports the RCB algorithm only.

    ## Input/output of Zoltan2Interface ##

    ### User parameters of Zoltan2Interface ###
    Parameter | type | default | master.xml | validated | requested | description
    ----------|------|---------|:----------:|:---------:|:---------:|------------
    | A                                      | Factory | null  |   | * | * | Generating factory of the matrix A used during the prolongator smoothing process |
    | Coordinates                            | Factory | null  |   | * | * | Factory generating coordinates vector used for rebalancing (RCB algorithm)
    | rowWeight                              | int     | 0     |   | * |   | Default weight to rows (total weight = nnz + rowWeight")
    | ParameterList                          | ParamterList | null |  | * |  | Zoltan2 parameters
    | number of partitions                   | GO      | - |  |  |  | Short-cut parameter set by RepartitionFactory. Avoid repartitioning algorithms if only one partition is necessary (see details below)

    The * in the @c master.xml column denotes that the parameter is defined in the @c master.xml file.<br>
    The * in the @c validated column means that the parameter is declared in the list of valid input parameters (see Zoltan2Interface::GetValidParameters).<br>
    The * in the @c requested column states that the data is requested as input with all dependencies (see Zoltan2Interface::DeclareInput).

    ### Variables provided by Zoltan2Interface ###

    After Zoltan2Interface::Build the following data is available (if requested)

    Parameter | generated by | description
    ----------|--------------|------------
    | Partition       | Zoltan2Interface   | GOVector based on the Row map of A (DOF-based) containing the process id the DOF should be living in after rebalancing/repartitioning

    The "Partition" vector is used as input for the RepartitionFactory class.
    If Re-partitioning/rebalancing is necessary it uses the "Partition" variable to create the corresponding Xpetra::Import object which then is used
    by the RebalanceFactory classes (e.g., RebalanceAcFactory, RebalanceTransferFactory,...) to rebalance the coarse level operators.

    The RepartitionHeuristicFactory calculates how many partitions are to be built when performing rebalancing.
    It stores the result in the "number of partitions" variable on the current level (type = GO).
    If it is "number of partitions=1" we skip the Zoltan2 call and just create an dummy "Partition" vector containing zeros only.
    If no repartitioning is necessary (i.e., just keep the current partitioning) we return "Partition = Teuchos::null".
    If "number of partitions" > 1, the algorithm tries to find the requested number of partitions.

  */

  //FIXME: this class should not be templated
  template <class Scalar,
            class LocalOrdinal = typename Xpetra::Matrix<Scalar>::local_ordinal_type,
            class GlobalOrdinal = typename Xpetra::Matrix<Scalar, LocalOrdinal>::global_ordinal_type,
            class Node = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
  class Zoltan2Interface : public SingleLevelFactoryBase {
#undef MUELU_ZOLTAN2INTERFACE_SHORT
#include "MueLu_UseShortNames.hpp"

  public:

    //! @name Constructors/Destructors
    //@{

    //! Constructor
    Zoltan2Interface();

    //! Destructor
    virtual ~Zoltan2Interface() { }
    //@}

    RCP<const ParameterList> GetValidParameterList() const;

    //! @name Input
    //@{
    void DeclareInput(Level& currentLevel) const;
    //@}

    //! @name Build methods.
    //@{
    void Build(Level& currentLevel) const;

    //@}

  private:
    RCP<ParameterList> defaultZoltan2Params;

  };  //class Zoltan2Interface

#ifdef HAVE_MUELU_EPETRA

#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
    (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))

#if defined(HAVE_MUELU_ZOLTAN)
  // Stub partial specialization of Zoltan2Interface for EpetraNode
  // Use ZoltanInterface instead of Zoltan2Interface
  template<>
  class Zoltan2Interface<double,int,int,Xpetra::EpetraNode> : public SingleLevelFactoryBase {
    typedef double              Scalar;
    typedef int                 LocalOrdinal;
    typedef int                 GlobalOrdinal;
    typedef Xpetra::EpetraNode  Node;
#undef MUELU_ZOLTAN2INTERFACE_SHORT
#include "MueLu_UseShortNames.hpp"
  public:
    Zoltan2Interface() {
      level_           = rcp(new Level());
      zoltanInterface_ = rcp(new ZoltanInterface());

      level_->SetLevelID(1);
    }
    virtual ~Zoltan2Interface() {
      zoltanInterface_ = Teuchos::null;
      level_           = Teuchos::null;
    }
    RCP<const ParameterList> GetValidParameterList() const {
      RCP<ParameterList> validParamList = rcp(new ParameterList());

      validParamList->set< RCP<const FactoryBase> >   ("A",             Teuchos::null, "Factory of the matrix A");
      validParamList->set< RCP<const FactoryBase> >   ("Coordinates",   Teuchos::null, "Factory of the coordinates");
      validParamList->set< RCP<const FactoryBase> >   ("number of partitions", Teuchos::null, "Instance of RepartitionHeuristicFactory.");
      validParamList->set< RCP<const ParameterList> > ("ParameterList", Teuchos::null, "Zoltan2 parameters");

      return validParamList;
    };
    void DeclareInput(Level& currentLevel) const {
      Input(currentLevel, "A");
      Input(currentLevel, "number of partitions");
      Input(currentLevel, "Coordinates");
    };
    void Build(Level& currentLevel) const {
      this->GetOStream(Warnings0) << "Tpetra does not support <double,int,int,EpetraNode> instantiation, "
          "switching Zoltan2Interface to ZoltanInterface" << std::endl;

      typedef Xpetra::MultiVector<double, LocalOrdinal, GlobalOrdinal, Node> dMultiVector;

      // Put the data into a fake level
      level_->Set("A",                    Get<RCP<Matrix> >      (currentLevel, "A"));
      level_->Set("Coordinates",          Get<RCP<dMultiVector> >(currentLevel, "Coordinates"));
      level_->Set("number of partitions", currentLevel.Get<GO>("number of partitions"));

      level_->Request("Partition", zoltanInterface_.get());
      zoltanInterface_->Build(*level_);

      RCP<Xpetra::Vector<GO,LO,GO,NO> > decomposition;
      level_->Get("Partition", decomposition, zoltanInterface_.get());
      Set(currentLevel, "Partition", decomposition);
    };

  private:
    RCP<Level>           level_;            // fake Level
    RCP<ZoltanInterface> zoltanInterface_;
  };
#else
  // Stub partial specialization of Zoltan2Interface for EpetraNode
  template<>
  class Zoltan2Interface<double,int,int,Xpetra::EpetraNode> : public SingleLevelFactoryBase {
  public:
    Zoltan2Interface() { throw Exceptions::RuntimeError("Tpetra does not support <double,int,int,EpetraNode> instantiation"); }
    virtual ~Zoltan2Interface() { }
    RCP<const ParameterList> GetValidParameterList() const { return Teuchos::null; };
    void DeclareInput(Level& level) const {};
    void Build(Level &level) const {};
  };
#endif // HAVE_MUELU_ZOLTAN

#endif

#endif // HAVE_MUELU_EPETRA

} //namespace MueLu

#define MUELU_ZOLTAN2INTERFACE_SHORT
#endif //if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)

#endif // MUELU_ZOLTAN2INTERFACE_DECL_HPP