Zoltan2
Zoltan2_GraphModel.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_GRAPHMODEL_HPP_
51 #define _ZOLTAN2_GRAPHMODEL_HPP_
52 
53 #include <Zoltan2_Model.hpp>
54 #include <Zoltan2_ModelHelpers.hpp>
55 #include <Zoltan2_InputTraits.hpp>
57 #include <Zoltan2_GraphAdapter.hpp>
60 #include <Zoltan2_MeshAdapter.hpp>
61 #include <Zoltan2_StridedData.hpp>
62 #include <unordered_map>
63 
64 namespace Zoltan2 {
65 
67 
79 template <typename Adapter>
80 class GraphModel : public Model<Adapter>
81 {
82 public:
83 
84 #ifndef DOXYGEN_SHOULD_SKIP_THIS
85  typedef typename Adapter::scalar_t scalar_t;
86  typedef typename Adapter::gno_t gno_t;
87  typedef typename Adapter::lno_t lno_t;
88  typedef typename Adapter::node_t node_t;
89  typedef typename Adapter::user_t user_t;
90  typedef typename Adapter::userCoord_t userCoord_t;
91  typedef StridedData<lno_t, scalar_t> input_t;
92 #endif
93 
96 
108  GraphModel(const RCP<const MatrixAdapter<user_t,userCoord_t> > &ia,
109  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
110  modelFlag_t &modelFlags);
111 
112  GraphModel(const RCP<const GraphAdapter<user_t,userCoord_t> > &ia,
113  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
114  modelFlag_t &modelFlags);
115 
116  GraphModel(const RCP<const MeshAdapter<user_t> > &ia,
117  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
118  modelFlag_t &modelflags);
119 
120  GraphModel(const RCP<const VectorAdapter<userCoord_t> > &ia,
121  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
122  modelFlag_t &flags)
123  {
124  throw std::runtime_error("cannot build GraphModel from VectorAdapter");
125  }
126 
127  GraphModel(const RCP<const IdentifierAdapter<user_t> > &ia,
128  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
129  modelFlag_t &flags)
130  {
131  throw std::runtime_error("cannot build GraphModel from IdentifierAdapter");
132  }
133 
136  const RCP<const Comm<int> > getComm() { return comm_; }
137 
140  size_t getLocalNumVertices() const { return nLocalVertices_; }
141 
144  size_t getGlobalNumVertices() const { return nGlobalVertices_; }
145 
149  size_t getLocalNumEdges() const { return nLocalEdges_; }
150 
154  size_t getGlobalNumEdges() const { return nGlobalEdges_; }
155 
158  int getNumWeightsPerVertex() const { return nWeightsPerVertex_; }
159 
162  int getNumWeightsPerEdge() const { return nWeightsPerEdge_; }
163 
166  int getCoordinateDim() const { return vCoordDim_; }
167 
177  ArrayView<const gno_t> &Ids,
178  ArrayView<input_t> &wgts) const
179  {
180  Ids = vGids_.view(0, nLocalVertices_);
181  wgts = vWeights_.view(0, nWeightsPerVertex_);
182  return nLocalVertices_;
183  }
184 
191  size_t getVertexCoords(ArrayView<input_t> &xyz) const
192  {
193  xyz = vCoords_.view(0, vCoordDim_);
194  return nLocalVertices_;
195  }
196 
208  // Implied Vertex LNOs from getVertexList are used as indices to offsets
209  // array.
210  // Vertex GNOs are returned as neighbors in edgeIds.
211 
212  size_t getEdgeList( ArrayView<const gno_t> &edgeIds,
213  ArrayView<const lno_t> &offsets,
214  ArrayView<input_t> &wgts) const
215  {
216  edgeIds = eGids_.view(0, nLocalEdges_);
217  offsets = eOffsets_.view(0, nLocalVertices_+1);
218  wgts = eWeights_.view(0, nWeightsPerEdge_);
219  return nLocalEdges_;
220  }
221 
223  // The Model interface.
225 
226  size_t getLocalNumObjects() const { return nLocalVertices_; }
227 
228  size_t getGlobalNumObjects() const { return nGlobalVertices_; }
229 
230 private:
231 
232  void shared_constructor(const RCP<const Adapter>&ia, modelFlag_t &modelFlags);
233 
234  template <typename AdapterWithCoords>
235  void shared_GetVertexCoords(const AdapterWithCoords *ia);
236 
237  void print(); // For debugging
238 
239  const RCP<const Environment > env_;
240  const RCP<const Comm<int> > comm_;
241 
242  bool localGraph_; // Flag indicating whether this graph is
243  // LOCAL with respect to the process;
244  // if !localGraph_, graph is GLOBAL with respect to
245  // the communicator.
246 
247 
248  size_t nLocalVertices_; // # local vertices in built graph
249  size_t nGlobalVertices_; // # global vertices in built graph
250  ArrayRCP<gno_t> vGids_; // vertices of graph built in model;
251  // may be same as adapter's input
252  // or may be renumbered 0 to (N-1).
253 
254  int nWeightsPerVertex_;
255  ArrayRCP<input_t> vWeights_;
256 
257  int vCoordDim_;
258  ArrayRCP<input_t> vCoords_;
259 
260  // Note: in some cases, size of these arrays
261  // may be larger than nLocalEdges_. So do not use .size().
262  // Use nLocalEdges_, nGlobalEdges_
263 
264  size_t nLocalEdges_; // # local edges in built graph
265  size_t nGlobalEdges_; // # global edges in built graph
266  ArrayRCP<gno_t> eGids_; // edges of graph built in model
267  ArrayRCP<lno_t> eOffsets_; // edge offsets build in model
268  // May be same as adapter's input
269  // or may differ
270  // due to renumbering, self-edge
271  // removal, or local graph.
272 
273  int nWeightsPerEdge_;
274  ArrayRCP<input_t> eWeights_; // edge weights in built graph
275  // May be same as adapter's input
276  // or may differ due to self-edge
277  // removal, or local graph.
278 
279 };
280 
281 
283 // GraphModel from MatrixAdapter
284 template <typename Adapter>
286  const RCP<const MatrixAdapter<user_t,userCoord_t> > &ia,
287  const RCP<const Environment> &env,
288  const RCP<const Comm<int> > &comm,
289  modelFlag_t &modelFlags):
290  env_(env),
291  comm_(comm),
292  localGraph_(false),
293  nLocalVertices_(0),
294  nGlobalVertices_(0),
295  vGids_(),
296  nWeightsPerVertex_(0),
297  vWeights_(),
298  vCoordDim_(0),
299  vCoords_(),
300  nLocalEdges_(0),
301  nGlobalEdges_(0),
302  eGids_(),
303  eOffsets_(),
304  nWeightsPerEdge_(0),
305  eWeights_()
306 {
307  // Model creation flags
308  localGraph_ = modelFlags.test(BUILD_LOCAL_GRAPH);
309 
310  bool symTranspose = modelFlags.test(SYMMETRIZE_INPUT_TRANSPOSE);
311  bool symBipartite = modelFlags.test(SYMMETRIZE_INPUT_BIPARTITE);
312  bool vertexCols = modelFlags.test(VERTICES_ARE_MATRIX_COLUMNS);
313  bool vertexNz = modelFlags.test(VERTICES_ARE_MATRIX_NONZEROS);
314 
315  if (symTranspose || symBipartite || vertexCols || vertexNz){
316  throw std::runtime_error("graph build option not yet implemented");
317  }
318 
319  // Get the matrix from the input adapter
320  gno_t const *vtxIds=NULL, *nborIds=NULL;
321  lno_t const *offsets=NULL;
322  try{
323  nLocalVertices_ = ia->getLocalNumIDs();
324  ia->getIDsView(vtxIds);
325  }
327  try{
328  if (ia->CRSViewAvailable()) {
329  ia->getCRSView(offsets, nborIds);
330  }
331  else {
332  // TODO: Add support for CCS matrix layout
333  throw std::runtime_error("Only MatrixAdapter::getCRSView is supported "
334  "in graph model");
335  }
336  }
338 
339  // Save the pointers from the input adapter
340  nLocalEdges_ = offsets[nLocalVertices_];
341  vGids_ = arcp_const_cast<gno_t>(
342  arcp<const gno_t>(vtxIds, 0, nLocalVertices_, false));
343  eGids_ = arcp_const_cast<gno_t>(
344  arcp<const gno_t>(nborIds, 0, nLocalEdges_, false));
345  eOffsets_ = arcp_const_cast<lno_t>(
346  arcp<const lno_t>(offsets, 0, nLocalVertices_+1, false));
347 
348  // Edge weights
349  nWeightsPerEdge_ = 0; // no edge weights from a matrix yet.
350  // TODO: use matrix values as edge weights
351 
352  // Do constructor common to all adapters
353  shared_constructor(ia, modelFlags);
354 
355  // Get vertex coordinates, if available
356  if (ia->coordinatesAvailable()) {
357  typedef VectorAdapter<userCoord_t> adapterWithCoords_t;
358  shared_GetVertexCoords<adapterWithCoords_t>(ia->getCoordinateInput());
359  }
360  print();
361 }
362 
363 
365 // GraphModel from GraphAdapter
366 template <typename Adapter>
368  const RCP<const GraphAdapter<user_t,userCoord_t> > &ia,
369  const RCP<const Environment> &env,
370  const RCP<const Comm<int> > &comm,
371  modelFlag_t &modelFlags):
372  env_(env),
373  comm_(comm),
374  localGraph_(false),
375  nLocalVertices_(0),
376  nGlobalVertices_(0),
377  vGids_(),
378  nWeightsPerVertex_(0),
379  vWeights_(),
380  vCoordDim_(0),
381  vCoords_(),
382  nLocalEdges_(0),
383  nGlobalEdges_(0),
384  eGids_(),
385  eOffsets_(),
386  nWeightsPerEdge_(0),
387  eWeights_()
388 {
389  // Model creation flags
390  localGraph_ = modelFlags.test(BUILD_LOCAL_GRAPH);
391 
392  // This GraphModel is built with vertices == GRAPH_VERTEX from GraphAdapter.
393  // It is not ready to use vertices == GRAPH_EDGE from GraphAdapter.
394  env_->localInputAssertion(__FILE__, __LINE__,
395  "GraphModel from GraphAdapter is implemented only for "
396  "Graph Vertices as primary object, not for Graph Edges",
397  ia->getPrimaryEntityType() == Zoltan2::GRAPH_VERTEX, BASIC_ASSERTION);
398 
399  // Get the graph from the input adapter
400 
401  gno_t const *vtxIds=NULL, *nborIds=NULL;
402  lno_t const *offsets=NULL;
403  try{
404  nLocalVertices_ = ia->getLocalNumVertices();
405  ia->getVertexIDsView(vtxIds);
406  ia->getEdgesView(offsets, nborIds);
407  }
409 
410  // Save the pointers from the input adapter
411  nLocalEdges_ = offsets[nLocalVertices_];
412  vGids_ = arcp_const_cast<gno_t>(
413  arcp<const gno_t>(vtxIds, 0, nLocalVertices_, false));
414  eGids_ = arcp_const_cast<gno_t>(
415  arcp<const gno_t>(nborIds, 0, nLocalEdges_, false));
416  eOffsets_ = arcp_const_cast<lno_t>(
417  arcp<const lno_t>(offsets, 0, nLocalVertices_+1, false));
418 
419  // Edge weights
420  nWeightsPerEdge_ = ia->getNumWeightsPerEdge();
421  if (nWeightsPerEdge_ > 0){
422  input_t *wgts = new input_t [nWeightsPerEdge_];
423  eWeights_ = arcp(wgts, 0, nWeightsPerEdge_, true);
424  }
425 
426  for (int w=0; w < nWeightsPerEdge_; w++){
427  const scalar_t *ewgts=NULL;
428  int stride=0;
429 
430  ia->getEdgeWeightsView(ewgts, stride, w);
431 
432  ArrayRCP<const scalar_t> wgtArray(ewgts, 0, nLocalEdges_, false);
433  eWeights_[w] = input_t(wgtArray, stride);
434  }
435 
436  // Do constructor common to all adapters
437  shared_constructor(ia, modelFlags);
438 
439  // Get vertex coordinates, if available
440  if (ia->coordinatesAvailable()) {
441  typedef VectorAdapter<userCoord_t> adapterWithCoords_t;
442  shared_GetVertexCoords<adapterWithCoords_t>(ia->getCoordinateInput());
443  }
444  print();
445 }
446 
448 // GraphModel from MeshAdapter
449 template <typename Adapter>
451  const RCP<const MeshAdapter<user_t> > &ia,
452  const RCP<const Environment> &env,
453  const RCP<const Comm<int> > &comm,
454  modelFlag_t &modelFlags):
455  env_(env),
456  comm_(comm),
457  localGraph_(false),
458  nLocalVertices_(0),
459  nGlobalVertices_(0),
460  vGids_(),
461  nWeightsPerVertex_(0),
462  vWeights_(),
463  vCoordDim_(0),
464  vCoords_(),
465  nLocalEdges_(0),
466  nGlobalEdges_(0),
467  eGids_(),
468  eOffsets_(),
469  nWeightsPerEdge_(0),
470  eWeights_()
471 {
472  env_->timerStart(MACRO_TIMERS, "GraphModel constructed from MeshAdapter");
473 
474  // Model creation flags
475  localGraph_ = modelFlags.test(BUILD_LOCAL_GRAPH);
476 
477  // This GraphModel is built with vertices == ia->getPrimaryEntityType()
478  // from MeshAdapter.
479 
480  // Get the graph from the input adapter
481 
482  Zoltan2::MeshEntityType primaryEType = ia->getPrimaryEntityType();
483  Zoltan2::MeshEntityType secondAdjEType = ia->getSecondAdjacencyEntityType();
484 
485  // Get the IDs of the primary entity type; these are graph vertices
486 
487  gno_t const *vtxIds=NULL;
488  try {
489  nLocalVertices_ = ia->getLocalNumOf(primaryEType);
490  ia->getIDsViewOf(primaryEType, vtxIds);
491  }
493 
494  vGids_ = arcp_const_cast<gno_t>(
495  arcp<const gno_t>(vtxIds, 0, nLocalVertices_, false));
496 
497  // Get the second adjacencies to construct edges of the dual graph.
498  gno_t const *nborIds=NULL;
499  lno_t const *offsets=NULL;
500 
501  if (!ia->avail2ndAdjs(primaryEType, secondAdjEType)) {
502  // KDDKDD TODO Want to do this differently for local and global graphs?
503  // KDDKDD TODO Currently getting global 2nd Adjs and filtering them for
504  // KDDKDD TODO local graphs. That approach is consistent with other
505  // KDDKDD TODO adapters, but is more expensive -- why build them just to
506  // KDDKDD TODO throw them away? Instead, perhaps should build
507  // KDDKDD TODO only local adjacencies.
508  // KDDKDD TODO Does it suffice to pass a serial comm for local graph?
509  try {
510  get2ndAdjsViewFromAdjs(ia, comm_, primaryEType, secondAdjEType, offsets,
511  nborIds);
512  }
514  }
515  else { // avail2ndAdjs
516  // Get the edges
517  try {
518  ia->get2ndAdjsView(primaryEType, secondAdjEType, offsets, nborIds);
519  }
521  }
522 
523  // Save the pointers from the input adapter
524  nLocalEdges_ = offsets[nLocalVertices_];
525  eGids_ = arcp_const_cast<gno_t>(
526  arcp<const gno_t>(nborIds, 0, nLocalEdges_, false));
527  eOffsets_ = arcp_const_cast<lno_t>(
528  arcp<const lno_t>(offsets, 0, nLocalVertices_+1, false));
529 
530  // Edge weights
531  // Cannot specify edge weights if Zoltan2 computes the second adjacencies;
532  // there's no way to know the correct order for the adjacencies and weights.
533  // InputAdapter must provide 2nd adjs in order for edge weights to be used.
534  if (ia->avail2ndAdjs(primaryEType, secondAdjEType)) {
535  nWeightsPerEdge_ = ia->getNumWeightsPer2ndAdj(primaryEType, secondAdjEType);
536  if (nWeightsPerEdge_ > 0){
537  input_t *wgts = new input_t [nWeightsPerEdge_];
538  eWeights_ = arcp(wgts, 0, nWeightsPerEdge_, true);
539  }
540 
541  for (int w=0; w < nWeightsPerEdge_; w++){
542  const scalar_t *ewgts=NULL;
543  int stride=0;
544 
545  ia->get2ndAdjWeightsView(primaryEType, secondAdjEType,
546  ewgts, stride, w);
547 
548  ArrayRCP<const scalar_t> wgtArray(ewgts, 0,
549  nLocalEdges_*stride, false);
550  eWeights_[w] = input_t(wgtArray, stride);
551  }
552  }
553 
554  // Do constructor common to all adapters
555  shared_constructor(ia, modelFlags);
556 
557  // Get vertex coordinates
558  typedef MeshAdapter<user_t> adapterWithCoords_t;
559  shared_GetVertexCoords<adapterWithCoords_t>(&(*ia));
560 
561  env_->timerStop(MACRO_TIMERS, "GraphModel constructed from MeshAdapter");
562  print();
563 }
564 
565 
567 template <typename Adapter>
569  const RCP<const Adapter> &ia,
570  modelFlag_t &modelFlags)
571 {
572  // Model creation flags
573  bool consecutiveIdsRequired = modelFlags.test(GENERATE_CONSECUTIVE_IDS);
574  bool removeSelfEdges = modelFlags.test(REMOVE_SELF_EDGES);
575  bool subsetGraph = modelFlags.test(BUILD_SUBSET_GRAPH);
576 
577  // May modify the graph provided from input adapter; save pointers to
578  // the input adapter's data.
579  size_t adapterNLocalEdges = nLocalEdges_;
580  ArrayRCP<gno_t> adapterVGids = vGids_; // vertices of graph from adapter
581  ArrayRCP<gno_t> adapterEGids = eGids_; // edges of graph from adapter
582  ArrayRCP<lno_t> adapterEOffsets = eOffsets_; // edge offsets from adapter
583  ArrayRCP<input_t> adapterEWeights = eWeights_; // edge weights from adapter
584 
585  if (localGraph_) {
586  // Local graph is requested.
587  // Renumber vertices 0 to nLocalVertices_-1
588  // Filter out off-process edges
589  // Renumber edge neighbors to be in range [0,nLocalVertices_-1]
590 
591  // Allocate new space for local graph info
592  // Note that eGids_ and eWeights_[w] may be larger than needed;
593  // we would have to pre-count the number of local edges to avoid overalloc
594  vGids_ = arcp(new gno_t[nLocalVertices_],
595  0, nLocalVertices_, true);
596  eGids_ = arcp(new gno_t[adapterNLocalEdges],
597  0, adapterNLocalEdges, true);
598  eOffsets_ = arcp(new lno_t[nLocalVertices_+1],
599  0, nLocalVertices_+1, true);
600 
601  scalar_t **tmpEWeights = NULL;
602  if (nWeightsPerEdge_ > 0){
603  eWeights_ = arcp(new input_t[nWeightsPerEdge_], 0,
604  nWeightsPerEdge_, true);
605  // Need to use temporary array because StridedData has const data
606  // so we can't write to it.
607  tmpEWeights = new scalar_t*[nWeightsPerEdge_];
608  for (int w = 0; w < nWeightsPerEdge_; w++)
609  tmpEWeights[w] = new scalar_t[adapterNLocalEdges];
610  }
611 
612  // Build map between global and local vertex numbers
613  std::unordered_map<gno_t, lno_t> globalToLocal(nLocalVertices_);
614  for (size_t i = 0; i < nLocalVertices_; i++)
615  globalToLocal[adapterVGids[i]] = i;
616 
617  // Loop over edges; keep only those that are local (i.e., on-rank)
618  eOffsets_[0] = 0;
619  lno_t ecnt = 0;
620  for (size_t i = 0; i < nLocalVertices_; i++) {
621  vGids_[i] = gno_t(i);
622  for (lno_t j = adapterEOffsets[i]; j < adapterEOffsets[i+1]; j++) {
623 
624  if (removeSelfEdges && (adapterEGids[j] == adapterVGids[i]))
625  continue; // Skipping self edge
626 
627  // Determine whether neighbor vertex is local
628  typename std::unordered_map<gno_t, lno_t>::iterator localidx;
629  if ((localidx = globalToLocal.find(adapterEGids[j])) !=
630  globalToLocal.end()) {
631  // neighbor vertex is local
632  // Keep the edge and its weights
633  eGids_[ecnt] = localidx->second;
634  for (int w = 0; w < nWeightsPerEdge_; w++)
635  tmpEWeights[w][ecnt] = adapterEWeights[w][j];
636 
637  ecnt++;
638  }
639  }
640  eOffsets_[i+1] = ecnt;
641  }
642  nLocalEdges_ = eOffsets_[nLocalVertices_];
643  if (nWeightsPerEdge_) {
644  for (int w = 0; w < nWeightsPerEdge_; w++) {
645  ArrayRCP<const scalar_t> wgtArray(tmpEWeights[w],
646  0, adapterNLocalEdges, true);
647  eWeights_[w] = input_t(wgtArray, 0);
648  }
649  delete [] tmpEWeights;
650  }
651  } // localGraph_
652 
653  else if (consecutiveIdsRequired || removeSelfEdges || subsetGraph) {
654  // Build a Global graph
655  // If we are here, we expect SOMETHING in the graph to change from input:
656  // removing self edges, or converting to consecutive IDs, or subsetting
657  // the graph.
658 
659 
660  // Determine vertex GIDs for the global GraphModel
661  Teuchos::ArrayRCP<size_t> vtxDist; // TODO keep vtxDist & add to interface?
662  if (consecutiveIdsRequired) {
663  // Allocate new memory for vertices for consecutiveIds
664  vGids_ = arcp(new gno_t[nLocalVertices_], 0, nLocalVertices_, true);
665 
666  // Build vtxDist array with starting vGid on each rank
667  int np = comm_->getSize();
668  vtxDist = arcp(new size_t[np+1], 0, np+1, true);
669  vtxDist[0] = 0;
670  Teuchos::gatherAll(*comm_, 1, &nLocalVertices_, np, &vtxDist[1]);
671  for (int i = 0; i < np; i++)
672  vtxDist[i+1] += vtxDist[i];
673  }
674 
675  // Allocate new memory for edges and offsets, as needed
676  // Note that eGids_ may or may not be larger than needed;
677  // would have to pre-count number of edges kept otherwise
678  eGids_ = arcp(new gno_t[adapterNLocalEdges],
679  0, adapterNLocalEdges, true);
680 
681  scalar_t **tmpEWeights = NULL;
682  if (subsetGraph || removeSelfEdges) {
683  // May change number of edges and, thus, the offsets
684  eOffsets_ = arcp(new lno_t[nLocalVertices_+1],
685  0, nLocalVertices_+1, true);
686  eOffsets_[0] = 0;
687 
688  // Need to copy weights if remove edges
689  if (nWeightsPerEdge_ > 0){
690  eWeights_ = arcp(new input_t[nWeightsPerEdge_], 0,
691  nWeightsPerEdge_, true);
692  // Need to use temporary array because StridedData has const data
693  // so we can't write to it.
694  tmpEWeights = new scalar_t*[nWeightsPerEdge_];
695  for (int w = 0; w < nWeightsPerEdge_; w++)
696  tmpEWeights[w] = new scalar_t[adapterNLocalEdges];
697  }
698  }
699 
700  // If needed, determine the owning ranks and its local index off-proc
701  Teuchos::ArrayRCP<int> edgeRemoteRanks;
702  Teuchos::ArrayRCP<lno_t> edgeRemoteLids;
703  std::unordered_map<gno_t, size_t> edgeRemoteUniqueMap;
704 
705  if (subsetGraph || consecutiveIdsRequired) {
706  gno_t dummy = Teuchos::OrdinalTraits<gno_t>::invalid();
707  Tpetra::Map<lno_t,gno_t> vtxMap(dummy, adapterVGids(), 0, comm_);
708 
709  // Need to filter requested edges to make a unique list,
710  // as Tpetra::Map does not return correct info for duplicated entries
711  // (See bug 6412)
712  // The local filter may be more efficient anyway -- fewer communicated
713  // values in the Tpetra directory
714  Teuchos::ArrayRCP<gno_t> edgeRemoteUniqueGids =
715  arcp(new gno_t[adapterNLocalEdges], 0, adapterNLocalEdges, true);
716 
717  size_t nEdgeUnique = 0;
718  for (size_t i = 0; i < adapterNLocalEdges; i++) {
719  if (edgeRemoteUniqueMap.find(adapterEGids[i]) ==
720  edgeRemoteUniqueMap.end()) {
721  edgeRemoteUniqueGids[nEdgeUnique] = adapterEGids[i];
722  edgeRemoteUniqueMap[adapterEGids[i]] = nEdgeUnique;
723  nEdgeUnique++;
724  }
725  }
726 
727  edgeRemoteRanks = arcp(new int[nEdgeUnique], 0, nEdgeUnique, true);
728  edgeRemoteLids = arcp(new lno_t[nEdgeUnique], 0, nEdgeUnique, true);
729  vtxMap.getRemoteIndexList(edgeRemoteUniqueGids(0, nEdgeUnique),
730  edgeRemoteRanks(), edgeRemoteLids());
731  }
732 
733  // Renumber and/or filter the edges and vertices
734  lno_t ecnt = 0;
735  int me = comm_->getRank();
736  for (size_t i = 0; i < nLocalVertices_; i++) {
737 
738  if (consecutiveIdsRequired)
739  vGids_[i] = vtxDist[me] + i;
740 
741  for (lno_t j = adapterEOffsets[i]; j < adapterEOffsets[i+1]; j++) {
742 
743  if (removeSelfEdges && (adapterVGids[i] == adapterEGids[j]))
744  continue; // Skipping self edge
745 
746  size_t remoteIdx = edgeRemoteUniqueMap[adapterEGids[j]];
747 
748  if (subsetGraph && (edgeRemoteRanks[remoteIdx] == -1))
749  continue; // Skipping edge with neighbor vertex that was not found
750  // in communicator
751 
752  if (consecutiveIdsRequired)
753  // Renumber edge using local number on remote rank plus first
754  // vtx number for remote rank
755  eGids_[ecnt] = edgeRemoteLids[remoteIdx]
756  + vtxDist[edgeRemoteRanks[remoteIdx]];
757  else
758  eGids_[ecnt] = adapterEGids[j];
759 
760  if (subsetGraph || removeSelfEdges) {
761  // Need to copy weights only if number of edges might change
762  for (int w = 0; w < nWeightsPerEdge_; w++)
763  tmpEWeights[w][ecnt] = adapterEWeights[w][j];
764  }
765 
766  ecnt++;
767  }
768  if (subsetGraph || removeSelfEdges)
769  eOffsets_[i+1] = ecnt;
770  }
771  nLocalEdges_ = ecnt;
772  if (nWeightsPerEdge_ && (subsetGraph || removeSelfEdges)) {
773  for (int w = 0; w < nWeightsPerEdge_; w++) {
774  ArrayRCP<const scalar_t> wgtArray(tmpEWeights[w],
775  0, nLocalEdges_, true);
776  eWeights_[w] = input_t(wgtArray, 1);
777  }
778  delete [] tmpEWeights;
779  }
780  }
781 
782  // Vertex weights
783  nWeightsPerVertex_ = ia->getNumWeightsPerID();
784 
785  if (nWeightsPerVertex_ > 0){
786  input_t *weightInfo = new input_t [nWeightsPerVertex_];
787  env_->localMemoryAssertion(__FILE__, __LINE__, nWeightsPerVertex_,
788  weightInfo);
789 
790  for (int idx=0; idx < nWeightsPerVertex_; idx++){
791  bool useNumNZ = ia->useDegreeAsWeight(idx);
792  if (useNumNZ){
793  scalar_t *wgts = new scalar_t [nLocalVertices_];
794  env_->localMemoryAssertion(__FILE__, __LINE__, nLocalVertices_, wgts);
795  ArrayRCP<const scalar_t> wgtArray = arcp(wgts,
796  0, nLocalVertices_, true);
797  for (size_t i=0; i < nLocalVertices_; i++)
798  wgts[i] = eOffsets_[i+1] - eOffsets_[i];
799  weightInfo[idx] = input_t(wgtArray, 1);
800  }
801  else{
802  const scalar_t *weights=NULL;
803  int stride=0;
804  ia->getWeightsView(weights, stride, idx);
805  ArrayRCP<const scalar_t> wgtArray = arcp(weights, 0,
806  stride*nLocalVertices_,
807  false);
808  weightInfo[idx] = input_t(wgtArray, stride);
809  }
810  }
811 
812  vWeights_ = arcp<input_t>(weightInfo, 0, nWeightsPerVertex_, true);
813  }
814 
815  // Compute global values
816  if (localGraph_) {
817  nGlobalVertices_ = nLocalVertices_;
818  nGlobalEdges_ = nLocalEdges_;
819  }
820  else {
821  reduceAll<int, size_t>(*comm_, Teuchos::REDUCE_SUM, 1,
822  &nLocalVertices_, &nGlobalVertices_);
823  reduceAll<int, size_t>(*comm_, Teuchos::REDUCE_SUM, 1,
824  &nLocalEdges_, &nGlobalEdges_);
825  }
826 
827  env_->memory("After construction of graph model");
828 }
829 
831 
832 template <typename Adapter>
833 template <typename AdapterWithCoords>
834 void GraphModel<Adapter>::shared_GetVertexCoords(const AdapterWithCoords *ia)
835 {
836  // get pointers to vertex coordinates from input adapter
837 
838  vCoordDim_ = ia->getDimension();
839 
840  if (vCoordDim_ > 0){
841  input_t *coordInfo = new input_t [vCoordDim_];
842  env_->localMemoryAssertion(__FILE__, __LINE__, vCoordDim_, coordInfo);
843 
844  for (int dim=0; dim < vCoordDim_; dim++){
845  const scalar_t *coords=NULL;
846  int stride=0;
847  ia->getCoordinatesView(coords, stride, dim);
848  ArrayRCP<const scalar_t> coordArray = arcp(coords, 0,
849  stride*nLocalVertices_,
850  false);
851  coordInfo[dim] = input_t(coordArray, stride);
852  }
853 
854  vCoords_ = arcp<input_t>(coordInfo, 0, vCoordDim_, true);
855  }
856 }
857 
859  template <typename Adapter>
861 {
862 // if (env_->getDebugLevel() < VERBOSE_DETAILED_STATUS)
863 // return;
864 
865  std::ostream *os = env_->getDebugOStream();
866 
867  int me = comm_->getRank();
868 
869  *os << me
870  << " " << (localGraph_ ? "LOCAL GRAPH " : "GLOBAL GRAPH ")
871  << " Nvtx " << nLocalVertices_
872  << " Nedge " << nLocalEdges_
873  << " NVWgt " << nWeightsPerVertex_
874  << " NEWgt " << nWeightsPerEdge_
875  << " CDim " << vCoordDim_
876  << std::endl;
877 
878  for (size_t i = 0; i < nLocalVertices_; i++) {
879  *os << me << " " << i << " GID " << vGids_[i] << ": ";
880  for (lno_t j = eOffsets_[i]; j < eOffsets_[i+1]; j++)
881  *os << eGids_[j] << " " ;
882  *os << std::endl;
883  }
884 
885  if (nWeightsPerVertex_) {
886  for (size_t i = 0; i < nLocalVertices_; i++) {
887  *os << me << " " << i << " VWGTS " << vGids_[i] << ": ";
888  for (int j = 0; j < nWeightsPerVertex_; j++)
889  *os << vWeights_[j][i] << " ";
890  *os << std::endl;
891  }
892  }
893 
894  if (nWeightsPerEdge_) {
895  for (size_t i = 0; i < nLocalVertices_; i++) {
896  *os << me << " " << i << " EWGTS " << vGids_[i] << ": ";
897  for (lno_t j = eOffsets_[i]; j < eOffsets_[i+1]; j++) {
898  *os << eGids_[j] << " (";
899  for (int w = 0; w < nWeightsPerEdge_; w++)
900  *os << eWeights_[w][j] << " ";
901  *os << ") ";
902  }
903  *os << std::endl;
904  }
905  }
906 
907  if (vCoordDim_) {
908  for (size_t i = 0; i < nLocalVertices_; i++) {
909  *os << me << " " << i << " COORDS " << vGids_[i] << ": ";
910  for (int j = 0; j < vCoordDim_; j++)
911  *os << vCoords_[j][i] << " ";
912  *os << std::endl;
913  }
914  }
915  else
916  *os << me << " NO COORDINATES AVAIL " << std::endl;
917 }
918 
919 } // namespace Zoltan2
920 
921 
922 #endif
923 
use columns as graph vertices
size_t getLocalNumEdges() const
Returns the number of edges on this process. In global or subset graphs, includes off-process edges...
size_t getVertexList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
Sets pointers to this process&#39; vertex Ids and their weights.
Time an algorithm (or other entity) as a whole.
IdentifierAdapter defines the interface for identifiers.
ignore invalid neighbors
int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
size_t getGlobalNumObjects() const
Return the global number of objects.
fast typical checks for valid arguments
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
GraphModel(const RCP< const VectorAdapter< userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
MatrixAdapter defines the adapter interface for matrices.
Defines the Model interface.
GraphAdapter defines the interface for graph-based user data.
algorithm requires consecutive ids
Defines the MeshAdapter interface.
MeshAdapter defines the interface for mesh input.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
model must symmetrize input
Defines the IdentifierAdapter interface.
Defines the VectorAdapter interface.
GraphModel(const RCP< const MatrixAdapter< user_t, userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &modelFlags)
Constructor.
model must symmetrize input
size_t getGlobalNumVertices() const
Returns the global number vertices.
Defines helper functions for use in the models.
algorithm requires no self edges
int getCoordinateDim() const
Returns the dimension (0 to 3) of vertex coordinates.
size_t getVertexCoords(ArrayView< input_t > &xyz) const
Sets pointers to this process&#39; vertex coordinates, if available. Order of coordinate info matches tha...
void get2ndAdjsViewFromAdjs(const Teuchos::RCP< const MeshAdapter< User > > &ia, const RCP< const Comm< int > > comm, Zoltan2::MeshEntityType sourcetarget, Zoltan2::MeshEntityType through, const typename MeshAdapter< User >::lno_t *&offsets, const typename MeshAdapter< User >::gno_t *&adjacencyIds)
use nonzeros as graph vertices
size_t getEdgeList(ArrayView< const gno_t > &edgeIds, ArrayView< const lno_t > &offsets, ArrayView< input_t > &wgts) const
Sets pointers to this process&#39; edge (neighbor) global Ids, including off-process edges.
size_t getGlobalNumEdges() const
Returns the global number edges. For local graphs, the number of global edges is the number of local ...
VectorAdapter defines the interface for vector input.
The StridedData class manages lists of weights or coordinates.
Traits for application input objects.
size_t getLocalNumVertices() const
Returns the number vertices on this process.
GraphModel defines the interface required for graph models.
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
size_t getLocalNumObjects() const
Return the local number of objects.
Defines the MatrixAdapter interface.
The base class for all model classes.
int getNumWeightsPerEdge() const
Returns the number (0 or greater) of weights per edge.
Defines the GraphAdapter interface.
GraphModel(const RCP< const IdentifierAdapter< user_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
const RCP< const Comm< int > > getComm()
Return the communicator used by the model.
model represents graph within only one rank
This file defines the StridedData class.