Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_LinePartitioner_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_LINE_PARTITIONER_DEF_HPP
11#define IFPACK2_LINE_PARTITIONER_DEF_HPP
12
13#include "Tpetra_CrsGraph.hpp"
14#include "Tpetra_Util.hpp"
15
16namespace Ifpack2 {
17
18template<class T>
19inline typename Teuchos::ScalarTraits<T>::magnitudeType square(T x) {
20 return Teuchos::ScalarTraits<T>::magnitude(x) * Teuchos::ScalarTraits<T>::magnitude(x);
21}
22
23//==============================================================================
24// Constructor
25template<class GraphType,class Scalar>
27LinePartitioner (const Teuchos::RCP<const row_graph_type>& graph) :
28 OverlappingPartitioner<GraphType> (graph), NumEqns_(1), threshold_(0.0)
29{}
30
31
32template<class GraphType,class Scalar>
34
35
36template<class GraphType,class Scalar>
37void
39setPartitionParameters(Teuchos::ParameterList& List) {
40 threshold_ = List.get("partitioner: line detection threshold",threshold_);
41 TEUCHOS_TEST_FOR_EXCEPTION(threshold_ < 0.0 || threshold_ > 1.0,
42 std::runtime_error,"Ifpack2::LinePartitioner: threshold not valid");
43
44 NumEqns_ = List.get("partitioner: PDE equations",NumEqns_);
45 TEUCHOS_TEST_FOR_EXCEPTION(NumEqns_<1,std::runtime_error,"Ifpack2::LinePartitioner: NumEqns not valid");
46
47 coord_ = List.get("partitioner: coordinates",coord_);
48 TEUCHOS_TEST_FOR_EXCEPTION(coord_.is_null(),std::runtime_error,"Ifpack2::LinePartitioner: coordinates not defined");
49}
50
51
52template<class GraphType,class Scalar>
54 const local_ordinal_type invalid = Teuchos::OrdinalTraits<local_ordinal_type>::invalid();
55
56 // Sanity Checks
57 TEUCHOS_TEST_FOR_EXCEPTION(coord_.is_null(),std::runtime_error,"Ifpack2::LinePartitioner: coordinates not defined");
58 TEUCHOS_TEST_FOR_EXCEPTION((size_t)this->Partition_.size() != this->Graph_->getLocalNumRows(),std::runtime_error,"Ifpack2::LinePartitioner: partition size error");
59
60 // Short circuit
61 if(this->Partition_.size() == 0) {this->NumLocalParts_ = 0; return;}
62
63 // Set partitions to invalid to initialize algorithm
64 for(size_t i=0; i<this->Graph_->getLocalNumRows(); i++)
65 this->Partition_[i] = invalid;
66
67 // Use the auto partitioner
68 this->NumLocalParts_ = this->Compute_Blocks_AutoLine(this->Partition_());
69
70 // Resize Parts_
71 this->Parts_.resize(this->NumLocalParts_);
72}
73
74
75// ============================================================================
76template<class GraphType,class Scalar>
77int LinePartitioner<GraphType,Scalar>::Compute_Blocks_AutoLine(Teuchos::ArrayView<local_ordinal_type> blockIndices) const {
78 typedef local_ordinal_type LO;
79 typedef magnitude_type MT;
80 const LO invalid = Teuchos::OrdinalTraits<LO>::invalid();
81 const MT zero = Teuchos::ScalarTraits<MT>::zero();
82
83 Teuchos::ArrayRCP<const MT> xvalsRCP, yvalsRCP, zvalsRCP;
84 Teuchos::ArrayView<const MT> xvals, yvals, zvals;
85 xvalsRCP = coord_->getData(0); xvals = xvalsRCP();
86 if(coord_->getNumVectors() > 1) { yvalsRCP = coord_->getData(1); yvals = yvalsRCP(); }
87 if(coord_->getNumVectors() > 2) { zvalsRCP = coord_->getData(2); zvals = zvalsRCP(); }
88
89 double tol = threshold_;
90 size_t N = this->Graph_->getLocalNumRows();
91 size_t allocated_space = this->Graph_->getLocalMaxNumRowEntries();
92
93 nonconst_local_inds_host_view_type cols("cols",allocated_space);
94 Teuchos::Array<LO> indices(allocated_space);
95 Teuchos::Array<MT> dist(allocated_space);
96
97 Teuchos::Array<LO> itemp(2*allocated_space);
98 Teuchos::Array<MT> dtemp(allocated_space);
99
100 LO num_lines = 0;
101
102 for(LO i=0; i<(LO)N; i+=NumEqns_) {
103 size_t nz=0;
104 // Short circuit if I've already been blocked
105 if(blockIndices[i] != invalid) continue;
106
107 // Get neighbors and sort by distance
108 this->Graph_->getLocalRowCopy(i,cols,nz);
109 MT x0 = (!xvals.is_null()) ? xvals[i/NumEqns_] : zero;
110 MT y0 = (!yvals.is_null()) ? yvals[i/NumEqns_] : zero;
111 MT z0 = (!zvals.is_null()) ? zvals[i/NumEqns_] : zero;
112
113 LO neighbor_len=0;
114 for(size_t j=0; j<nz; j+=NumEqns_) {
115 MT mydist = zero;
116 LO nn = cols[j] / NumEqns_;
117 if(cols[j] >=(LO)N) continue; // Check for off-proc entries
118 if(!xvals.is_null()) mydist += square<MT>(x0 - xvals[nn]);
119 if(!yvals.is_null()) mydist += square<MT>(y0 - yvals[nn]);
120 if(!zvals.is_null()) mydist += square<MT>(z0 - zvals[nn]);
121 dist[neighbor_len] = Teuchos::ScalarTraits<MT>::squareroot(mydist);
122 indices[neighbor_len]=cols[j];
123 neighbor_len++;
124 }
125
126 Teuchos::ArrayView<MT> dist_view = dist(0,neighbor_len);
127 Tpetra::sort2(dist_view.begin(),dist_view.end(),indices.begin());
128
129 // Number myself
130 for(LO k=0; k<NumEqns_; k++)
131 blockIndices[i + k] = num_lines;
132
133 // Fire off a neighbor line search (nearest neighbor)
134 if(neighbor_len > 2 && dist[1]/dist[neighbor_len-1] < tol) {
135 local_automatic_line_search(NumEqns_,blockIndices,i,indices[1],num_lines,tol,itemp,dtemp);
136 }
137 // Fire off a neighbor line search (second nearest neighbor)
138 if(neighbor_len > 3 && dist[2]/dist[neighbor_len-1] < tol) {
139 local_automatic_line_search(NumEqns_,blockIndices,i,indices[2],num_lines,tol,itemp,dtemp);
140 }
141
142 num_lines++;
143 }
144 return num_lines;
145}
146// ============================================================================
147template<class GraphType,class Scalar>
148void LinePartitioner<GraphType,Scalar>::local_automatic_line_search(int NumEqns, Teuchos::ArrayView <local_ordinal_type> blockIndices, local_ordinal_type last, local_ordinal_type next, local_ordinal_type LineID, double tol, Teuchos::Array<local_ordinal_type> itemp, Teuchos::Array<magnitude_type> dtemp) const {
149 typedef local_ordinal_type LO;
150 typedef magnitude_type MT;
151 const LO invalid = Teuchos::OrdinalTraits<LO>::invalid();
152 const MT zero = Teuchos::ScalarTraits<MT>::zero();
153
154 Teuchos::ArrayRCP<const MT> xvalsRCP, yvalsRCP, zvalsRCP;
155 Teuchos::ArrayView<const MT> xvals, yvals, zvals;
156 xvalsRCP = coord_->getData(0); xvals = xvalsRCP();
157 if(coord_->getNumVectors() > 1) { yvalsRCP = coord_->getData(1); yvals = yvalsRCP(); }
158 if(coord_->getNumVectors() > 2) { zvalsRCP = coord_->getData(2); zvals = zvalsRCP(); }
159
160 size_t N = this->Graph_->getLocalNumRows();
161 size_t allocated_space = this->Graph_->getLocalMaxNumRowEntries();
162
163 nonconst_local_inds_host_view_type cols(itemp.data(),allocated_space);
164 Teuchos::ArrayView<LO> indices = itemp.view(allocated_space,allocated_space);
165 Teuchos::ArrayView<MT> dist= dtemp();
166
167 while (blockIndices[next] == invalid) {
168 // Get the next row
169 size_t nz=0;
170 LO neighbors_in_line=0;
171
172 this->Graph_->getLocalRowCopy(next,cols,nz);
173 MT x0 = (!xvals.is_null()) ? xvals[next/NumEqns_] : zero;
174 MT y0 = (!yvals.is_null()) ? yvals[next/NumEqns_] : zero;
175 MT z0 = (!zvals.is_null()) ? zvals[next/NumEqns_] : zero;
176
177 // Calculate neighbor distances & sort
178 LO neighbor_len=0;
179 for(size_t i=0; i<nz; i+=NumEqns) {
180 MT mydist = zero;
181 if(cols[i] >=(LO)N) continue; // Check for off-proc entries
182 LO nn = cols[i] / NumEqns;
183 if(blockIndices[nn]==LineID) neighbors_in_line++;
184 if(!xvals.is_null()) mydist += square<MT>(x0 - xvals[nn]);
185 if(!yvals.is_null()) mydist += square<MT>(y0 - yvals[nn]);
186 if(!zvals.is_null()) mydist += square<MT>(z0 - zvals[nn]);
187 dist[neighbor_len] = Teuchos::ScalarTraits<MT>::squareroot(mydist);
188 indices[neighbor_len]=cols[i];
189 neighbor_len++;
190 }
191 // If more than one of my neighbors is already in this line. I
192 // can't be because I'd create a cycle
193 if(neighbors_in_line > 1) break;
194
195 // Otherwise add me to the line
196 for(LO k=0; k<NumEqns; k++)
197 blockIndices[next + k] = LineID;
198
199 // Try to find the next guy in the line (only check the closest two that aren't element 0 (diagonal))
200 Teuchos::ArrayView<MT> dist_view = dist(0,neighbor_len);
201 Tpetra::sort2(dist_view.begin(),dist_view.end(),indices.begin());
202
203 if(neighbor_len > 2 && indices[1] != last && blockIndices[indices[1]] == -1 && dist[1]/dist[neighbor_len-1] < tol) {
204 last=next;
205 next=indices[1];
206 }
207 else if(neighbor_len > 3 && indices[2] != last && blockIndices[indices[2]] == -1 && dist[2]/dist[neighbor_len-1] < tol) {
208 last=next;
209 next=indices[2];
210 }
211 else {
212 // I have no further neighbors in this line
213 break;
214 }
215 }
216}
217
218
219
220}// namespace Ifpack2
221
222#define IFPACK2_LINEPARTITIONER_INSTANT(S,LO,GO,N) \
223 template class Ifpack2::LinePartitioner<Tpetra::CrsGraph< LO, GO, N >,S >; \
224 template class Ifpack2::LinePartitioner<Tpetra::RowGraph< LO, GO, N >,S >;
225
226#endif // IFPACK2_LINEPARTITIONER_DEF_HPP
LinePartitioner(const Teuchos::RCP< const row_graph_type > &graph)
Constructor.
Definition Ifpack2_LinePartitioner_def.hpp:27
void setPartitionParameters(Teuchos::ParameterList &List)
Set the partitioner's parameters (none for linear partitioning).
Definition Ifpack2_LinePartitioner_def.hpp:39
virtual ~LinePartitioner()
Destructor.
Definition Ifpack2_LinePartitioner_def.hpp:33
void computePartitions()
Compute the partitions.
Definition Ifpack2_LinePartitioner_def.hpp:53
Teuchos::Array< local_ordinal_type > Partition_
Mapping from local row to partition number.
Definition Ifpack2_OverlappingPartitioner_decl.hpp:145
OverlappingPartitioner(const Teuchos::RCP< const row_graph_type > &graph)
Constructor.
Definition Ifpack2_OverlappingPartitioner_def.hpp:24
Teuchos::RCP< const row_graph_type > Graph_
The graph to be partitioned.
Definition Ifpack2_OverlappingPartitioner_decl.hpp:155
Teuchos::Array< Teuchos::ArrayRCP< local_ordinal_type > > Parts_
Mapping from partition to all rows it contains.
Definition Ifpack2_OverlappingPartitioner_decl.hpp:152
int NumLocalParts_
Number of local subgraphs.
Definition Ifpack2_OverlappingPartitioner_decl.hpp:138
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41