30 template <
typename Derived1,
typename Derived2>
32 Eigen::PlainObjectBase<Derived1> &leftMatrix,
33 Eigen::PlainObjectBase<Derived2> &rightMatrix,
34 Eigen::PlainObjectBase<Derived2> &result,
37 bool cyclicComm =
true,
38 bool dotProductComputation =
true)
41 PRECICE_ASSERT(result.cols() == rightMatrix.cols(), result.cols(), rightMatrix.cols());
42 PRECICE_ASSERT(leftMatrix.cols() == rightMatrix.rows(), leftMatrix.cols(), rightMatrix.rows());
46 result.noalias() = leftMatrix * rightMatrix;
67 if (dotProductComputation)
87 template <
typename Derived1,
typename Derived2,
typename Derived3>
89 const Eigen::MatrixBase<Derived1> &leftMatrix,
90 const Eigen::MatrixBase<Derived2> &rightMatrix,
91 Eigen::PlainObjectBase<Derived3> &result,
96 PRECICE_ASSERT(leftMatrix.cols() == rightMatrix.rows(), leftMatrix.cols(), rightMatrix.rows());
100 Eigen::MatrixXd localResult(result.rows(), result.cols());
101 localResult.noalias() = leftMatrix * rightMatrix;
105 result = localResult;
115 template <
typename Derived1,
typename Derived2>
117 Eigen::PlainObjectBase<Derived1> &leftMatrix,
118 Eigen::PlainObjectBase<Derived2> &rightMatrix,
119 Eigen::PlainObjectBase<Derived2> &result,
141 int rows_rcv = (prevProc > 0) ? offsets[prevProc + 1] - offsets[prevProc] : offsets[1];
143 Eigen::MatrixXd leftMatrix_rcv(rows_rcv, q);
149 if (leftMatrix.size() > 0)
153 if (leftMatrix_rcv.size() > 0)
158 Eigen::MatrixXd diagBlock(leftMatrix.rows(), leftMatrix.rows());
159 diagBlock.noalias() = leftMatrix * rightMatrix;
163 PRECICE_ASSERT(result.cols() == diagBlock.cols(), result.cols(), diagBlock.cols());
164 result.block(off, 0, diagBlock.rows(), diagBlock.cols()) = diagBlock;
172 if (requestSend !=
nullptr)
174 if (requestRcv !=
nullptr)
178 Eigen::MatrixXd leftMatrix_copy(leftMatrix_rcv);
182 if (leftMatrix_copy.size() > 0)
191 int rows_rcv_nextCycle = (sourceProc_nextCycle > 0) ? offsets[sourceProc_nextCycle + 1] - offsets[sourceProc_nextCycle] : offsets[1];
192 rows_rcv = (sourceProc > 0) ? offsets[sourceProc + 1] - offsets[sourceProc] : offsets[1];
193 leftMatrix_rcv = Eigen::MatrixXd::Zero(rows_rcv_nextCycle, q);
197 if (leftMatrix_rcv.size() > 0)
201 if (requestSend !=
nullptr)
204 Eigen::MatrixXd block(rows_rcv, rightMatrix.cols());
205 block.noalias() = leftMatrix_copy * rightMatrix;
211 off = offsets[sourceProc];
212 PRECICE_ASSERT(result.cols() == block.cols(), result.cols(), block.cols());
213 result.block(off, 0, block.rows(), block.cols()) = block;
218 template <
typename Derived1,
typename Derived2>
220 Eigen::PlainObjectBase<Derived1> &leftMatrix,
221 Eigen::PlainObjectBase<Derived2> &rightMatrix,
222 Eigen::PlainObjectBase<Derived2> &result,
227 for (
int i = 0; i < leftMatrix.rows(); i++) {
234 while (i >= offsets[rank + 1])
237 Eigen::VectorXd lMRow = leftMatrix.row(i);
239 for (
int j = 0; j < r; j++) {
241 Eigen::VectorXd rMCol = rightMatrix.col(j);
247 local_row = i - offsets[rank];
248 result(local_row, j) = res_ij;
255 template <
typename Derived1,
typename Derived2>
257 Eigen::PlainObjectBase<Derived1> &leftMatrix,
258 Eigen::PlainObjectBase<Derived2> &rightMatrix,
259 Eigen::PlainObjectBase<Derived2> &result,
266 PRECICE_ASSERT(
static_cast<int>(leftMatrix.IsRowMajor) ==
static_cast<int>(rightMatrix.IsRowMajor),
267 static_cast<int>(leftMatrix.IsRowMajor),
static_cast<int>(rightMatrix.IsRowMajor));
271 Eigen::MatrixXd block = Eigen::MatrixXd::Zero(p, r);
272 block.noalias() = leftMatrix * rightMatrix;
279 Eigen::MatrixXd summarizedBlocks = Eigen::MatrixXd::Zero(p, r);
284 if (result.size() > 0)
291 result = summarizedBlocks.block(0, 0, offsets[1], r);
294 int off = offsets[secondaryRank];
295 int send_rows = offsets[secondaryRank + 1] - offsets[secondaryRank];
297 if (summarizedBlocks.block(off, 0, send_rows, r).size() > 0) {
300 Eigen::MatrixXd sendBlock = summarizedBlocks.block(off, 0, send_rows, r);
#define PRECICE_TRACE(...)
#define PRECICE_ASSERT(...)
void _multiplyNM_dotProduct(Eigen::PlainObjectBase< Derived1 > &leftMatrix, Eigen::PlainObjectBase< Derived2 > &rightMatrix, Eigen::PlainObjectBase< Derived2 > &result, const std::vector< int > &offsets, int p, int q, int r)
void closeCircularCommunication()
void initialize(const bool needCyclicComm)
Initializes the acceleration.
void multiply(const Eigen::MatrixBase< Derived1 > &leftMatrix, const Eigen::MatrixBase< Derived2 > &rightMatrix, Eigen::PlainObjectBase< Derived3 > &result, int p, int q, int r)
: Method computes the matrix-matrix/matrix-vector product of a (p x q) matrix that is distributed col...
void _multiply_cyclic(Eigen::PlainObjectBase< Derived1 > &leftMatrix, Eigen::PlainObjectBase< Derived2 > &rightMatrix, Eigen::PlainObjectBase< Derived2 > &result, const std::vector< int > &offsets, int p, int q, int r)
void multiply(Eigen::PlainObjectBase< Derived1 > &leftMatrix, Eigen::PlainObjectBase< Derived2 > &rightMatrix, Eigen::PlainObjectBase< Derived2 > &result, const std::vector< int > &offsets, int p, int q, int r, bool cyclicComm=true, bool dotProductComputation=true)
void _multiplyNM_block(Eigen::PlainObjectBase< Derived1 > &leftMatrix, Eigen::PlainObjectBase< Derived2 > &rightMatrix, Eigen::PlainObjectBase< Derived2 > &result, const std::vector< int > &offsets, int p, int q, int r)
Multiplies matrices based on a SAXPY-like block-wise computation with a rectangular result matrix of ...
com::PtrCommunication _cyclicCommLeft
Communication between neighboring ranks, backwards.
com::PtrCommunication _cyclicCommRight
Communication between neighboring ranks, forward.
~ParallelMatrixOperations()
void establishCircularCommunication()
This class provides a lightweight logger.
static void allreduceSum(precice::span< const double > sendData, precice::span< double > rcvData)
static int getSize()
Number of ranks. This includes ranks from both participants, e.g. minimal size is 2.
static Rank getRank()
Current rank.
static double dot(const Eigen::VectorXd &vec1, const Eigen::VectorXd &vec2)
static bool isPrimary()
True if this process is running the primary rank.
static auto allSecondaryRanks()
Returns an iterable range over salve ranks [1, _size)
static bool isParallel()
True if this process is running in parallel.
static bool isSecondary()
True if this process is running a secondary rank.
static com::PtrCommunication & getCommunication()
Intra-participant communication.
static void reduceSum(precice::span< const double > sendData, precice::span< double > rcvData)