Zoltan2
Zoltan2_EvaluatePartitionFactory.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
50 #ifndef ZOLTAN2_EVALUATE_PARTITION_FACTORY_HPP
51 #define ZOLTAN2_EVALUATE_PARTITION_FACTORY_HPP
55  public:
56 
64  static EvaluatePartition<basic_id_t>*newEvaluatePartition
65  (partitioning_problem_t *problem, const std::string &adapter_name,
66  base_adapter_t *input, ParameterList *params)
67  {
68  RCP<const Teuchos::Comm<int> > CommT = problem->getComm();
69 
70  if (adapter_name == "BasicIdentifier") {
71  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
73  (reinterpret_cast<const basic_id_t *>(input),
74  params, CommT, reinterpret_cast
76  (&problem->getSolution())));
77  } else if (adapter_name == "XpetraMultiVector") {
78  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
80  (reinterpret_cast<const xpetra_mv_adapter *>(input),
81  params, CommT, reinterpret_cast
83  (&problem->getSolution())));
84  } else if (adapter_name == "XpetraCrsGraph") {
85  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
87  (reinterpret_cast<const xcrsGraph_adapter *>(input),
88  params, CommT, reinterpret_cast
90  (&problem->getSolution())));
91  } else if (adapter_name == "XpetraCrsMatrix") {
92  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
94  (reinterpret_cast<const xcrsMatrix_adapter *>(input),
95  params, CommT, reinterpret_cast
97  (&problem->getSolution())));
98  } else if (adapter_name == "BasicVector") {
99  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
101  (reinterpret_cast<const basic_vector_adapter *>(input),
102  params, CommT, reinterpret_cast
104  (&problem->getSolution())));
105  } else if (adapter_name == "PamgenMesh") {
106  return reinterpret_cast<Zoltan2::EvaluatePartition<basic_id_t> *>
108  (reinterpret_cast<const pamgen_adapter_t *>(input),
109  params, CommT, reinterpret_cast
111  (&problem->getSolution())));
112  }
113  return nullptr; // adapter type not known
114  }
115  };
116 }
117 #endif // ZOLTAN2_EVALUATE_PARTITION_FACTORY_HPP
static EvaluatePartition< basic_id_t > * newEvaluatePartition(partitioning_problem_t *problem, const std::string &adapter_name, base_adapter_t *input, ParameterList *params)
Zoltan2::EvaluatePartition factory method.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Provides access for Zoltan2 to Xpetra::CrsGraph data.
This class represents a collection of global Identifiers and their associated weights, if any.
A PartitioningSolution is a solution to a partitioning problem.
RCP< const Comm< int > > getComm()
Return the communicator used by the problem.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
An adapter for Xpetra::MultiVector.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
BaseAdapter defines methods required by all Adapters.
A class that computes and returns quality metrics.
brief EvaluatePartitionFActory class contains 1 static factory method