File: ParSparseMatrix.cpp

package info (click to toggle)
lammps 20220106.git7586adbb6a%2Bds1-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 348,064 kB
  • sloc: cpp: 831,421; python: 24,896; xml: 14,949; f90: 10,845; ansic: 7,967; sh: 4,226; perl: 4,064; fortran: 2,424; makefile: 1,501; objc: 238; lisp: 163; csh: 16; awk: 14; tcl: 6
file content (302 lines) | stat: -rw-r--r-- 9,847 bytes parent folder | download | duplicates (2)
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
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
#include "ParSparseMatrix.h"
#include <fstream>
#ifdef TIMING_ON
double time_diff(timespec &start, timespec &end)
{
  return (double)(1e9 * (end.tv_sec - start.tv_sec) +
                  end.tv_nsec - start.tv_nsec) / 1e9;
}
#endif
using namespace MPI_Wrappers;
using namespace std;
namespace ATC_matrix {

// All the same constructors as for SparseMatrix
ParSparseMatrix<double>::ParSparseMatrix(MPI_Comm comm, INDEX rows, INDEX cols)
  : SparseMatrix<double>(rows, cols), _comm(comm) {}

ParSparseMatrix<double>::ParSparseMatrix(MPI_Comm comm,
    const SparseMatrix<double> &c) :
    SparseMatrix<double>(c), _comm(comm) {}

ParSparseMatrix<double>::ParSparseMatrix(MPI_Comm comm,
    INDEX* rows, INDEX* cols, double* vals, INDEX size,
    INDEX nRows, INDEX nCols, INDEX nRowsCRS)
  : SparseMatrix<double>(rows, cols, vals, size, nRows,
    nCols, nRowsCRS), _comm(comm) {}

//============================================================
void ParSparseMatrix<double>::MultMv(const Vector<double>& v,
    DenseVector<double>& c) const
{
  int numProcs = MPI_Wrappers::size(_comm);
#ifdef DISABLE_PAR_HEURISTICS
  // Use much more lenient heuristics to exercise parallel code
  if (numProcs == 1 ||  _size < 300) {
#else
  // These are simple heuristics to perform multiplication in serial if
  //   parallel will be slower. They were determined experimentally.
  if ( numProcs == 1 ||
      (_size < 50000  || _size > 10000000) ||
     ((_size < 150000 || _size > 5000000) && numProcs > 8) ||
     ((_size < 500000 || _size > 2500000) && numProcs > 16 ) ||
     (numProcs > 32)) {
#endif
    SparseMatrix<double>::MultMv(v, c);
    return;
  }


  SparseMatrix<double>::compress(*this);
  GCK(*this, v, this->nCols() != v.size(), "ParSparseMatrix * Vector")

  SparseMatrix<double> A_local;

  // Split the sparse matrix. partition() takes a ParSparMat, so we cast.
  partition(*static_cast<ParSparseMatrix<double>*>(&A_local));

  // actually do multiplication - end up with partial result vector
  // on each processor
#ifdef TIMING_ON
  timespec before, after;
//  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &before);
#endif
  DenseVector<double> c_local = A_local * v;
#ifdef TIMING_ON
//  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &after);
  cout << "P" << MPI_Wrappers::rank(MPI_COMM_WORLD) << " " << time_diff(before,after) << " mat.vec time\n";
  //LammpsInterface::instance()->all_print((after-before),"mat.vec time");
  barrier(MPI_COMM_WORLD);
#endif

  // destroy A_local intelligently
  static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

  // Add all the result vectors together on each processor.
#ifdef TIMING_ON
  barrier(MPI_COMM_WORLD);
//barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &before);
#endif
  allsum(_comm, c_local.ptr(), c.ptr(), c_local.size());
#ifdef TIMING_ON
//barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &after);
  cout << "P" << MPI_Wrappers::rank(MPI_COMM_WORLD) << " " << time_diff(before,after) << " allsum time\n";
  //LammpsInterface::instance()->print_msg_once((after-before),"allsum time");
#endif
}

DenseVector<double> ParSparseMatrix<double>::transMat(
    const Vector<double>& v) const {
  SparseMatrix<double>::compress(*this);
  GCK(*this, v, this->nRows() != v.size(), "ParSparseMatrix transpose * Vector")

  DenseVector<double> c(nCols(), true);

  SparseMatrix<double> A_local;
  partition(*static_cast<ParSparseMatrix<double>*>(&A_local));

  // actually do multiplication - end up with partial result vector
  // on each processor
  DenseVector<double> c_local = A_local.transMat(v);

  static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

  // Add all the result vectors together on each processor.
  allsum(_comm, c_local.ptr(), c.ptr(), c_local.size());

  return c;
}

void ParSparseMatrix<double>::MultAB(const Matrix<double>& B,
    DenseMatrix<double>& C) const {
  SparseMatrix<double>::compress(*this);
  GCK(*this, B, this->nCols() != B.nRows(), "ParSparseMatrix * Matrix")

  SparseMatrix<double> A_local;
  partition(*static_cast<ParSparseMatrix<double>*>(&A_local));

  // actually do multiplication - end up with partial result matrix
  // on each processor

#ifdef TIMING_ON
  timespec before, after;
  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &before);
#endif
  DenseMatrix<double> C_local = A_local * B;
#ifdef TIMING_ON
  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &after);
  cout << "P" << MPI_Wrappers::rank(MPI_COMM_WORLD) << " " << time_diff(after,before) << " mat.vec time\n";
  //LammpsInterface::instance()->all_print((after-before),"mat.vec time");
#endif

  static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

  // Add all the result vectors together on each processor.
#ifdef TIMING_ON
  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &before);
#endif
  allsum(_comm, C_local.ptr(), C.ptr(), C_local.size());
#ifdef TIMING_ON
  barrier(MPI_COMM_WORLD);
  clock_gettime(CLOCK_MONOTONIC, &after);
  cout << "P" << MPI_Wrappers::rank(MPI_COMM_WORLD) << " " << time_diff(after,before) << " allsum time\n";
  //LammpsInterface::instance()->print_msg_once((after-before),"allsum time");
#endif
}

DenseMatrix<double> ParSparseMatrix<double>::transMat(
    const DenseMatrix<double>& B) const {
  SparseMatrix<double>::compress(*this);
  GCK(*this, B, this->nRows() != B.nRows(), "ParSparseMatrix transpose * Matrix")

  DenseMatrix<double> C(nCols(), B.nCols(), true);

  SparseMatrix<double> A_local;
  partition(*static_cast<ParSparseMatrix<double>*>(&A_local));

  // actually do multiplication - end up with partial result matrix
  // on each processor
  DenseMatrix<double> C_local = A_local.transMat(B);

  static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

  // Add all the result vectors together on each processor.
  allsum(_comm, C_local.ptr(), C.ptr(), C_local.size());

  return C;
}

/*
 The two commented-out functions both need to return SparseMatrices. It's hard
 to combine sparse matrices between processors, so this has not yet been completed.

 void ParMultAB(const SparseMatrix<double> &B, SparseMatrix<double> &C) const
 {
 //SparseMatrix<T>::compress(*this);
 GCK(*this, B, this->nCols()!=B.nRows(), "ParSparseMatrix * SparseMatrix")

 ParSparseMatrix<double> A_local(this->_comm);
 this->partition(A_local);

 // actually do multiplication - end up with partial result matrix
 // on each processor

 SparseMatrix<double> C_local = ((SparseMatrix<double>)A_local) * B;

 // destroy newA intelligently
 static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

 // Add all the result vectors together on each processor.
 sumSparse(C_local, C);
 }*/

DenseMatrix<double> ParSparseMatrix<double>::transMat(
    const SparseMatrix<double>& B) const {
  SparseMatrix<double>::compress(*this);
  GCK(*this, B, this->nRows() != B.nRows(), "ParSparseMatrix transpose * SparseMatrix")

  DenseMatrix<double> C(nCols(), B.nCols(), true);

  SparseMatrix<double> A_local;
  partition(*static_cast<ParSparseMatrix<double>*>(&A_local));

  // actually do multiplication - end up with partial result matrix
  // on each processor
  DenseMatrix<double> C_local = A_local.transMat(B);

  static_cast<ParSparseMatrix<double>*>(&A_local)->finalize();

  // Add all the result vectors together on each processor.
  allsum(_comm, C_local.ptr(), C.ptr(), C_local.size());

  return C;
}

/*void ParMultAB(const DiagonalMatrix<double> &B, SparseMatrix<double> &C) const
 {
 //SparseMatrix<T>::compress(*this);
 GCK(*this, B, this->nCols()!=B.nRows(), "ParSparseMatrix * DiagonalMatrix")

 ParSparseMatrix<double> A_local(this->_comm);
 this->partition(A_local);

 // actually do multiplication - end up with partial result matrix
 // on each processor

 SparseMatrix<double> C_local = ((SparseMatrix<double>)A_local) * B;

 // destroy newA intelligently
 A_local._val = nullptr;
 A_local._ja = nullptr;

 // Add all the result vectors together on each processor.
 sumSparse(C_local, C);
 }*/

void ParSparseMatrix<double>::partition(
    ParSparseMatrix<double>& A_local) const {
  // create new sparse matrix on each processor, with same size and
  // a disjoint subset of A's elements.
  //
  // Ex: on two processors,
  //
  // |0 1 0|             |0 1 0|               |0 0 0|
  // |2 6 0| splits into |2 0 0| on proc 1 and |0 6 0| on proc 2
  // |0 0 3|             |0 0 0|               |0 0 3|
  //
  // We compute the subproducts individually on each processor, then
  // sum up all the vectors to get our final result.
  //

  // decide which elements will be in each submatrix
  INDEX startIndex = (MPI_Wrappers::rank(_comm) * size()) / MPI_Wrappers::size(_comm);
  INDEX endIndex = ((MPI_Wrappers::rank(_comm) + 1) * size()) / MPI_Wrappers::size(_comm);

  // update number of elements
  A_local._nRows = _nRows;
  A_local._nCols = _nCols;
  A_local._size = endIndex - startIndex;
  A_local._nRowsCRS = _nRowsCRS;
  // use pointer arithmetic to:
  // set newA's _val (to inside A's _val)
  A_local._val = _val + startIndex;
  // set newA's _ja (to inside A's _ja)
  A_local._ja = _ja + startIndex;
  // set newA's _ia (from scratch)
  A_local._ia = new INDEX[nRowsCRS() + 1];
  INDEX numRows = nRowsCRS();
  if (A_local._size > 0) {
    for (INDEX i = 0; i < numRows + 1; i++) {
      A_local._ia[i] = std::min(std::max((_ia[i] - startIndex), 0),
          endIndex - startIndex);
    }
  } else {
    A_local._nRowsCRS = 0;
  }
}

// Prepare an A_local matrix for deletion after it has been loaded with
//   data members from another matrix.
void ParSparseMatrix<double>::finalize() {
  _val = nullptr;
  _ja = nullptr;
}

void ParSparseMatrix<double>::operator=(const SparseMatrix<double> &source)
{
  copy(source);
}


/*void sumSparse(SparseMatrix<double> &C_local, SparseMatrix<double> &C)
 {
 }*/

}