50 #ifndef _ZOLTAN2_MATRIXPARTITIONINGPROBLEM_HPP_
51 #define _ZOLTAN2_MATRIXPARTITIONINGPROBLEM_HPP_
104 template<
typename Adapter>
110 typedef typename Adapter::gno_t
gno_t;
111 typedef typename Adapter::lno_t
lno_t;
121 const RCP<
const Teuchos::Comm<int> > &comm):
122 Problem<Adapter>(A,p,comm), solution_(),
132 #ifdef HAVE_ZOLTAN2_MPI
133 typedef Teuchos::OpaqueWrapper<MPI_Comm> mpiWrapper_t;
138 rcp<const Comm<int> >(new Teuchos::MpiComm<int>(
139 Teuchos::opaqueWrapper(mpicomm))))
182 void solve(
bool updateInputData=
true );
190 return *(solution_.getRawPtr());
208 Array<std::string> algorithm_names(1);
209 algorithm_names[0] =
"2D Cartesian";
210 RCP<Teuchos::StringValidator> algorithm_Validator = Teuchos::rcp(
211 new Teuchos::StringValidator( algorithm_names ));
212 pl.set(
"algorithm",
"2D Cartesian",
"partitioning algorithm",
213 algorithm_Validator);
229 pl.set(
"imbalance_tolerance", 1.1,
"imbalance tolerance, ratio of "
233 RCP<Teuchos::EnhancedNumberValidator<int>> num_global_parts_Validator =
234 Teuchos::rcp(
new Teuchos::EnhancedNumberValidator<int>(
235 1, Teuchos::EnhancedNumberTraits<int>::max()) );
236 pl.set(
"num_global_parts", 1,
"global number of parts to compute "
237 "(0 means use the number of processes)", num_global_parts_Validator);
240 RCP<Teuchos::EnhancedNumberValidator<int>> num_local_parts_Validator =
241 Teuchos::rcp(
new Teuchos::EnhancedNumberValidator<int>(
242 0, Teuchos::EnhancedNumberTraits<int>::max()) );
243 pl.set(
"num_local_parts", 0,
"number of parts to compute for this "
244 "process (num_global_parts == sum of all num_local_parts)",
245 num_local_parts_Validator);
269 void initializeProblem();
271 void createPartitioningProblem(
bool newData);
273 RCP<MatrixPartitioningSolution<Adapter> > solution_;
279 std::string algName_;
287 template <
typename Adapter>
288 void MatrixPartitioningProblem<Adapter>::initializeProblem()
292 this->env_->debug(
DETAILED_STATUS,
"MatrixPartitioningProblem::initializeProblem");
294 inputType_ = this->inputAdapter_->adapterType();
299 std::cerr <<
"Error: only matrix adapter type supported" << std::endl;
303 this->env_->memory(
"After initializeProblem");
309 template <
typename Adapter>
312 std::cout <<
"MatrixPartitioningProblem solve " << std::endl;
321 createPartitioningProblem(updateInputData);
345 this->baseInputAdapter_));
354 this->env_->timerStart(
MACRO_TIMERS,
"create solution");
359 this->envConst_, this->comm_,
364 solution_ = rcp(soln);
367 this->env_->memory(
"After creating Solution");
373 this->algorithm_->partitionMatrix(solution_);
383 template <
typename Adapter>
387 "MatrixPartitioningProblem::createPartitioningProblem");
389 using Teuchos::ParameterList;
412 const Teuchos::ParameterEntry *pe;
413 std::string defString(
"default");
426 std::string algorithm(defString);
427 pe = pl.getEntryPtr(
"algorithm");
429 algorithm = pe->getValue<std::string>(&algorithm);
437 if (algorithm != defString)
441 if (algorithm == std::string(
"2D Cartesian"))
443 algName_ = algorithm;
448 throw std::logic_error(
"parameter list algorithm is invalid");
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Time an algorithm (or other entity) as a whole.
MatrixPartitioningProblem(Adapter *A, ParameterList *p, const RCP< const Teuchos::Comm< int > > &comm)
Constructor where Teuchos communicator is specified.
~MatrixPartitioningProblem()
Destructor.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this Problem.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
MatrixPartitioningProblem sets up partitioning problems for the user.
Adapter::base_adapter_t base_adapter_t
MatrixPartitioningProblem(Adapter *A, ParameterList *p)
Constructor where communicator is the Teuchos default.
static RCP< Teuchos::AnyNumberParameterEntryValidator > getAnyDoubleValidator()
Exists to make setting up validators less cluttered.
sub-steps, each method's entry and exit
SparseMatrixAdapter_t::part_t part_t
Teuchos::ParameterList & getParametersNonConst()
Returns a reference to a non-const copy of the parameters.
A PartitioningSolution is a solution to a partitioning problem.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
BaseAdapterType
An enum to identify general types of adapters.
Problem base class from which other classes (PartitioningProblem, ColoringProblem, OrderingProblem, MatchingProblem, etc.) derive.
Defines the Problem base class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Adapter::scalar_t scalar_t
A PartitioningSolution is a solution to a partitioning problem.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t