Actual source code: mpimatmatmult.c
1: /*
2: Defines matrix-matrix product routines for pairs of MPIAIJ matrices
3: C = A * B
4: */
5: #include <../src/mat/impls/aij/seq/aij.h>
6: #include <../src/mat/utils/freespace.h>
7: #include <../src/mat/impls/aij/mpi/mpiaij.h>
8: #include <petscbt.h>
9: #include <../src/mat/impls/dense/mpi/mpidense.h>
10: #include <petsc/private/vecimpl.h>
11: #include <petsc/private/sfimpl.h>
13: #if defined(PETSC_HAVE_HYPRE)
14: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
15: #endif
17: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt_MPIAIJ_MPIAIJ(Mat C)
18: {
19: Mat_Product *product = C->product;
20: Mat B = product->B;
22: PetscFunctionBegin;
23: PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &product->B));
24: PetscCall(MatDestroy(&B));
25: PetscCall(MatProductSymbolic_AB_MPIAIJ_MPIAIJ(C));
26: PetscFunctionReturn(PETSC_SUCCESS);
27: }
29: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat C)
30: {
31: Mat_Product *product = C->product;
32: Mat A = product->A, B = product->B;
33: MatProductAlgorithm alg = product->alg;
34: PetscReal fill = product->fill;
35: PetscBool flg;
37: PetscFunctionBegin;
38: /* scalable */
39: PetscCall(PetscStrcmp(alg, "scalable", &flg));
40: if (flg) {
41: PetscCall(MatMatMultSymbolic_MPIAIJ_MPIAIJ(A, B, fill, C));
42: PetscFunctionReturn(PETSC_SUCCESS);
43: }
45: /* nonscalable */
46: PetscCall(PetscStrcmp(alg, "nonscalable", &flg));
47: if (flg) {
48: PetscCall(MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A, B, fill, C));
49: PetscFunctionReturn(PETSC_SUCCESS);
50: }
52: /* seqmpi */
53: PetscCall(PetscStrcmp(alg, "seqmpi", &flg));
54: if (flg) {
55: PetscCall(MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A, B, fill, C));
56: PetscFunctionReturn(PETSC_SUCCESS);
57: }
59: /* backend general code */
60: PetscCall(PetscStrcmp(alg, "backend", &flg));
61: if (flg) {
62: PetscCall(MatProductSymbolic_MPIAIJBACKEND(C));
63: PetscFunctionReturn(PETSC_SUCCESS);
64: }
66: #if defined(PETSC_HAVE_HYPRE)
67: PetscCall(PetscStrcmp(alg, "hypre", &flg));
68: if (flg) {
69: PetscCall(MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A, B, fill, C));
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
72: #endif
73: SETERRQ(PetscObjectComm((PetscObject)C), PETSC_ERR_SUP, "Mat Product Algorithm is not supported");
74: }
76: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *data)
77: {
78: Mat_APMPI *ptap = (Mat_APMPI *)data;
80: PetscFunctionBegin;
81: PetscCall(PetscFree2(ptap->startsj_s, ptap->startsj_r));
82: PetscCall(PetscFree(ptap->bufa));
83: PetscCall(MatDestroy(&ptap->P_loc));
84: PetscCall(MatDestroy(&ptap->P_oth));
85: PetscCall(MatDestroy(&ptap->Pt));
86: PetscCall(PetscFree(ptap->api));
87: PetscCall(PetscFree(ptap->apj));
88: PetscCall(PetscFree(ptap->apa));
89: PetscCall(PetscFree(ptap));
90: PetscFunctionReturn(PETSC_SUCCESS);
91: }
93: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A, Mat P, Mat C)
94: {
95: Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *c = (Mat_MPIAIJ *)C->data;
96: Mat_SeqAIJ *ad = (Mat_SeqAIJ *)a->A->data, *ao = (Mat_SeqAIJ *)a->B->data;
97: Mat_SeqAIJ *cd = (Mat_SeqAIJ *)c->A->data, *co = (Mat_SeqAIJ *)c->B->data;
98: PetscScalar *cda = cd->a, *coa = co->a;
99: Mat_SeqAIJ *p_loc, *p_oth;
100: PetscScalar *apa, *ca;
101: PetscInt cm = C->rmap->n;
102: Mat_APMPI *ptap;
103: PetscInt *api, *apj, *apJ, i, k;
104: PetscInt cstart = C->cmap->rstart;
105: PetscInt cdnz, conz, k0, k1;
106: const PetscScalar *dummy;
107: MPI_Comm comm;
108: PetscMPIInt size;
110: PetscFunctionBegin;
111: MatCheckProduct(C, 3);
112: ptap = (Mat_APMPI *)C->product->data;
113: PetscCheck(ptap, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtAP cannot be computed. Missing data");
114: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
115: PetscCallMPI(MPI_Comm_size(comm, &size));
116: PetscCheck(ptap->P_oth || size <= 1, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "AP cannot be reused. Do not call MatProductClear()");
118: /* flag CPU mask for C */
119: #if defined(PETSC_HAVE_DEVICE)
120: if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
121: if (c->A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->A->offloadmask = PETSC_OFFLOAD_CPU;
122: if (c->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->B->offloadmask = PETSC_OFFLOAD_CPU;
123: #endif
125: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
126: /* update numerical values of P_oth and P_loc */
127: PetscCall(MatGetBrowsOfAoCols_MPIAIJ(A, P, MAT_REUSE_MATRIX, &ptap->startsj_s, &ptap->startsj_r, &ptap->bufa, &ptap->P_oth));
128: PetscCall(MatMPIAIJGetLocalMat(P, MAT_REUSE_MATRIX, &ptap->P_loc));
130: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
131: /* get data from symbolic products */
132: p_loc = (Mat_SeqAIJ *)ptap->P_loc->data;
133: p_oth = NULL;
134: if (size > 1) p_oth = (Mat_SeqAIJ *)ptap->P_oth->data;
136: /* get apa for storing dense row A[i,:]*P */
137: apa = ptap->apa;
139: api = ptap->api;
140: apj = ptap->apj;
141: /* trigger copy to CPU */
142: PetscCall(MatSeqAIJGetArrayRead(a->A, &dummy));
143: PetscCall(MatSeqAIJRestoreArrayRead(a->A, &dummy));
144: PetscCall(MatSeqAIJGetArrayRead(a->B, &dummy));
145: PetscCall(MatSeqAIJRestoreArrayRead(a->B, &dummy));
146: for (i = 0; i < cm; i++) {
147: /* compute apa = A[i,:]*P */
148: AProw_nonscalable(i, ad, ao, p_loc, p_oth, apa);
150: /* set values in C */
151: apJ = apj + api[i];
152: cdnz = cd->i[i + 1] - cd->i[i];
153: conz = co->i[i + 1] - co->i[i];
155: /* 1st off-diagonal part of C */
156: ca = PetscSafePointerPlusOffset(coa, co->i[i]);
157: k = 0;
158: for (k0 = 0; k0 < conz; k0++) {
159: if (apJ[k] >= cstart) break;
160: ca[k0] = apa[apJ[k]];
161: apa[apJ[k++]] = 0.0;
162: }
164: /* diagonal part of C */
165: ca = PetscSafePointerPlusOffset(cda, cd->i[i]);
166: for (k1 = 0; k1 < cdnz; k1++) {
167: ca[k1] = apa[apJ[k]];
168: apa[apJ[k++]] = 0.0;
169: }
171: /* 2nd off-diagonal part of C */
172: ca = PetscSafePointerPlusOffset(coa, co->i[i]);
173: for (; k0 < conz; k0++) {
174: ca[k0] = apa[apJ[k]];
175: apa[apJ[k++]] = 0.0;
176: }
177: }
178: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
179: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A, Mat P, PetscReal fill, Mat C)
184: {
185: MPI_Comm comm;
186: PetscMPIInt size;
187: Mat_APMPI *ptap;
188: PetscFreeSpaceList free_space = NULL, current_space = NULL;
189: Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
190: Mat_SeqAIJ *ad = (Mat_SeqAIJ *)a->A->data, *ao = (Mat_SeqAIJ *)a->B->data, *p_loc, *p_oth;
191: PetscInt *pi_loc, *pj_loc, *pi_oth, *pj_oth, *dnz, *onz;
192: PetscInt *adi = ad->i, *adj = ad->j, *aoi = ao->i, *aoj = ao->j, rstart = A->rmap->rstart;
193: PetscInt *lnk, i, pnz, row, *api, *apj, *Jptr, apnz, nspacedouble = 0, j, nzi;
194: PetscInt am = A->rmap->n, pN = P->cmap->N, pn = P->cmap->n, pm = P->rmap->n;
195: PetscBT lnkbt;
196: PetscReal afill;
197: MatType mtype;
199: PetscFunctionBegin;
200: MatCheckProduct(C, 4);
201: PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty");
202: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
203: PetscCallMPI(MPI_Comm_size(comm, &size));
205: /* create struct Mat_APMPI and attached it to C later */
206: PetscCall(PetscNew(&ptap));
208: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
209: PetscCall(MatGetBrowsOfAoCols_MPIAIJ(A, P, MAT_INITIAL_MATRIX, &ptap->startsj_s, &ptap->startsj_r, &ptap->bufa, &ptap->P_oth));
211: /* get P_loc by taking all local rows of P */
212: PetscCall(MatMPIAIJGetLocalMat(P, MAT_INITIAL_MATRIX, &ptap->P_loc));
214: p_loc = (Mat_SeqAIJ *)ptap->P_loc->data;
215: pi_loc = p_loc->i;
216: pj_loc = p_loc->j;
217: if (size > 1) {
218: p_oth = (Mat_SeqAIJ *)ptap->P_oth->data;
219: pi_oth = p_oth->i;
220: pj_oth = p_oth->j;
221: } else {
222: p_oth = NULL;
223: pi_oth = NULL;
224: pj_oth = NULL;
225: }
227: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
228: PetscCall(PetscMalloc1(am + 2, &api));
229: ptap->api = api;
230: api[0] = 0;
232: /* create and initialize a linked list */
233: PetscCall(PetscLLCondensedCreate(pN, pN, &lnk, &lnkbt));
235: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
236: PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(adi[am], PetscIntSumTruncate(aoi[am], pi_loc[pm]))), &free_space));
237: current_space = free_space;
239: MatPreallocateBegin(comm, am, pn, dnz, onz);
240: for (i = 0; i < am; i++) {
241: /* diagonal portion of A */
242: nzi = adi[i + 1] - adi[i];
243: for (j = 0; j < nzi; j++) {
244: row = *adj++;
245: pnz = pi_loc[row + 1] - pi_loc[row];
246: Jptr = pj_loc + pi_loc[row];
247: /* add non-zero cols of P into the sorted linked list lnk */
248: PetscCall(PetscLLCondensedAddSorted(pnz, Jptr, lnk, lnkbt));
249: }
250: /* off-diagonal portion of A */
251: nzi = aoi[i + 1] - aoi[i];
252: for (j = 0; j < nzi; j++) {
253: row = *aoj++;
254: pnz = pi_oth[row + 1] - pi_oth[row];
255: Jptr = pj_oth + pi_oth[row];
256: PetscCall(PetscLLCondensedAddSorted(pnz, Jptr, lnk, lnkbt));
257: }
258: /* add possible missing diagonal entry */
259: if (C->force_diagonals) {
260: j = i + rstart; /* column index */
261: PetscCall(PetscLLCondensedAddSorted(1, &j, lnk, lnkbt));
262: }
264: apnz = lnk[0];
265: api[i + 1] = api[i] + apnz;
267: /* if free space is not available, double the total space in the list */
268: if (current_space->local_remaining < apnz) {
269: PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(apnz, current_space->total_array_size), ¤t_space));
270: nspacedouble++;
271: }
273: /* Copy data into free space, then initialize lnk */
274: PetscCall(PetscLLCondensedClean(pN, apnz, current_space->array, lnk, lnkbt));
275: PetscCall(MatPreallocateSet(i + rstart, apnz, current_space->array, dnz, onz));
277: current_space->array += apnz;
278: current_space->local_used += apnz;
279: current_space->local_remaining -= apnz;
280: }
282: /* Allocate space for apj, initialize apj, and */
283: /* destroy list of free space and other temporary array(s) */
284: PetscCall(PetscMalloc1(api[am] + 1, &ptap->apj));
285: apj = ptap->apj;
286: PetscCall(PetscFreeSpaceContiguous(&free_space, ptap->apj));
287: PetscCall(PetscLLDestroy(lnk, lnkbt));
289: /* malloc apa to store dense row A[i,:]*P */
290: PetscCall(PetscCalloc1(pN, &ptap->apa));
292: /* set and assemble symbolic parallel matrix C */
293: PetscCall(MatSetSizes(C, am, pn, PETSC_DETERMINE, PETSC_DETERMINE));
294: PetscCall(MatSetBlockSizesFromMats(C, A, P));
296: PetscCall(MatGetType(A, &mtype));
297: PetscCall(MatSetType(C, mtype));
298: PetscCall(MatMPIAIJSetPreallocation(C, 0, dnz, 0, onz));
299: MatPreallocateEnd(dnz, onz);
301: PetscCall(MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api));
302: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
303: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
304: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
306: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
307: C->ops->productnumeric = MatProductNumeric_AB;
309: /* attach the supporting struct to C for reuse */
310: C->product->data = ptap;
311: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
313: /* set MatInfo */
314: afill = (PetscReal)api[am] / (adi[am] + aoi[am] + pi_loc[pm] + 1) + 1.e-5;
315: if (afill < 1.0) afill = 1.0;
316: C->info.mallocs = nspacedouble;
317: C->info.fill_ratio_given = fill;
318: C->info.fill_ratio_needed = afill;
320: #if defined(PETSC_USE_INFO)
321: if (api[am]) {
322: PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
323: PetscCall(PetscInfo(C, "Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n", (double)afill));
324: } else {
325: PetscCall(PetscInfo(C, "Empty matrix product\n"));
326: }
327: #endif
328: PetscFunctionReturn(PETSC_SUCCESS);
329: }
331: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat, Mat, PetscReal, Mat);
332: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat, Mat, Mat);
334: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AB(Mat C)
335: {
336: Mat_Product *product = C->product;
337: Mat A = product->A, B = product->B;
339: PetscFunctionBegin;
340: if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend)
341: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, (%" PetscInt_FMT ", %" PetscInt_FMT ") != (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->cmap->rstart, A->cmap->rend, B->rmap->rstart, B->rmap->rend);
343: C->ops->matmultsymbolic = MatMatMultSymbolic_MPIAIJ_MPIDense;
344: C->ops->productsymbolic = MatProductSymbolic_AB;
345: PetscFunctionReturn(PETSC_SUCCESS);
346: }
348: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(Mat C)
349: {
350: Mat_Product *product = C->product;
351: Mat A = product->A, B = product->B;
353: PetscFunctionBegin;
354: if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
355: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, (%" PetscInt_FMT ", %" PetscInt_FMT ") != (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->rstart, A->rmap->rend, B->rmap->rstart, B->rmap->rend);
357: C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense;
358: C->ops->productsymbolic = MatProductSymbolic_AtB;
359: PetscFunctionReturn(PETSC_SUCCESS);
360: }
362: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat C)
363: {
364: Mat_Product *product = C->product;
366: PetscFunctionBegin;
367: switch (product->type) {
368: case MATPRODUCT_AB:
369: PetscCall(MatProductSetFromOptions_MPIAIJ_MPIDense_AB(C));
370: break;
371: case MATPRODUCT_AtB:
372: PetscCall(MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(C));
373: break;
374: default:
375: break;
376: }
377: PetscFunctionReturn(PETSC_SUCCESS);
378: }
380: typedef struct {
381: Mat workB, workB1;
382: MPI_Request *rwaits, *swaits;
383: PetscInt nsends, nrecvs;
384: MPI_Datatype *stype, *rtype;
385: PetscInt blda;
386: } MPIAIJ_MPIDense;
388: static PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
389: {
390: MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense *)ctx;
391: PetscInt i;
393: PetscFunctionBegin;
394: PetscCall(MatDestroy(&contents->workB));
395: PetscCall(MatDestroy(&contents->workB1));
396: for (i = 0; i < contents->nsends; i++) PetscCallMPI(MPI_Type_free(&contents->stype[i]));
397: for (i = 0; i < contents->nrecvs; i++) PetscCallMPI(MPI_Type_free(&contents->rtype[i]));
398: PetscCall(PetscFree4(contents->stype, contents->rtype, contents->rwaits, contents->swaits));
399: PetscCall(PetscFree(contents));
400: PetscFunctionReturn(PETSC_SUCCESS);
401: }
403: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A, Mat B, PetscReal fill, Mat C)
404: {
405: Mat_MPIAIJ *aij = (Mat_MPIAIJ *)A->data;
406: PetscInt nz = aij->B->cmap->n, nsends, nrecvs, i, nrows_to, j, blda, m, M, n, N;
407: MPIAIJ_MPIDense *contents;
408: VecScatter ctx = aij->Mvctx;
409: PetscInt Am = A->rmap->n, Bm = B->rmap->n, BN = B->cmap->N, Bbn, Bbn1, bs, nrows_from, numBb;
410: MPI_Comm comm;
411: MPI_Datatype type1, *stype, *rtype;
412: const PetscInt *sindices, *sstarts, *rstarts;
413: PetscMPIInt *disp;
414: PetscBool cisdense;
416: PetscFunctionBegin;
417: MatCheckProduct(C, 4);
418: PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
419: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
420: PetscCall(PetscObjectBaseTypeCompare((PetscObject)C, MATMPIDENSE, &cisdense));
421: if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
422: PetscCall(MatGetLocalSize(C, &m, &n));
423: PetscCall(MatGetSize(C, &M, &N));
424: if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) PetscCall(MatSetSizes(C, Am, B->cmap->n, A->rmap->N, BN));
425: PetscCall(MatSetBlockSizesFromMats(C, A, B));
426: PetscCall(MatSetUp(C));
427: PetscCall(MatDenseGetLDA(B, &blda));
428: PetscCall(PetscNew(&contents));
430: PetscCall(VecScatterGetRemote_Private(ctx, PETSC_TRUE /*send*/, &nsends, &sstarts, &sindices, NULL, NULL));
431: PetscCall(VecScatterGetRemoteOrdered_Private(ctx, PETSC_FALSE /*recv*/, &nrecvs, &rstarts, NULL, NULL, NULL));
433: /* Create column block of B and C for memory scalability when BN is too large */
434: /* Estimate Bbn, column size of Bb */
435: if (nz) {
436: Bbn1 = 2 * Am * BN / nz;
437: if (!Bbn1) Bbn1 = 1;
438: } else Bbn1 = BN;
440: bs = PetscAbs(B->cmap->bs);
441: Bbn1 = Bbn1 / bs * bs; /* Bbn1 is a multiple of bs */
442: if (Bbn1 > BN) Bbn1 = BN;
443: PetscCall(MPIU_Allreduce(&Bbn1, &Bbn, 1, MPIU_INT, MPI_MAX, comm));
445: /* Enable runtime option for Bbn */
446: PetscOptionsBegin(comm, ((PetscObject)C)->prefix, "MatMatMult", "Mat");
447: PetscCall(PetscOptionsInt("-matmatmult_Bbn", "Number of columns in Bb", "MatMatMult", Bbn, &Bbn, NULL));
448: PetscOptionsEnd();
449: Bbn = PetscMin(Bbn, BN);
451: if (Bbn > 0 && Bbn < BN) {
452: numBb = BN / Bbn;
453: Bbn1 = BN - numBb * Bbn;
454: } else numBb = 0;
456: if (numBb) {
457: PetscCall(PetscInfo(C, "use Bb, BN=%" PetscInt_FMT ", Bbn=%" PetscInt_FMT "; numBb=%" PetscInt_FMT "\n", BN, Bbn, numBb));
458: if (Bbn1) { /* Create workB1 for the remaining columns */
459: PetscCall(PetscInfo(C, "use Bb1, BN=%" PetscInt_FMT ", Bbn1=%" PetscInt_FMT "\n", BN, Bbn1));
460: /* Create work matrix used to store off processor rows of B needed for local product */
461: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nz, Bbn1, NULL, &contents->workB1));
462: } else contents->workB1 = NULL;
463: }
465: /* Create work matrix used to store off processor rows of B needed for local product */
466: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nz, Bbn, NULL, &contents->workB));
468: /* Use MPI derived data type to reduce memory required by the send/recv buffers */
469: PetscCall(PetscMalloc4(nsends, &stype, nrecvs, &rtype, nrecvs, &contents->rwaits, nsends, &contents->swaits));
470: contents->stype = stype;
471: contents->nsends = nsends;
473: contents->rtype = rtype;
474: contents->nrecvs = nrecvs;
475: contents->blda = blda;
477: PetscCall(PetscMalloc1(Bm + 1, &disp));
478: for (i = 0; i < nsends; i++) {
479: nrows_to = sstarts[i + 1] - sstarts[i];
480: for (j = 0; j < nrows_to; j++) disp[j] = sindices[sstarts[i] + j]; /* rowB to be sent */
481: PetscCallMPI(MPI_Type_create_indexed_block(nrows_to, 1, disp, MPIU_SCALAR, &type1));
482: PetscCallMPI(MPI_Type_create_resized(type1, 0, blda * sizeof(PetscScalar), &stype[i]));
483: PetscCallMPI(MPI_Type_commit(&stype[i]));
484: PetscCallMPI(MPI_Type_free(&type1));
485: }
487: for (i = 0; i < nrecvs; i++) {
488: /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */
489: nrows_from = rstarts[i + 1] - rstarts[i];
490: disp[0] = 0;
491: PetscCallMPI(MPI_Type_create_indexed_block(1, nrows_from, disp, MPIU_SCALAR, &type1));
492: PetscCallMPI(MPI_Type_create_resized(type1, 0, nz * sizeof(PetscScalar), &rtype[i]));
493: PetscCallMPI(MPI_Type_commit(&rtype[i]));
494: PetscCallMPI(MPI_Type_free(&type1));
495: }
497: PetscCall(PetscFree(disp));
498: PetscCall(VecScatterRestoreRemote_Private(ctx, PETSC_TRUE /*send*/, &nsends, &sstarts, &sindices, NULL, NULL));
499: PetscCall(VecScatterRestoreRemoteOrdered_Private(ctx, PETSC_FALSE /*recv*/, &nrecvs, &rstarts, NULL, NULL, NULL));
500: PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
501: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
502: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
504: C->product->data = contents;
505: C->product->destroy = MatMPIAIJ_MPIDenseDestroy;
506: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
507: PetscFunctionReturn(PETSC_SUCCESS);
508: }
510: PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat, Mat, Mat, const PetscBool);
512: /*
513: Performs an efficient scatter on the rows of B needed by this process; this is
514: a modification of the VecScatterBegin_() routines.
516: Input: If Bbidx = 0, uses B = Bb, else B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense()
517: */
519: static PetscErrorCode MatMPIDenseScatter(Mat A, Mat B, PetscInt Bbidx, Mat C, Mat *outworkB)
520: {
521: Mat_MPIAIJ *aij = (Mat_MPIAIJ *)A->data;
522: const PetscScalar *b;
523: PetscScalar *rvalues;
524: VecScatter ctx = aij->Mvctx;
525: const PetscInt *sindices, *sstarts, *rstarts;
526: const PetscMPIInt *sprocs, *rprocs;
527: PetscInt i, nsends, nrecvs;
528: MPI_Request *swaits, *rwaits;
529: MPI_Comm comm;
530: PetscMPIInt tag = ((PetscObject)ctx)->tag, ncols = B->cmap->N, nrows = aij->B->cmap->n, nsends_mpi, nrecvs_mpi;
531: MPIAIJ_MPIDense *contents;
532: Mat workB;
533: MPI_Datatype *stype, *rtype;
534: PetscInt blda;
536: PetscFunctionBegin;
537: MatCheckProduct(C, 4);
538: PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
539: contents = (MPIAIJ_MPIDense *)C->product->data;
540: PetscCall(VecScatterGetRemote_Private(ctx, PETSC_TRUE /*send*/, &nsends, &sstarts, &sindices, &sprocs, NULL /*bs*/));
541: PetscCall(VecScatterGetRemoteOrdered_Private(ctx, PETSC_FALSE /*recv*/, &nrecvs, &rstarts, NULL, &rprocs, NULL /*bs*/));
542: PetscCall(PetscMPIIntCast(nsends, &nsends_mpi));
543: PetscCall(PetscMPIIntCast(nrecvs, &nrecvs_mpi));
544: if (Bbidx == 0) workB = *outworkB = contents->workB;
545: else workB = *outworkB = contents->workB1;
546: PetscCheck(nrows == workB->rmap->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of rows of workB %" PetscInt_FMT " not equal to columns of aij->B %d", workB->cmap->n, nrows);
547: swaits = contents->swaits;
548: rwaits = contents->rwaits;
550: PetscCall(MatDenseGetArrayRead(B, &b));
551: PetscCall(MatDenseGetLDA(B, &blda));
552: PetscCheck(blda == contents->blda, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot reuse an input matrix with lda %" PetscInt_FMT " != %" PetscInt_FMT, blda, contents->blda);
553: PetscCall(MatDenseGetArray(workB, &rvalues));
555: /* Post recv, use MPI derived data type to save memory */
556: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
557: rtype = contents->rtype;
558: for (i = 0; i < nrecvs; i++) PetscCallMPI(MPI_Irecv(rvalues + (rstarts[i] - rstarts[0]), ncols, rtype[i], rprocs[i], tag, comm, rwaits + i));
560: stype = contents->stype;
561: for (i = 0; i < nsends; i++) PetscCallMPI(MPI_Isend(b, ncols, stype[i], sprocs[i], tag, comm, swaits + i));
563: if (nrecvs) PetscCallMPI(MPI_Waitall(nrecvs_mpi, rwaits, MPI_STATUSES_IGNORE));
564: if (nsends) PetscCallMPI(MPI_Waitall(nsends_mpi, swaits, MPI_STATUSES_IGNORE));
566: PetscCall(VecScatterRestoreRemote_Private(ctx, PETSC_TRUE /*send*/, &nsends, &sstarts, &sindices, &sprocs, NULL));
567: PetscCall(VecScatterRestoreRemoteOrdered_Private(ctx, PETSC_FALSE /*recv*/, &nrecvs, &rstarts, NULL, &rprocs, NULL));
568: PetscCall(MatDenseRestoreArrayRead(B, &b));
569: PetscCall(MatDenseRestoreArray(workB, &rvalues));
570: PetscFunctionReturn(PETSC_SUCCESS);
571: }
573: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A, Mat B, Mat C)
574: {
575: Mat_MPIAIJ *aij = (Mat_MPIAIJ *)A->data;
576: Mat_MPIDense *bdense = (Mat_MPIDense *)B->data;
577: Mat_MPIDense *cdense = (Mat_MPIDense *)C->data;
578: Mat workB;
579: MPIAIJ_MPIDense *contents;
581: PetscFunctionBegin;
582: MatCheckProduct(C, 3);
583: PetscCheck(C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
584: contents = (MPIAIJ_MPIDense *)C->product->data;
585: /* diagonal block of A times all local rows of B */
586: /* TODO: this calls a symbolic multiplication every time, which could be avoided */
587: PetscCall(MatMatMult(aij->A, bdense->A, MAT_REUSE_MATRIX, PETSC_DEFAULT, &cdense->A));
588: if (contents->workB->cmap->n == B->cmap->N) {
589: /* get off processor parts of B needed to complete C=A*B */
590: PetscCall(MatMPIDenseScatter(A, B, 0, C, &workB));
592: /* off-diagonal block of A times nonlocal rows of B */
593: PetscCall(MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B, workB, cdense->A, PETSC_TRUE));
594: } else {
595: Mat Bb, Cb;
596: PetscInt BN = B->cmap->N, n = contents->workB->cmap->n, i;
597: PetscBool ccpu;
599: PetscCheck(n > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Column block size %" PetscInt_FMT " must be positive", n);
600: /* Prevent from unneeded copies back and forth from the GPU
601: when getting and restoring the submatrix
602: We need a proper GPU code for AIJ * dense in parallel */
603: PetscCall(MatBoundToCPU(C, &ccpu));
604: PetscCall(MatBindToCPU(C, PETSC_TRUE));
605: for (i = 0; i < BN; i += n) {
606: PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, i, PetscMin(i + n, BN), &Bb));
607: PetscCall(MatDenseGetSubMatrix(C, PETSC_DECIDE, PETSC_DECIDE, i, PetscMin(i + n, BN), &Cb));
609: /* get off processor parts of B needed to complete C=A*B */
610: PetscCall(MatMPIDenseScatter(A, Bb, (i + n) > BN, C, &workB));
612: /* off-diagonal block of A times nonlocal rows of B */
613: cdense = (Mat_MPIDense *)Cb->data;
614: PetscCall(MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B, workB, cdense->A, PETSC_TRUE));
615: PetscCall(MatDenseRestoreSubMatrix(B, &Bb));
616: PetscCall(MatDenseRestoreSubMatrix(C, &Cb));
617: }
618: PetscCall(MatBindToCPU(C, ccpu));
619: }
620: PetscFunctionReturn(PETSC_SUCCESS);
621: }
623: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A, Mat P, Mat C)
624: {
625: Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data, *c = (Mat_MPIAIJ *)C->data;
626: Mat_SeqAIJ *ad = (Mat_SeqAIJ *)a->A->data, *ao = (Mat_SeqAIJ *)a->B->data;
627: Mat_SeqAIJ *cd = (Mat_SeqAIJ *)c->A->data, *co = (Mat_SeqAIJ *)c->B->data;
628: PetscInt *adi = ad->i, *adj, *aoi = ao->i, *aoj;
629: PetscScalar *ada, *aoa, *cda = cd->a, *coa = co->a;
630: Mat_SeqAIJ *p_loc, *p_oth;
631: PetscInt *pi_loc, *pj_loc, *pi_oth, *pj_oth, *pj;
632: PetscScalar *pa_loc, *pa_oth, *pa, valtmp, *ca;
633: PetscInt cm = C->rmap->n, anz, pnz;
634: Mat_APMPI *ptap;
635: PetscScalar *apa_sparse;
636: const PetscScalar *dummy;
637: PetscInt *api, *apj, *apJ, i, j, k, row;
638: PetscInt cstart = C->cmap->rstart;
639: PetscInt cdnz, conz, k0, k1, nextp;
640: MPI_Comm comm;
641: PetscMPIInt size;
643: PetscFunctionBegin;
644: MatCheckProduct(C, 3);
645: ptap = (Mat_APMPI *)C->product->data;
646: PetscCheck(ptap, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtAP cannot be computed. Missing data");
647: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
648: PetscCallMPI(MPI_Comm_size(comm, &size));
649: PetscCheck(ptap->P_oth || size <= 1, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "AP cannot be reused. Do not call MatProductClear()");
651: /* flag CPU mask for C */
652: #if defined(PETSC_HAVE_DEVICE)
653: if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
654: if (c->A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->A->offloadmask = PETSC_OFFLOAD_CPU;
655: if (c->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->B->offloadmask = PETSC_OFFLOAD_CPU;
656: #endif
657: apa_sparse = ptap->apa;
659: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
660: /* update numerical values of P_oth and P_loc */
661: PetscCall(MatGetBrowsOfAoCols_MPIAIJ(A, P, MAT_REUSE_MATRIX, &ptap->startsj_s, &ptap->startsj_r, &ptap->bufa, &ptap->P_oth));
662: PetscCall(MatMPIAIJGetLocalMat(P, MAT_REUSE_MATRIX, &ptap->P_loc));
664: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
665: /* get data from symbolic products */
666: p_loc = (Mat_SeqAIJ *)ptap->P_loc->data;
667: pi_loc = p_loc->i;
668: pj_loc = p_loc->j;
669: pa_loc = p_loc->a;
670: if (size > 1) {
671: p_oth = (Mat_SeqAIJ *)ptap->P_oth->data;
672: pi_oth = p_oth->i;
673: pj_oth = p_oth->j;
674: pa_oth = p_oth->a;
675: } else {
676: p_oth = NULL;
677: pi_oth = NULL;
678: pj_oth = NULL;
679: pa_oth = NULL;
680: }
682: /* trigger copy to CPU */
683: PetscCall(MatSeqAIJGetArrayRead(a->A, &dummy));
684: PetscCall(MatSeqAIJRestoreArrayRead(a->A, &dummy));
685: PetscCall(MatSeqAIJGetArrayRead(a->B, &dummy));
686: PetscCall(MatSeqAIJRestoreArrayRead(a->B, &dummy));
687: api = ptap->api;
688: apj = ptap->apj;
689: for (i = 0; i < cm; i++) {
690: apJ = apj + api[i];
692: /* diagonal portion of A */
693: anz = adi[i + 1] - adi[i];
694: adj = ad->j + adi[i];
695: ada = ad->a + adi[i];
696: for (j = 0; j < anz; j++) {
697: row = adj[j];
698: pnz = pi_loc[row + 1] - pi_loc[row];
699: pj = pj_loc + pi_loc[row];
700: pa = pa_loc + pi_loc[row];
701: /* perform sparse axpy */
702: valtmp = ada[j];
703: nextp = 0;
704: for (k = 0; nextp < pnz; k++) {
705: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
706: apa_sparse[k] += valtmp * pa[nextp++];
707: }
708: }
709: PetscCall(PetscLogFlops(2.0 * pnz));
710: }
712: /* off-diagonal portion of A */
713: anz = aoi[i + 1] - aoi[i];
714: aoj = PetscSafePointerPlusOffset(ao->j, aoi[i]);
715: aoa = PetscSafePointerPlusOffset(ao->a, aoi[i]);
716: for (j = 0; j < anz; j++) {
717: row = aoj[j];
718: pnz = pi_oth[row + 1] - pi_oth[row];
719: pj = pj_oth + pi_oth[row];
720: pa = pa_oth + pi_oth[row];
721: /* perform sparse axpy */
722: valtmp = aoa[j];
723: nextp = 0;
724: for (k = 0; nextp < pnz; k++) {
725: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
726: apa_sparse[k] += valtmp * pa[nextp++];
727: }
728: }
729: PetscCall(PetscLogFlops(2.0 * pnz));
730: }
732: /* set values in C */
733: cdnz = cd->i[i + 1] - cd->i[i];
734: conz = co->i[i + 1] - co->i[i];
736: /* 1st off-diagonal part of C */
737: ca = PetscSafePointerPlusOffset(coa, co->i[i]);
738: k = 0;
739: for (k0 = 0; k0 < conz; k0++) {
740: if (apJ[k] >= cstart) break;
741: ca[k0] = apa_sparse[k];
742: apa_sparse[k] = 0.0;
743: k++;
744: }
746: /* diagonal part of C */
747: ca = cda + cd->i[i];
748: for (k1 = 0; k1 < cdnz; k1++) {
749: ca[k1] = apa_sparse[k];
750: apa_sparse[k] = 0.0;
751: k++;
752: }
754: /* 2nd off-diagonal part of C */
755: ca = PetscSafePointerPlusOffset(coa, co->i[i]);
756: for (; k0 < conz; k0++) {
757: ca[k0] = apa_sparse[k];
758: apa_sparse[k] = 0.0;
759: k++;
760: }
761: }
762: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
763: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
764: PetscFunctionReturn(PETSC_SUCCESS);
765: }
767: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
768: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A, Mat P, PetscReal fill, Mat C)
769: {
770: MPI_Comm comm;
771: PetscMPIInt size;
772: Mat_APMPI *ptap;
773: PetscFreeSpaceList free_space = NULL, current_space = NULL;
774: Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
775: Mat_SeqAIJ *ad = (Mat_SeqAIJ *)a->A->data, *ao = (Mat_SeqAIJ *)a->B->data, *p_loc, *p_oth;
776: PetscInt *pi_loc, *pj_loc, *pi_oth, *pj_oth, *dnz, *onz;
777: PetscInt *adi = ad->i, *adj = ad->j, *aoi = ao->i, *aoj = ao->j, rstart = A->rmap->rstart;
778: PetscInt i, pnz, row, *api, *apj, *Jptr, apnz, nspacedouble = 0, j, nzi, *lnk, apnz_max = 1;
779: PetscInt am = A->rmap->n, pn = P->cmap->n, pm = P->rmap->n, lsize = pn + 20;
780: PetscReal afill;
781: MatType mtype;
783: PetscFunctionBegin;
784: MatCheckProduct(C, 4);
785: PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty");
786: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
787: PetscCallMPI(MPI_Comm_size(comm, &size));
789: /* create struct Mat_APMPI and attached it to C later */
790: PetscCall(PetscNew(&ptap));
792: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
793: PetscCall(MatGetBrowsOfAoCols_MPIAIJ(A, P, MAT_INITIAL_MATRIX, &ptap->startsj_s, &ptap->startsj_r, &ptap->bufa, &ptap->P_oth));
795: /* get P_loc by taking all local rows of P */
796: PetscCall(MatMPIAIJGetLocalMat(P, MAT_INITIAL_MATRIX, &ptap->P_loc));
798: p_loc = (Mat_SeqAIJ *)ptap->P_loc->data;
799: pi_loc = p_loc->i;
800: pj_loc = p_loc->j;
801: if (size > 1) {
802: p_oth = (Mat_SeqAIJ *)ptap->P_oth->data;
803: pi_oth = p_oth->i;
804: pj_oth = p_oth->j;
805: } else {
806: p_oth = NULL;
807: pi_oth = NULL;
808: pj_oth = NULL;
809: }
811: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
812: PetscCall(PetscMalloc1(am + 2, &api));
813: ptap->api = api;
814: api[0] = 0;
816: PetscCall(PetscLLCondensedCreate_Scalable(lsize, &lnk));
818: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
819: PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(adi[am], PetscIntSumTruncate(aoi[am], pi_loc[pm]))), &free_space));
820: current_space = free_space;
821: MatPreallocateBegin(comm, am, pn, dnz, onz);
822: for (i = 0; i < am; i++) {
823: /* diagonal portion of A */
824: nzi = adi[i + 1] - adi[i];
825: for (j = 0; j < nzi; j++) {
826: row = *adj++;
827: pnz = pi_loc[row + 1] - pi_loc[row];
828: Jptr = pj_loc + pi_loc[row];
829: /* Expand list if it is not long enough */
830: if (pnz + apnz_max > lsize) {
831: lsize = pnz + apnz_max;
832: PetscCall(PetscLLCondensedExpand_Scalable(lsize, &lnk));
833: }
834: /* add non-zero cols of P into the sorted linked list lnk */
835: PetscCall(PetscLLCondensedAddSorted_Scalable(pnz, Jptr, lnk));
836: apnz = *lnk; /* The first element in the list is the number of items in the list */
837: api[i + 1] = api[i] + apnz;
838: if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
839: }
840: /* off-diagonal portion of A */
841: nzi = aoi[i + 1] - aoi[i];
842: for (j = 0; j < nzi; j++) {
843: row = *aoj++;
844: pnz = pi_oth[row + 1] - pi_oth[row];
845: Jptr = pj_oth + pi_oth[row];
846: /* Expand list if it is not long enough */
847: if (pnz + apnz_max > lsize) {
848: lsize = pnz + apnz_max;
849: PetscCall(PetscLLCondensedExpand_Scalable(lsize, &lnk));
850: }
851: /* add non-zero cols of P into the sorted linked list lnk */
852: PetscCall(PetscLLCondensedAddSorted_Scalable(pnz, Jptr, lnk));
853: apnz = *lnk; /* The first element in the list is the number of items in the list */
854: api[i + 1] = api[i] + apnz;
855: if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
856: }
858: /* add missing diagonal entry */
859: if (C->force_diagonals) {
860: j = i + rstart; /* column index */
861: PetscCall(PetscLLCondensedAddSorted_Scalable(1, &j, lnk));
862: }
864: apnz = *lnk;
865: api[i + 1] = api[i] + apnz;
866: if (apnz > apnz_max) apnz_max = apnz;
868: /* if free space is not available, double the total space in the list */
869: if (current_space->local_remaining < apnz) {
870: PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(apnz, current_space->total_array_size), ¤t_space));
871: nspacedouble++;
872: }
874: /* Copy data into free space, then initialize lnk */
875: PetscCall(PetscLLCondensedClean_Scalable(apnz, current_space->array, lnk));
876: PetscCall(MatPreallocateSet(i + rstart, apnz, current_space->array, dnz, onz));
878: current_space->array += apnz;
879: current_space->local_used += apnz;
880: current_space->local_remaining -= apnz;
881: }
883: /* Allocate space for apj, initialize apj, and */
884: /* destroy list of free space and other temporary array(s) */
885: PetscCall(PetscMalloc1(api[am] + 1, &ptap->apj));
886: apj = ptap->apj;
887: PetscCall(PetscFreeSpaceContiguous(&free_space, ptap->apj));
888: PetscCall(PetscLLCondensedDestroy_Scalable(lnk));
890: /* create and assemble symbolic parallel matrix C */
891: PetscCall(MatSetSizes(C, am, pn, PETSC_DETERMINE, PETSC_DETERMINE));
892: PetscCall(MatSetBlockSizesFromMats(C, A, P));
893: PetscCall(MatGetType(A, &mtype));
894: PetscCall(MatSetType(C, mtype));
895: PetscCall(MatMPIAIJSetPreallocation(C, 0, dnz, 0, onz));
896: MatPreallocateEnd(dnz, onz);
898: /* malloc apa for assembly C */
899: PetscCall(PetscCalloc1(apnz_max, &ptap->apa));
901: PetscCall(MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api));
902: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
903: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
904: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
906: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
907: C->ops->productnumeric = MatProductNumeric_AB;
909: /* attach the supporting struct to C for reuse */
910: C->product->data = ptap;
911: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
913: /* set MatInfo */
914: afill = (PetscReal)api[am] / (adi[am] + aoi[am] + pi_loc[pm] + 1) + 1.e-5;
915: if (afill < 1.0) afill = 1.0;
916: C->info.mallocs = nspacedouble;
917: C->info.fill_ratio_given = fill;
918: C->info.fill_ratio_needed = afill;
920: #if defined(PETSC_USE_INFO)
921: if (api[am]) {
922: PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
923: PetscCall(PetscInfo(C, "Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n", (double)afill));
924: } else {
925: PetscCall(PetscInfo(C, "Empty matrix product\n"));
926: }
927: #endif
928: PetscFunctionReturn(PETSC_SUCCESS);
929: }
931: /* This function is needed for the seqMPI matrix-matrix multiplication. */
932: /* Three input arrays are merged to one output array. The size of the */
933: /* output array is also output. Duplicate entries only show up once. */
934: static void Merge3SortedArrays(PetscInt size1, PetscInt *in1, PetscInt size2, PetscInt *in2, PetscInt size3, PetscInt *in3, PetscInt *size4, PetscInt *out)
935: {
936: int i = 0, j = 0, k = 0, l = 0;
938: /* Traverse all three arrays */
939: while (i < size1 && j < size2 && k < size3) {
940: if (in1[i] < in2[j] && in1[i] < in3[k]) {
941: out[l++] = in1[i++];
942: } else if (in2[j] < in1[i] && in2[j] < in3[k]) {
943: out[l++] = in2[j++];
944: } else if (in3[k] < in1[i] && in3[k] < in2[j]) {
945: out[l++] = in3[k++];
946: } else if (in1[i] == in2[j] && in1[i] < in3[k]) {
947: out[l++] = in1[i];
948: i++, j++;
949: } else if (in1[i] == in3[k] && in1[i] < in2[j]) {
950: out[l++] = in1[i];
951: i++, k++;
952: } else if (in3[k] == in2[j] && in2[j] < in1[i]) {
953: out[l++] = in2[j];
954: k++, j++;
955: } else if (in1[i] == in2[j] && in1[i] == in3[k]) {
956: out[l++] = in1[i];
957: i++, j++, k++;
958: }
959: }
961: /* Traverse two remaining arrays */
962: while (i < size1 && j < size2) {
963: if (in1[i] < in2[j]) {
964: out[l++] = in1[i++];
965: } else if (in1[i] > in2[j]) {
966: out[l++] = in2[j++];
967: } else {
968: out[l++] = in1[i];
969: i++, j++;
970: }
971: }
973: while (i < size1 && k < size3) {
974: if (in1[i] < in3[k]) {
975: out[l++] = in1[i++];
976: } else if (in1[i] > in3[k]) {
977: out[l++] = in3[k++];
978: } else {
979: out[l++] = in1[i];
980: i++, k++;
981: }
982: }
984: while (k < size3 && j < size2) {
985: if (in3[k] < in2[j]) {
986: out[l++] = in3[k++];
987: } else if (in3[k] > in2[j]) {
988: out[l++] = in2[j++];
989: } else {
990: out[l++] = in3[k];
991: k++, j++;
992: }
993: }
995: /* Traverse one remaining array */
996: while (i < size1) out[l++] = in1[i++];
997: while (j < size2) out[l++] = in2[j++];
998: while (k < size3) out[l++] = in3[k++];
1000: *size4 = l;
1001: }
1003: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */
1004: /* adds up the products. Two of these three multiplications are performed with existing (sequential) */
1005: /* matrix-matrix multiplications. */
1006: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C)
1007: {
1008: MPI_Comm comm;
1009: PetscMPIInt size;
1010: Mat_APMPI *ptap;
1011: PetscFreeSpaceList free_space_diag = NULL, current_space = NULL;
1012: Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
1013: Mat_SeqAIJ *ad = (Mat_SeqAIJ *)a->A->data, *ao = (Mat_SeqAIJ *)a->B->data, *p_loc;
1014: Mat_MPIAIJ *p = (Mat_MPIAIJ *)P->data;
1015: Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq;
1016: PetscInt adponz, adpdnz;
1017: PetscInt *pi_loc, *dnz, *onz;
1018: PetscInt *adi = ad->i, *adj = ad->j, *aoi = ao->i, rstart = A->rmap->rstart;
1019: PetscInt *lnk, i, i1 = 0, pnz, row, *adpoi, *adpoj, *api, *adpoJ, *aopJ, *apJ, *Jptr, aopnz, nspacedouble = 0, j, nzi, *apj, apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
1020: PetscInt am = A->rmap->n, pN = P->cmap->N, pn = P->cmap->n, pm = P->rmap->n, p_colstart, p_colend;
1021: PetscBT lnkbt;
1022: PetscReal afill;
1023: PetscMPIInt rank;
1024: Mat adpd, aopoth;
1025: MatType mtype;
1026: const char *prefix;
1028: PetscFunctionBegin;
1029: MatCheckProduct(C, 4);
1030: PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty");
1031: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1032: PetscCallMPI(MPI_Comm_size(comm, &size));
1033: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1034: PetscCall(MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend));
1036: /* create struct Mat_APMPI and attached it to C later */
1037: PetscCall(PetscNew(&ptap));
1039: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1040: PetscCall(MatGetBrowsOfAoCols_MPIAIJ(A, P, MAT_INITIAL_MATRIX, &ptap->startsj_s, &ptap->startsj_r, &ptap->bufa, &ptap->P_oth));
1042: /* get P_loc by taking all local rows of P */
1043: PetscCall(MatMPIAIJGetLocalMat(P, MAT_INITIAL_MATRIX, &ptap->P_loc));
1045: p_loc = (Mat_SeqAIJ *)ptap->P_loc->data;
1046: pi_loc = p_loc->i;
1048: /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
1049: PetscCall(PetscMalloc1(am + 2, &api));
1050: PetscCall(PetscMalloc1(am + 2, &adpoi));
1052: adpoi[0] = 0;
1053: ptap->api = api;
1054: api[0] = 0;
1056: /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1057: PetscCall(PetscLLCondensedCreate(pN, pN, &lnk, &lnkbt));
1058: MatPreallocateBegin(comm, am, pn, dnz, onz);
1060: /* Symbolic calc of A_loc_diag * P_loc_diag */
1061: PetscCall(MatGetOptionsPrefix(A, &prefix));
1062: PetscCall(MatProductCreate(a->A, p->A, NULL, &adpd));
1063: PetscCall(MatGetOptionsPrefix(A, &prefix));
1064: PetscCall(MatSetOptionsPrefix(adpd, prefix));
1065: PetscCall(MatAppendOptionsPrefix(adpd, "inner_diag_"));
1067: PetscCall(MatProductSetType(adpd, MATPRODUCT_AB));
1068: PetscCall(MatProductSetAlgorithm(adpd, "sorted"));
1069: PetscCall(MatProductSetFill(adpd, fill));
1070: PetscCall(MatProductSetFromOptions(adpd));
1072: adpd->force_diagonals = C->force_diagonals;
1073: PetscCall(MatProductSymbolic(adpd));
1075: adpd_seq = (Mat_SeqAIJ *)((adpd)->data);
1076: adpdi = adpd_seq->i;
1077: adpdj = adpd_seq->j;
1078: p_off = (Mat_SeqAIJ *)p->B->data;
1079: poff_i = p_off->i;
1080: poff_j = p_off->j;
1082: /* j_temp stores indices of a result row before they are added to the linked list */
1083: PetscCall(PetscMalloc1(pN + 2, &j_temp));
1085: /* Symbolic calc of the A_diag * p_loc_off */
1086: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1087: PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(adi[am], PetscIntSumTruncate(aoi[am], pi_loc[pm]))), &free_space_diag));
1088: current_space = free_space_diag;
1090: for (i = 0; i < am; i++) {
1091: /* A_diag * P_loc_off */
1092: nzi = adi[i + 1] - adi[i];
1093: for (j = 0; j < nzi; j++) {
1094: row = *adj++;
1095: pnz = poff_i[row + 1] - poff_i[row];
1096: Jptr = poff_j + poff_i[row];
1097: for (i1 = 0; i1 < pnz; i1++) j_temp[i1] = p->garray[Jptr[i1]];
1098: /* add non-zero cols of P into the sorted linked list lnk */
1099: PetscCall(PetscLLCondensedAddSorted(pnz, j_temp, lnk, lnkbt));
1100: }
1102: adponz = lnk[0];
1103: adpoi[i + 1] = adpoi[i] + adponz;
1105: /* if free space is not available, double the total space in the list */
1106: if (current_space->local_remaining < adponz) {
1107: PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(adponz, current_space->total_array_size), ¤t_space));
1108: nspacedouble++;
1109: }
1111: /* Copy data into free space, then initialize lnk */
1112: PetscCall(PetscLLCondensedClean(pN, adponz, current_space->array, lnk, lnkbt));
1114: current_space->array += adponz;
1115: current_space->local_used += adponz;
1116: current_space->local_remaining -= adponz;
1117: }
1119: /* Symbolic calc of A_off * P_oth */
1120: PetscCall(MatSetOptionsPrefix(a->B, prefix));
1121: PetscCall(MatAppendOptionsPrefix(a->B, "inner_offdiag_"));
1122: PetscCall(MatCreate(PETSC_COMM_SELF, &aopoth));
1123: PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth));
1124: aopoth_seq = (Mat_SeqAIJ *)((aopoth)->data);
1125: aopothi = aopoth_seq->i;
1126: aopothj = aopoth_seq->j;
1128: /* Allocate space for apj, adpj, aopj, ... */
1129: /* destroy lists of free space and other temporary array(s) */
1131: PetscCall(PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am] + 2, &ptap->apj));
1132: PetscCall(PetscMalloc1(adpoi[am] + 2, &adpoj));
1134: /* Copy from linked list to j-array */
1135: PetscCall(PetscFreeSpaceContiguous(&free_space_diag, adpoj));
1136: PetscCall(PetscLLDestroy(lnk, lnkbt));
1138: adpoJ = adpoj;
1139: adpdJ = adpdj;
1140: aopJ = aopothj;
1141: apj = ptap->apj;
1142: apJ = apj; /* still empty */
1144: /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1145: /* A_diag * P_loc_diag to get A*P */
1146: for (i = 0; i < am; i++) {
1147: aopnz = aopothi[i + 1] - aopothi[i];
1148: adponz = adpoi[i + 1] - adpoi[i];
1149: adpdnz = adpdi[i + 1] - adpdi[i];
1151: /* Correct indices from A_diag*P_diag */
1152: for (i1 = 0; i1 < adpdnz; i1++) adpdJ[i1] += p_colstart;
1153: /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1154: Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1155: PetscCall(MatPreallocateSet(i + rstart, apnz, apJ, dnz, onz));
1157: aopJ += aopnz;
1158: adpoJ += adponz;
1159: adpdJ += adpdnz;
1160: apJ += apnz;
1161: api[i + 1] = api[i] + apnz;
1162: }
1164: /* malloc apa to store dense row A[i,:]*P */
1165: PetscCall(PetscCalloc1(pN + 2, &ptap->apa));
1167: /* create and assemble symbolic parallel matrix C */
1168: PetscCall(MatSetSizes(C, am, pn, PETSC_DETERMINE, PETSC_DETERMINE));
1169: PetscCall(MatSetBlockSizesFromMats(C, A, P));
1170: PetscCall(MatGetType(A, &mtype));
1171: PetscCall(MatSetType(C, mtype));
1172: PetscCall(MatMPIAIJSetPreallocation(C, 0, dnz, 0, onz));
1173: MatPreallocateEnd(dnz, onz);
1175: PetscCall(MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api));
1176: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1177: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1178: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1180: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1181: C->ops->productnumeric = MatProductNumeric_AB;
1183: /* attach the supporting struct to C for reuse */
1184: C->product->data = ptap;
1185: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
1187: /* set MatInfo */
1188: afill = (PetscReal)api[am] / (adi[am] + aoi[am] + pi_loc[pm] + 1) + 1.e-5;
1189: if (afill < 1.0) afill = 1.0;
1190: C->info.mallocs = nspacedouble;
1191: C->info.fill_ratio_given = fill;
1192: C->info.fill_ratio_needed = afill;
1194: #if defined(PETSC_USE_INFO)
1195: if (api[am]) {
1196: PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
1197: PetscCall(PetscInfo(C, "Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n", (double)afill));
1198: } else {
1199: PetscCall(PetscInfo(C, "Empty matrix product\n"));
1200: }
1201: #endif
1203: PetscCall(MatDestroy(&aopoth));
1204: PetscCall(MatDestroy(&adpd));
1205: PetscCall(PetscFree(j_temp));
1206: PetscCall(PetscFree(adpoj));
1207: PetscCall(PetscFree(adpoi));
1208: PetscFunctionReturn(PETSC_SUCCESS);
1209: }
1211: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1212: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P, Mat A, Mat C)
1213: {
1214: Mat_APMPI *ptap;
1215: Mat Pt;
1217: PetscFunctionBegin;
1218: MatCheckProduct(C, 3);
1219: ptap = (Mat_APMPI *)C->product->data;
1220: PetscCheck(ptap, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtAP cannot be computed. Missing data");
1221: PetscCheck(ptap->Pt, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtA cannot be reused. Do not call MatProductClear()");
1223: Pt = ptap->Pt;
1224: PetscCall(MatTransposeSetPrecursor(P, Pt));
1225: PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &Pt));
1226: PetscCall(MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt, A, C));
1227: PetscFunctionReturn(PETSC_SUCCESS);
1228: }
1230: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1231: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P, Mat A, PetscReal fill, Mat C)
1232: {
1233: Mat_APMPI *ptap;
1234: Mat_MPIAIJ *p = (Mat_MPIAIJ *)P->data;
1235: MPI_Comm comm;
1236: PetscMPIInt size, rank;
1237: PetscFreeSpaceList free_space = NULL, current_space = NULL;
1238: PetscInt pn = P->cmap->n, aN = A->cmap->N, an = A->cmap->n;
1239: PetscInt *lnk, i, k, nsend, rstart;
1240: PetscBT lnkbt;
1241: PetscMPIInt tagi, tagj, *len_si, *len_s, *len_ri, nrecv;
1242: PETSC_UNUSED PetscMPIInt icompleted = 0;
1243: PetscInt **buf_rj, **buf_ri, **buf_ri_k, row, ncols, *cols;
1244: PetscInt len, proc, *dnz, *onz, *owners, nzi;
1245: PetscInt nrows, *buf_s, *buf_si, *buf_si_i, **nextrow, **nextci;
1246: MPI_Request *swaits, *rwaits;
1247: MPI_Status *sstatus, rstatus;
1248: PetscLayout rowmap;
1249: PetscInt *owners_co, *coi, *coj; /* i and j array of (p->B)^T*A*P - used in the communication */
1250: PetscMPIInt *len_r, *id_r; /* array of length of comm->size, store send/recv matrix values */
1251: PetscInt *Jptr, *prmap = p->garray, con, j, Crmax;
1252: Mat_SeqAIJ *a_loc, *c_loc, *c_oth;
1253: PetscHMapI ta;
1254: MatType mtype;
1255: const char *prefix;
1257: PetscFunctionBegin;
1258: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1259: PetscCallMPI(MPI_Comm_size(comm, &size));
1260: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1262: /* create symbolic parallel matrix C */
1263: PetscCall(MatGetType(A, &mtype));
1264: PetscCall(MatSetType(C, mtype));
1266: C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1268: /* create struct Mat_APMPI and attached it to C later */
1269: PetscCall(PetscNew(&ptap));
1271: /* (0) compute Rd = Pd^T, Ro = Po^T */
1272: PetscCall(MatTranspose(p->A, MAT_INITIAL_MATRIX, &ptap->Rd));
1273: PetscCall(MatTranspose(p->B, MAT_INITIAL_MATRIX, &ptap->Ro));
1275: /* (1) compute symbolic A_loc */
1276: PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &ptap->A_loc));
1278: /* (2-1) compute symbolic C_oth = Ro*A_loc */
1279: PetscCall(MatGetOptionsPrefix(A, &prefix));
1280: PetscCall(MatSetOptionsPrefix(ptap->Ro, prefix));
1281: PetscCall(MatAppendOptionsPrefix(ptap->Ro, "inner_offdiag_"));
1282: PetscCall(MatCreate(PETSC_COMM_SELF, &ptap->C_oth));
1283: PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro, ptap->A_loc, fill, ptap->C_oth));
1285: /* (3) send coj of C_oth to other processors */
1286: /* determine row ownership */
1287: PetscCall(PetscLayoutCreate(comm, &rowmap));
1288: rowmap->n = pn;
1289: rowmap->bs = 1;
1290: PetscCall(PetscLayoutSetUp(rowmap));
1291: owners = rowmap->range;
1293: /* determine the number of messages to send, their lengths */
1294: PetscCall(PetscMalloc4(size, &len_s, size, &len_si, size, &sstatus, size + 2, &owners_co));
1295: PetscCall(PetscArrayzero(len_s, size));
1296: PetscCall(PetscArrayzero(len_si, size));
1298: c_oth = (Mat_SeqAIJ *)ptap->C_oth->data;
1299: coi = c_oth->i;
1300: coj = c_oth->j;
1301: con = ptap->C_oth->rmap->n;
1302: proc = 0;
1303: for (i = 0; i < con; i++) {
1304: while (prmap[i] >= owners[proc + 1]) proc++;
1305: len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */
1306: len_s[proc] += coi[i + 1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1307: }
1309: len = 0; /* max length of buf_si[], see (4) */
1310: owners_co[0] = 0;
1311: nsend = 0;
1312: for (proc = 0; proc < size; proc++) {
1313: owners_co[proc + 1] = owners_co[proc] + len_si[proc];
1314: if (len_s[proc]) {
1315: nsend++;
1316: len_si[proc] = 2 * (len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1317: len += len_si[proc];
1318: }
1319: }
1321: /* determine the number and length of messages to receive for coi and coj */
1322: PetscCall(PetscGatherNumberOfMessages(comm, NULL, len_s, &nrecv));
1323: PetscCall(PetscGatherMessageLengths2(comm, nsend, nrecv, len_s, len_si, &id_r, &len_r, &len_ri));
1325: /* post the Irecv and Isend of coj */
1326: PetscCall(PetscCommGetNewTag(comm, &tagj));
1327: PetscCall(PetscPostIrecvInt(comm, tagj, nrecv, id_r, len_r, &buf_rj, &rwaits));
1328: PetscCall(PetscMalloc1(nsend + 1, &swaits));
1329: for (proc = 0, k = 0; proc < size; proc++) {
1330: if (!len_s[proc]) continue;
1331: i = owners_co[proc];
1332: PetscCallMPI(MPI_Isend(coj + coi[i], len_s[proc], MPIU_INT, proc, tagj, comm, swaits + k));
1333: k++;
1334: }
1336: /* (2-2) compute symbolic C_loc = Rd*A_loc */
1337: PetscCall(MatSetOptionsPrefix(ptap->Rd, prefix));
1338: PetscCall(MatAppendOptionsPrefix(ptap->Rd, "inner_diag_"));
1339: PetscCall(MatCreate(PETSC_COMM_SELF, &ptap->C_loc));
1340: PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd, ptap->A_loc, fill, ptap->C_loc));
1341: c_loc = (Mat_SeqAIJ *)ptap->C_loc->data;
1343: /* receives coj are complete */
1344: for (i = 0; i < nrecv; i++) PetscCallMPI(MPI_Waitany(nrecv, rwaits, &icompleted, &rstatus));
1345: PetscCall(PetscFree(rwaits));
1346: if (nsend) PetscCallMPI(MPI_Waitall(nsend, swaits, sstatus));
1348: /* add received column indices into ta to update Crmax */
1349: a_loc = (Mat_SeqAIJ *)ptap->A_loc->data;
1351: /* create and initialize a linked list */
1352: PetscCall(PetscHMapICreateWithSize(an, &ta)); /* for compute Crmax */
1353: MatRowMergeMax_SeqAIJ(a_loc, ptap->A_loc->rmap->N, ta);
1355: for (k = 0; k < nrecv; k++) { /* k-th received message */
1356: Jptr = buf_rj[k];
1357: for (j = 0; j < len_r[k]; j++) PetscCall(PetscHMapISet(ta, *(Jptr + j) + 1, 1));
1358: }
1359: PetscCall(PetscHMapIGetSize(ta, &Crmax));
1360: PetscCall(PetscHMapIDestroy(&ta));
1362: /* (4) send and recv coi */
1363: PetscCall(PetscCommGetNewTag(comm, &tagi));
1364: PetscCall(PetscPostIrecvInt(comm, tagi, nrecv, id_r, len_ri, &buf_ri, &rwaits));
1365: PetscCall(PetscMalloc1(len + 1, &buf_s));
1366: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1367: for (proc = 0, k = 0; proc < size; proc++) {
1368: if (!len_s[proc]) continue;
1369: /* form outgoing message for i-structure:
1370: buf_si[0]: nrows to be sent
1371: [1:nrows]: row index (global)
1372: [nrows+1:2*nrows+1]: i-structure index
1373: */
1374: nrows = len_si[proc] / 2 - 1; /* num of rows in Co to be sent to [proc] */
1375: buf_si_i = buf_si + nrows + 1;
1376: buf_si[0] = nrows;
1377: buf_si_i[0] = 0;
1378: nrows = 0;
1379: for (i = owners_co[proc]; i < owners_co[proc + 1]; i++) {
1380: nzi = coi[i + 1] - coi[i];
1381: buf_si_i[nrows + 1] = buf_si_i[nrows] + nzi; /* i-structure */
1382: buf_si[nrows + 1] = prmap[i] - owners[proc]; /* local row index */
1383: nrows++;
1384: }
1385: PetscCallMPI(MPI_Isend(buf_si, len_si[proc], MPIU_INT, proc, tagi, comm, swaits + k));
1386: k++;
1387: buf_si += len_si[proc];
1388: }
1389: for (i = 0; i < nrecv; i++) PetscCallMPI(MPI_Waitany(nrecv, rwaits, &icompleted, &rstatus));
1390: PetscCall(PetscFree(rwaits));
1391: if (nsend) PetscCallMPI(MPI_Waitall(nsend, swaits, sstatus));
1393: PetscCall(PetscFree4(len_s, len_si, sstatus, owners_co));
1394: PetscCall(PetscFree(len_ri));
1395: PetscCall(PetscFree(swaits));
1396: PetscCall(PetscFree(buf_s));
1398: /* (5) compute the local portion of C */
1399: /* set initial free space to be Crmax, sufficient for holding nonzeros in each row of C */
1400: PetscCall(PetscFreeSpaceGet(Crmax, &free_space));
1401: current_space = free_space;
1403: PetscCall(PetscMalloc3(nrecv, &buf_ri_k, nrecv, &nextrow, nrecv, &nextci));
1404: for (k = 0; k < nrecv; k++) {
1405: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1406: nrows = *buf_ri_k[k];
1407: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1408: nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure */
1409: }
1411: MatPreallocateBegin(comm, pn, an, dnz, onz);
1412: PetscCall(PetscLLCondensedCreate(Crmax, aN, &lnk, &lnkbt));
1413: for (i = 0; i < pn; i++) { /* for each local row of C */
1414: /* add C_loc into C */
1415: nzi = c_loc->i[i + 1] - c_loc->i[i];
1416: Jptr = c_loc->j + c_loc->i[i];
1417: PetscCall(PetscLLCondensedAddSorted(nzi, Jptr, lnk, lnkbt));
1419: /* add received col data into lnk */
1420: for (k = 0; k < nrecv; k++) { /* k-th received message */
1421: if (i == *nextrow[k]) { /* i-th row */
1422: nzi = *(nextci[k] + 1) - *nextci[k];
1423: Jptr = buf_rj[k] + *nextci[k];
1424: PetscCall(PetscLLCondensedAddSorted(nzi, Jptr, lnk, lnkbt));
1425: nextrow[k]++;
1426: nextci[k]++;
1427: }
1428: }
1430: /* add missing diagonal entry */
1431: if (C->force_diagonals) {
1432: k = i + owners[rank]; /* column index */
1433: PetscCall(PetscLLCondensedAddSorted(1, &k, lnk, lnkbt));
1434: }
1436: nzi = lnk[0];
1438: /* copy data into free space, then initialize lnk */
1439: PetscCall(PetscLLCondensedClean(aN, nzi, current_space->array, lnk, lnkbt));
1440: PetscCall(MatPreallocateSet(i + owners[rank], nzi, current_space->array, dnz, onz));
1441: }
1442: PetscCall(PetscFree3(buf_ri_k, nextrow, nextci));
1443: PetscCall(PetscLLDestroy(lnk, lnkbt));
1444: PetscCall(PetscFreeSpaceDestroy(free_space));
1446: /* local sizes and preallocation */
1447: PetscCall(MatSetSizes(C, pn, an, PETSC_DETERMINE, PETSC_DETERMINE));
1448: if (P->cmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(C->rmap, P->cmap->bs));
1449: if (A->cmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(C->cmap, A->cmap->bs));
1450: PetscCall(MatMPIAIJSetPreallocation(C, 0, dnz, 0, onz));
1451: MatPreallocateEnd(dnz, onz);
1453: /* add C_loc and C_oth to C */
1454: PetscCall(MatGetOwnershipRange(C, &rstart, NULL));
1455: for (i = 0; i < pn; i++) {
1456: ncols = c_loc->i[i + 1] - c_loc->i[i];
1457: cols = c_loc->j + c_loc->i[i];
1458: row = rstart + i;
1459: PetscCall(MatSetValues(C, 1, (const PetscInt *)&row, ncols, (const PetscInt *)cols, NULL, INSERT_VALUES));
1461: if (C->force_diagonals) PetscCall(MatSetValues(C, 1, (const PetscInt *)&row, 1, (const PetscInt *)&row, NULL, INSERT_VALUES));
1462: }
1463: for (i = 0; i < con; i++) {
1464: ncols = c_oth->i[i + 1] - c_oth->i[i];
1465: cols = c_oth->j + c_oth->i[i];
1466: row = prmap[i];
1467: PetscCall(MatSetValues(C, 1, (const PetscInt *)&row, ncols, (const PetscInt *)cols, NULL, INSERT_VALUES));
1468: }
1469: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1470: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1471: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1473: /* members in merge */
1474: PetscCall(PetscFree(id_r));
1475: PetscCall(PetscFree(len_r));
1476: PetscCall(PetscFree(buf_ri[0]));
1477: PetscCall(PetscFree(buf_ri));
1478: PetscCall(PetscFree(buf_rj[0]));
1479: PetscCall(PetscFree(buf_rj));
1480: PetscCall(PetscLayoutDestroy(&rowmap));
1482: /* attach the supporting struct to C for reuse */
1483: C->product->data = ptap;
1484: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
1485: PetscFunctionReturn(PETSC_SUCCESS);
1486: }
1488: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P, Mat A, Mat C)
1489: {
1490: Mat_MPIAIJ *p = (Mat_MPIAIJ *)P->data;
1491: Mat_SeqAIJ *c_seq;
1492: Mat_APMPI *ptap;
1493: Mat A_loc, C_loc, C_oth;
1494: PetscInt i, rstart, rend, cm, ncols, row;
1495: const PetscInt *cols;
1496: const PetscScalar *vals;
1498: PetscFunctionBegin;
1499: MatCheckProduct(C, 3);
1500: ptap = (Mat_APMPI *)C->product->data;
1501: PetscCheck(ptap, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtAP cannot be computed. Missing data");
1502: PetscCheck(ptap->A_loc, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtA cannot be reused. Do not call MatProductClear()");
1503: PetscCall(MatZeroEntries(C));
1505: /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1506: /* 1) get R = Pd^T, Ro = Po^T */
1507: PetscCall(MatTransposeSetPrecursor(p->A, ptap->Rd));
1508: PetscCall(MatTranspose(p->A, MAT_REUSE_MATRIX, &ptap->Rd));
1509: PetscCall(MatTransposeSetPrecursor(p->B, ptap->Ro));
1510: PetscCall(MatTranspose(p->B, MAT_REUSE_MATRIX, &ptap->Ro));
1512: /* 2) compute numeric A_loc */
1513: PetscCall(MatMPIAIJGetLocalMat(A, MAT_REUSE_MATRIX, &ptap->A_loc));
1515: /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1516: A_loc = ptap->A_loc;
1517: PetscCall(ptap->C_loc->ops->matmultnumeric(ptap->Rd, A_loc, ptap->C_loc));
1518: PetscCall(ptap->C_oth->ops->matmultnumeric(ptap->Ro, A_loc, ptap->C_oth));
1519: C_loc = ptap->C_loc;
1520: C_oth = ptap->C_oth;
1522: /* add C_loc and C_oth to C */
1523: PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
1525: /* C_loc -> C */
1526: cm = C_loc->rmap->N;
1527: c_seq = (Mat_SeqAIJ *)C_loc->data;
1528: cols = c_seq->j;
1529: vals = c_seq->a;
1530: for (i = 0; i < cm; i++) {
1531: ncols = c_seq->i[i + 1] - c_seq->i[i];
1532: row = rstart + i;
1533: PetscCall(MatSetValues(C, 1, &row, ncols, cols, vals, ADD_VALUES));
1534: cols += ncols;
1535: vals += ncols;
1536: }
1538: /* Co -> C, off-processor part */
1539: cm = C_oth->rmap->N;
1540: c_seq = (Mat_SeqAIJ *)C_oth->data;
1541: cols = c_seq->j;
1542: vals = c_seq->a;
1543: for (i = 0; i < cm; i++) {
1544: ncols = c_seq->i[i + 1] - c_seq->i[i];
1545: row = p->garray[i];
1546: PetscCall(MatSetValues(C, 1, &row, ncols, cols, vals, ADD_VALUES));
1547: cols += ncols;
1548: vals += ncols;
1549: }
1550: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1551: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1552: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1553: PetscFunctionReturn(PETSC_SUCCESS);
1554: }
1556: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P, Mat A, Mat C)
1557: {
1558: Mat_Merge_SeqsToMPI *merge;
1559: Mat_MPIAIJ *p = (Mat_MPIAIJ *)P->data;
1560: Mat_SeqAIJ *pd = (Mat_SeqAIJ *)p->A->data, *po = (Mat_SeqAIJ *)p->B->data;
1561: Mat_APMPI *ap;
1562: PetscInt *adj;
1563: PetscInt i, j, k, anz, pnz, row, *cj, nexta;
1564: MatScalar *ada, *ca, valtmp;
1565: PetscInt am = A->rmap->n, cm = C->rmap->n, pon = (p->B)->cmap->n;
1566: MPI_Comm comm;
1567: PetscMPIInt size, rank, taga, *len_s;
1568: PetscInt *owners, proc, nrows, **buf_ri_k, **nextrow, **nextci;
1569: PetscInt **buf_ri, **buf_rj;
1570: PetscInt cnz = 0, *bj_i, *bi, *bj, bnz, nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1571: MPI_Request *s_waits, *r_waits;
1572: MPI_Status *status;
1573: MatScalar **abuf_r, *ba_i, *pA, *coa, *ba;
1574: const PetscScalar *dummy;
1575: PetscInt *ai, *aj, *coi, *coj, *poJ, *pdJ;
1576: Mat A_loc;
1577: Mat_SeqAIJ *a_loc;
1579: PetscFunctionBegin;
1580: MatCheckProduct(C, 3);
1581: ap = (Mat_APMPI *)C->product->data;
1582: PetscCheck(ap, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtA cannot be computed. Missing data");
1583: PetscCheck(ap->A_loc, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "PtA cannot be reused. Do not call MatProductClear()");
1584: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
1585: PetscCallMPI(MPI_Comm_size(comm, &size));
1586: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1588: merge = ap->merge;
1590: /* 2) compute numeric C_seq = P_loc^T*A_loc */
1591: /* get data from symbolic products */
1592: coi = merge->coi;
1593: coj = merge->coj;
1594: PetscCall(PetscCalloc1(coi[pon] + 1, &coa));
1595: bi = merge->bi;
1596: bj = merge->bj;
1597: owners = merge->rowmap->range;
1598: PetscCall(PetscCalloc1(bi[cm] + 1, &ba));
1600: /* get A_loc by taking all local rows of A */
1601: A_loc = ap->A_loc;
1602: PetscCall(MatMPIAIJGetLocalMat(A, MAT_REUSE_MATRIX, &A_loc));
1603: a_loc = (Mat_SeqAIJ *)(A_loc)->data;
1604: ai = a_loc->i;
1605: aj = a_loc->j;
1607: /* trigger copy to CPU */
1608: PetscCall(MatSeqAIJGetArrayRead(p->A, &dummy));
1609: PetscCall(MatSeqAIJRestoreArrayRead(p->A, &dummy));
1610: PetscCall(MatSeqAIJGetArrayRead(p->B, &dummy));
1611: PetscCall(MatSeqAIJRestoreArrayRead(p->B, &dummy));
1612: for (i = 0; i < am; i++) {
1613: anz = ai[i + 1] - ai[i];
1614: adj = aj + ai[i];
1615: ada = a_loc->a + ai[i];
1617: /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1618: /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1619: pnz = po->i[i + 1] - po->i[i];
1620: poJ = po->j + po->i[i];
1621: pA = po->a + po->i[i];
1622: for (j = 0; j < pnz; j++) {
1623: row = poJ[j];
1624: cj = coj + coi[row];
1625: ca = coa + coi[row];
1626: /* perform sparse axpy */
1627: nexta = 0;
1628: valtmp = pA[j];
1629: for (k = 0; nexta < anz; k++) {
1630: if (cj[k] == adj[nexta]) {
1631: ca[k] += valtmp * ada[nexta];
1632: nexta++;
1633: }
1634: }
1635: PetscCall(PetscLogFlops(2.0 * anz));
1636: }
1638: /* put the value into Cd (diagonal part) */
1639: pnz = pd->i[i + 1] - pd->i[i];
1640: pdJ = pd->j + pd->i[i];
1641: pA = pd->a + pd->i[i];
1642: for (j = 0; j < pnz; j++) {
1643: row = pdJ[j];
1644: cj = bj + bi[row];
1645: ca = ba + bi[row];
1646: /* perform sparse axpy */
1647: nexta = 0;
1648: valtmp = pA[j];
1649: for (k = 0; nexta < anz; k++) {
1650: if (cj[k] == adj[nexta]) {
1651: ca[k] += valtmp * ada[nexta];
1652: nexta++;
1653: }
1654: }
1655: PetscCall(PetscLogFlops(2.0 * anz));
1656: }
1657: }
1659: /* 3) send and recv matrix values coa */
1660: buf_ri = merge->buf_ri;
1661: buf_rj = merge->buf_rj;
1662: len_s = merge->len_s;
1663: PetscCall(PetscCommGetNewTag(comm, &taga));
1664: PetscCall(PetscPostIrecvScalar(comm, taga, merge->nrecv, merge->id_r, merge->len_r, &abuf_r, &r_waits));
1666: PetscCall(PetscMalloc2(merge->nsend + 1, &s_waits, size, &status));
1667: for (proc = 0, k = 0; proc < size; proc++) {
1668: if (!len_s[proc]) continue;
1669: i = merge->owners_co[proc];
1670: PetscCallMPI(MPI_Isend(coa + coi[i], len_s[proc], MPIU_MATSCALAR, proc, taga, comm, s_waits + k));
1671: k++;
1672: }
1673: if (merge->nrecv) PetscCallMPI(MPI_Waitall(merge->nrecv, r_waits, status));
1674: if (merge->nsend) PetscCallMPI(MPI_Waitall(merge->nsend, s_waits, status));
1676: PetscCall(PetscFree2(s_waits, status));
1677: PetscCall(PetscFree(r_waits));
1678: PetscCall(PetscFree(coa));
1680: /* 4) insert local Cseq and received values into Cmpi */
1681: PetscCall(PetscMalloc3(merge->nrecv, &buf_ri_k, merge->nrecv, &nextrow, merge->nrecv, &nextci));
1682: for (k = 0; k < merge->nrecv; k++) {
1683: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1684: nrows = *buf_ri_k[k];
1685: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1686: nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure */
1687: }
1689: for (i = 0; i < cm; i++) {
1690: row = owners[rank] + i; /* global row index of C_seq */
1691: bj_i = bj + bi[i]; /* col indices of the i-th row of C */
1692: ba_i = ba + bi[i];
1693: bnz = bi[i + 1] - bi[i];
1694: /* add received vals into ba */
1695: for (k = 0; k < merge->nrecv; k++) { /* k-th received message */
1696: /* i-th row */
1697: if (i == *nextrow[k]) {
1698: cnz = *(nextci[k] + 1) - *nextci[k];
1699: cj = buf_rj[k] + *nextci[k];
1700: ca = abuf_r[k] + *nextci[k];
1701: nextcj = 0;
1702: for (j = 0; nextcj < cnz; j++) {
1703: if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1704: ba_i[j] += ca[nextcj++];
1705: }
1706: }
1707: nextrow[k]++;
1708: nextci[k]++;
1709: PetscCall(PetscLogFlops(2.0 * cnz));
1710: }
1711: }
1712: PetscCall(MatSetValues(C, 1, &row, bnz, bj_i, ba_i, INSERT_VALUES));
1713: }
1714: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1715: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1717: PetscCall(PetscFree(ba));
1718: PetscCall(PetscFree(abuf_r[0]));
1719: PetscCall(PetscFree(abuf_r));
1720: PetscCall(PetscFree3(buf_ri_k, nextrow, nextci));
1721: PetscFunctionReturn(PETSC_SUCCESS);
1722: }
1724: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P, Mat A, PetscReal fill, Mat C)
1725: {
1726: Mat A_loc;
1727: Mat_APMPI *ap;
1728: PetscFreeSpaceList free_space = NULL, current_space = NULL;
1729: Mat_MPIAIJ *p = (Mat_MPIAIJ *)P->data, *a = (Mat_MPIAIJ *)A->data;
1730: PetscInt *pdti, *pdtj, *poti, *potj, *ptJ;
1731: PetscInt nnz;
1732: PetscInt *lnk, *owners_co, *coi, *coj, i, k, pnz, row;
1733: PetscInt am = A->rmap->n, pn = P->cmap->n;
1734: MPI_Comm comm;
1735: PetscMPIInt size, rank, tagi, tagj, *len_si, *len_s, *len_ri;
1736: PetscInt **buf_rj, **buf_ri, **buf_ri_k;
1737: PetscInt len, proc, *dnz, *onz, *owners;
1738: PetscInt nzi, *bi, *bj;
1739: PetscInt nrows, *buf_s, *buf_si, *buf_si_i, **nextrow, **nextci;
1740: MPI_Request *swaits, *rwaits;
1741: MPI_Status *sstatus, rstatus;
1742: Mat_Merge_SeqsToMPI *merge;
1743: PetscInt *ai, *aj, *Jptr, anz, *prmap = p->garray, pon, nspacedouble = 0, j;
1744: PetscReal afill = 1.0, afill_tmp;
1745: PetscInt rstart = P->cmap->rstart, rmax, Armax;
1746: Mat_SeqAIJ *a_loc;
1747: PetscHMapI ta;
1748: MatType mtype;
1750: PetscFunctionBegin;
1751: PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1752: /* check if matrix local sizes are compatible */
1753: PetscCheck(A->rmap->rstart == P->rmap->rstart && A->rmap->rend == P->rmap->rend, comm, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, A (%" PetscInt_FMT ", %" PetscInt_FMT ") != P (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->rstart,
1754: A->rmap->rend, P->rmap->rstart, P->rmap->rend);
1756: PetscCallMPI(MPI_Comm_size(comm, &size));
1757: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1759: /* create struct Mat_APMPI and attached it to C later */
1760: PetscCall(PetscNew(&ap));
1762: /* get A_loc by taking all local rows of A */
1763: PetscCall(MatMPIAIJGetLocalMat(A, MAT_INITIAL_MATRIX, &A_loc));
1765: ap->A_loc = A_loc;
1766: a_loc = (Mat_SeqAIJ *)(A_loc)->data;
1767: ai = a_loc->i;
1768: aj = a_loc->j;
1770: /* determine symbolic Co=(p->B)^T*A - send to others */
1771: PetscCall(MatGetSymbolicTranspose_SeqAIJ(p->A, &pdti, &pdtj));
1772: PetscCall(MatGetSymbolicTranspose_SeqAIJ(p->B, &poti, &potj));
1773: pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1774: >= (num of nonzero rows of C_seq) - pn */
1775: PetscCall(PetscMalloc1(pon + 1, &coi));
1776: coi[0] = 0;
1778: /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1779: nnz = PetscRealIntMultTruncate(fill, PetscIntSumTruncate(poti[pon], ai[am]));
1780: PetscCall(PetscFreeSpaceGet(nnz, &free_space));
1781: current_space = free_space;
1783: /* create and initialize a linked list */
1784: PetscCall(PetscHMapICreateWithSize(A->cmap->n + a->B->cmap->N, &ta));
1785: MatRowMergeMax_SeqAIJ(a_loc, am, ta);
1786: PetscCall(PetscHMapIGetSize(ta, &Armax));
1788: PetscCall(PetscLLCondensedCreate_Scalable(Armax, &lnk));
1790: for (i = 0; i < pon; i++) {
1791: pnz = poti[i + 1] - poti[i];
1792: ptJ = potj + poti[i];
1793: for (j = 0; j < pnz; j++) {
1794: row = ptJ[j]; /* row of A_loc == col of Pot */
1795: anz = ai[row + 1] - ai[row];
1796: Jptr = aj + ai[row];
1797: /* add non-zero cols of AP into the sorted linked list lnk */
1798: PetscCall(PetscLLCondensedAddSorted_Scalable(anz, Jptr, lnk));
1799: }
1800: nnz = lnk[0];
1802: /* If free space is not available, double the total space in the list */
1803: if (current_space->local_remaining < nnz) {
1804: PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(nnz, current_space->total_array_size), ¤t_space));
1805: nspacedouble++;
1806: }
1808: /* Copy data into free space, and zero out denserows */
1809: PetscCall(PetscLLCondensedClean_Scalable(nnz, current_space->array, lnk));
1811: current_space->array += nnz;
1812: current_space->local_used += nnz;
1813: current_space->local_remaining -= nnz;
1815: coi[i + 1] = coi[i] + nnz;
1816: }
1818: PetscCall(PetscMalloc1(coi[pon] + 1, &coj));
1819: PetscCall(PetscFreeSpaceContiguous(&free_space, coj));
1820: PetscCall(PetscLLCondensedDestroy_Scalable(lnk)); /* must destroy to get a new one for C */
1822: afill_tmp = (PetscReal)coi[pon] / (poti[pon] + ai[am] + 1);
1823: if (afill_tmp > afill) afill = afill_tmp;
1825: /* send j-array (coj) of Co to other processors */
1826: /* determine row ownership */
1827: PetscCall(PetscNew(&merge));
1828: PetscCall(PetscLayoutCreate(comm, &merge->rowmap));
1830: merge->rowmap->n = pn;
1831: merge->rowmap->bs = 1;
1833: PetscCall(PetscLayoutSetUp(merge->rowmap));
1834: owners = merge->rowmap->range;
1836: /* determine the number of messages to send, their lengths */
1837: PetscCall(PetscCalloc1(size, &len_si));
1838: PetscCall(PetscCalloc1(size, &merge->len_s));
1840: len_s = merge->len_s;
1841: merge->nsend = 0;
1843: PetscCall(PetscMalloc1(size + 2, &owners_co));
1845: proc = 0;
1846: for (i = 0; i < pon; i++) {
1847: while (prmap[i] >= owners[proc + 1]) proc++;
1848: len_si[proc]++; /* num of rows in Co to be sent to [proc] */
1849: len_s[proc] += coi[i + 1] - coi[i];
1850: }
1852: len = 0; /* max length of buf_si[] */
1853: owners_co[0] = 0;
1854: for (proc = 0; proc < size; proc++) {
1855: owners_co[proc + 1] = owners_co[proc] + len_si[proc];
1856: if (len_si[proc]) {
1857: merge->nsend++;
1858: len_si[proc] = 2 * (len_si[proc] + 1);
1859: len += len_si[proc];
1860: }
1861: }
1863: /* determine the number and length of messages to receive for coi and coj */
1864: PetscCall(PetscGatherNumberOfMessages(comm, NULL, len_s, &merge->nrecv));
1865: PetscCall(PetscGatherMessageLengths2(comm, merge->nsend, merge->nrecv, len_s, len_si, &merge->id_r, &merge->len_r, &len_ri));
1867: /* post the Irecv and Isend of coj */
1868: PetscCall(PetscCommGetNewTag(comm, &tagj));
1869: PetscCall(PetscPostIrecvInt(comm, tagj, merge->nrecv, merge->id_r, merge->len_r, &buf_rj, &rwaits));
1870: PetscCall(PetscMalloc1(merge->nsend + 1, &swaits));
1871: for (proc = 0, k = 0; proc < size; proc++) {
1872: if (!len_s[proc]) continue;
1873: i = owners_co[proc];
1874: PetscCallMPI(MPI_Isend(coj + coi[i], len_s[proc], MPIU_INT, proc, tagj, comm, swaits + k));
1875: k++;
1876: }
1878: /* receives and sends of coj are complete */
1879: PetscCall(PetscMalloc1(size, &sstatus));
1880: for (i = 0; i < merge->nrecv; i++) {
1881: PETSC_UNUSED PetscMPIInt icompleted;
1882: PetscCallMPI(MPI_Waitany(merge->nrecv, rwaits, &icompleted, &rstatus));
1883: }
1884: PetscCall(PetscFree(rwaits));
1885: if (merge->nsend) PetscCallMPI(MPI_Waitall(merge->nsend, swaits, sstatus));
1887: /* add received column indices into table to update Armax */
1888: /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */
1889: for (k = 0; k < merge->nrecv; k++) { /* k-th received message */
1890: Jptr = buf_rj[k];
1891: for (j = 0; j < merge->len_r[k]; j++) PetscCall(PetscHMapISet(ta, *(Jptr + j) + 1, 1));
1892: }
1893: PetscCall(PetscHMapIGetSize(ta, &Armax));
1895: /* send and recv coi */
1896: PetscCall(PetscCommGetNewTag(comm, &tagi));
1897: PetscCall(PetscPostIrecvInt(comm, tagi, merge->nrecv, merge->id_r, len_ri, &buf_ri, &rwaits));
1898: PetscCall(PetscMalloc1(len + 1, &buf_s));
1899: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1900: for (proc = 0, k = 0; proc < size; proc++) {
1901: if (!len_s[proc]) continue;
1902: /* form outgoing message for i-structure:
1903: buf_si[0]: nrows to be sent
1904: [1:nrows]: row index (global)
1905: [nrows+1:2*nrows+1]: i-structure index
1906: */
1907: nrows = len_si[proc] / 2 - 1;
1908: buf_si_i = buf_si + nrows + 1;
1909: buf_si[0] = nrows;
1910: buf_si_i[0] = 0;
1911: nrows = 0;
1912: for (i = owners_co[proc]; i < owners_co[proc + 1]; i++) {
1913: nzi = coi[i + 1] - coi[i];
1914: buf_si_i[nrows + 1] = buf_si_i[nrows] + nzi; /* i-structure */
1915: buf_si[nrows + 1] = prmap[i] - owners[proc]; /* local row index */
1916: nrows++;
1917: }
1918: PetscCallMPI(MPI_Isend(buf_si, len_si[proc], MPIU_INT, proc, tagi, comm, swaits + k));
1919: k++;
1920: buf_si += len_si[proc];
1921: }
1922: i = merge->nrecv;
1923: while (i--) {
1924: PETSC_UNUSED PetscMPIInt icompleted;
1925: PetscCallMPI(MPI_Waitany(merge->nrecv, rwaits, &icompleted, &rstatus));
1926: }
1927: PetscCall(PetscFree(rwaits));
1928: if (merge->nsend) PetscCallMPI(MPI_Waitall(merge->nsend, swaits, sstatus));
1929: PetscCall(PetscFree(len_si));
1930: PetscCall(PetscFree(len_ri));
1931: PetscCall(PetscFree(swaits));
1932: PetscCall(PetscFree(sstatus));
1933: PetscCall(PetscFree(buf_s));
1935: /* compute the local portion of C (mpi mat) */
1936: /* allocate bi array and free space for accumulating nonzero column info */
1937: PetscCall(PetscMalloc1(pn + 1, &bi));
1938: bi[0] = 0;
1940: /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1941: nnz = PetscRealIntMultTruncate(fill, PetscIntSumTruncate(pdti[pn], PetscIntSumTruncate(poti[pon], ai[am])));
1942: PetscCall(PetscFreeSpaceGet(nnz, &free_space));
1943: current_space = free_space;
1945: PetscCall(PetscMalloc3(merge->nrecv, &buf_ri_k, merge->nrecv, &nextrow, merge->nrecv, &nextci));
1946: for (k = 0; k < merge->nrecv; k++) {
1947: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1948: nrows = *buf_ri_k[k];
1949: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1950: nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th received i-structure */
1951: }
1953: PetscCall(PetscLLCondensedCreate_Scalable(Armax, &lnk));
1954: MatPreallocateBegin(comm, pn, A->cmap->n, dnz, onz);
1955: rmax = 0;
1956: for (i = 0; i < pn; i++) {
1957: /* add pdt[i,:]*AP into lnk */
1958: pnz = pdti[i + 1] - pdti[i];
1959: ptJ = pdtj + pdti[i];
1960: for (j = 0; j < pnz; j++) {
1961: row = ptJ[j]; /* row of AP == col of Pt */
1962: anz = ai[row + 1] - ai[row];
1963: Jptr = aj + ai[row];
1964: /* add non-zero cols of AP into the sorted linked list lnk */
1965: PetscCall(PetscLLCondensedAddSorted_Scalable(anz, Jptr, lnk));
1966: }
1968: /* add received col data into lnk */
1969: for (k = 0; k < merge->nrecv; k++) { /* k-th received message */
1970: if (i == *nextrow[k]) { /* i-th row */
1971: nzi = *(nextci[k] + 1) - *nextci[k];
1972: Jptr = buf_rj[k] + *nextci[k];
1973: PetscCall(PetscLLCondensedAddSorted_Scalable(nzi, Jptr, lnk));
1974: nextrow[k]++;
1975: nextci[k]++;
1976: }
1977: }
1979: /* add missing diagonal entry */
1980: if (C->force_diagonals) {
1981: k = i + owners[rank]; /* column index */
1982: PetscCall(PetscLLCondensedAddSorted_Scalable(1, &k, lnk));
1983: }
1985: nnz = lnk[0];
1987: /* if free space is not available, make more free space */
1988: if (current_space->local_remaining < nnz) {
1989: PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(nnz, current_space->total_array_size), ¤t_space));
1990: nspacedouble++;
1991: }
1992: /* copy data into free space, then initialize lnk */
1993: PetscCall(PetscLLCondensedClean_Scalable(nnz, current_space->array, lnk));
1994: PetscCall(MatPreallocateSet(i + owners[rank], nnz, current_space->array, dnz, onz));
1996: current_space->array += nnz;
1997: current_space->local_used += nnz;
1998: current_space->local_remaining -= nnz;
2000: bi[i + 1] = bi[i] + nnz;
2001: if (nnz > rmax) rmax = nnz;
2002: }
2003: PetscCall(PetscFree3(buf_ri_k, nextrow, nextci));
2005: PetscCall(PetscMalloc1(bi[pn] + 1, &bj));
2006: PetscCall(PetscFreeSpaceContiguous(&free_space, bj));
2007: afill_tmp = (PetscReal)bi[pn] / (pdti[pn] + poti[pon] + ai[am] + 1);
2008: if (afill_tmp > afill) afill = afill_tmp;
2009: PetscCall(PetscLLCondensedDestroy_Scalable(lnk));
2010: PetscCall(PetscHMapIDestroy(&ta));
2011: PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(p->A, &pdti, &pdtj));
2012: PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(p->B, &poti, &potj));
2014: /* create symbolic parallel matrix C - why cannot be assembled in Numeric part */
2015: PetscCall(MatSetSizes(C, pn, A->cmap->n, PETSC_DETERMINE, PETSC_DETERMINE));
2016: PetscCall(MatSetBlockSizes(C, PetscAbs(P->cmap->bs), PetscAbs(A->cmap->bs)));
2017: PetscCall(MatGetType(A, &mtype));
2018: PetscCall(MatSetType(C, mtype));
2019: PetscCall(MatMPIAIJSetPreallocation(C, 0, dnz, 0, onz));
2020: MatPreallocateEnd(dnz, onz);
2021: PetscCall(MatSetBlockSize(C, 1));
2022: PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2023: for (i = 0; i < pn; i++) {
2024: row = i + rstart;
2025: nnz = bi[i + 1] - bi[i];
2026: Jptr = bj + bi[i];
2027: PetscCall(MatSetValues(C, 1, &row, nnz, Jptr, NULL, INSERT_VALUES));
2028: }
2029: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2030: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
2031: PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
2032: merge->bi = bi;
2033: merge->bj = bj;
2034: merge->coi = coi;
2035: merge->coj = coj;
2036: merge->buf_ri = buf_ri;
2037: merge->buf_rj = buf_rj;
2038: merge->owners_co = owners_co;
2040: /* attach the supporting struct to C for reuse */
2041: C->product->data = ap;
2042: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2043: ap->merge = merge;
2045: C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2047: #if defined(PETSC_USE_INFO)
2048: if (bi[pn] != 0) {
2049: PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
2050: PetscCall(PetscInfo(C, "Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n", (double)afill));
2051: } else {
2052: PetscCall(PetscInfo(C, "Empty matrix product\n"));
2053: }
2054: #endif
2055: PetscFunctionReturn(PETSC_SUCCESS);
2056: }
2058: static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C)
2059: {
2060: Mat_Product *product = C->product;
2061: Mat A = product->A, B = product->B;
2062: PetscReal fill = product->fill;
2063: PetscBool flg;
2065: PetscFunctionBegin;
2066: /* scalable */
2067: PetscCall(PetscStrcmp(product->alg, "scalable", &flg));
2068: if (flg) {
2069: PetscCall(MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A, B, fill, C));
2070: goto next;
2071: }
2073: /* nonscalable */
2074: PetscCall(PetscStrcmp(product->alg, "nonscalable", &flg));
2075: if (flg) {
2076: PetscCall(MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A, B, fill, C));
2077: goto next;
2078: }
2080: /* matmatmult */
2081: PetscCall(PetscStrcmp(product->alg, "at*b", &flg));
2082: if (flg) {
2083: Mat At;
2084: Mat_APMPI *ptap;
2086: PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &At));
2087: PetscCall(MatMatMultSymbolic_MPIAIJ_MPIAIJ(At, B, fill, C));
2088: ptap = (Mat_APMPI *)C->product->data;
2089: if (ptap) {
2090: ptap->Pt = At;
2091: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2092: }
2093: C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
2094: goto next;
2095: }
2097: /* backend general code */
2098: PetscCall(PetscStrcmp(product->alg, "backend", &flg));
2099: if (flg) {
2100: PetscCall(MatProductSymbolic_MPIAIJBACKEND(C));
2101: PetscFunctionReturn(PETSC_SUCCESS);
2102: }
2104: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProduct type is not supported");
2106: next:
2107: C->ops->productnumeric = MatProductNumeric_AtB;
2108: PetscFunctionReturn(PETSC_SUCCESS);
2109: }
2111: /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */
2112: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C)
2113: {
2114: Mat_Product *product = C->product;
2115: Mat A = product->A, B = product->B;
2116: #if defined(PETSC_HAVE_HYPRE)
2117: const char *algTypes[5] = {"scalable", "nonscalable", "seqmpi", "backend", "hypre"};
2118: PetscInt nalg = 5;
2119: #else
2120: const char *algTypes[4] = {
2121: "scalable",
2122: "nonscalable",
2123: "seqmpi",
2124: "backend",
2125: };
2126: PetscInt nalg = 4;
2127: #endif
2128: PetscInt alg = 1; /* set nonscalable algorithm as default */
2129: PetscBool flg;
2130: MPI_Comm comm;
2132: PetscFunctionBegin;
2133: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
2135: /* Set "nonscalable" as default algorithm */
2136: PetscCall(PetscStrcmp(C->product->alg, "default", &flg));
2137: if (flg) {
2138: PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2140: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2141: if (B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2142: MatInfo Ainfo, Binfo;
2143: PetscInt nz_local;
2144: PetscBool alg_scalable_loc = PETSC_FALSE, alg_scalable;
2146: PetscCall(MatGetInfo(A, MAT_LOCAL, &Ainfo));
2147: PetscCall(MatGetInfo(B, MAT_LOCAL, &Binfo));
2148: nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2150: if (B->cmap->N > product->fill * nz_local) alg_scalable_loc = PETSC_TRUE;
2151: PetscCall(MPIU_Allreduce(&alg_scalable_loc, &alg_scalable, 1, MPIU_BOOL, MPI_LOR, comm));
2153: if (alg_scalable) {
2154: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2155: PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2156: PetscCall(PetscInfo(B, "Use scalable algorithm, BN %" PetscInt_FMT ", fill*nz_allocated %g\n", B->cmap->N, (double)(product->fill * nz_local)));
2157: }
2158: }
2159: }
2161: /* Get runtime option */
2162: if (product->api_user) {
2163: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatMatMult", "Mat");
2164: PetscCall(PetscOptionsEList("-matmatmult_via", "Algorithmic approach", "MatMatMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2165: PetscOptionsEnd();
2166: } else {
2167: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_AB", "Mat");
2168: PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatMatMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2169: PetscOptionsEnd();
2170: }
2171: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2173: C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ;
2174: PetscFunctionReturn(PETSC_SUCCESS);
2175: }
2177: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABt(Mat C)
2178: {
2179: PetscFunctionBegin;
2180: PetscCall(MatProductSetFromOptions_MPIAIJ_AB(C));
2181: C->ops->productsymbolic = MatProductSymbolic_ABt_MPIAIJ_MPIAIJ;
2182: PetscFunctionReturn(PETSC_SUCCESS);
2183: }
2185: /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */
2186: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C)
2187: {
2188: Mat_Product *product = C->product;
2189: Mat A = product->A, B = product->B;
2190: const char *algTypes[4] = {"scalable", "nonscalable", "at*b", "backend"};
2191: PetscInt nalg = 4;
2192: PetscInt alg = 1; /* set default algorithm */
2193: PetscBool flg;
2194: MPI_Comm comm;
2196: PetscFunctionBegin;
2197: /* Check matrix local sizes */
2198: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
2199: PetscCheck(A->rmap->rstart == B->rmap->rstart && A->rmap->rend == B->rmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, A (%" PetscInt_FMT ", %" PetscInt_FMT ") != B (%" PetscInt_FMT ",%" PetscInt_FMT ")",
2200: A->rmap->rstart, A->rmap->rend, B->rmap->rstart, B->rmap->rend);
2202: /* Set default algorithm */
2203: PetscCall(PetscStrcmp(C->product->alg, "default", &flg));
2204: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2206: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2207: if (alg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2208: MatInfo Ainfo, Binfo;
2209: PetscInt nz_local;
2210: PetscBool alg_scalable_loc = PETSC_FALSE, alg_scalable;
2212: PetscCall(MatGetInfo(A, MAT_LOCAL, &Ainfo));
2213: PetscCall(MatGetInfo(B, MAT_LOCAL, &Binfo));
2214: nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2216: if (B->cmap->N > product->fill * nz_local) alg_scalable_loc = PETSC_TRUE;
2217: PetscCall(MPIU_Allreduce(&alg_scalable_loc, &alg_scalable, 1, MPIU_BOOL, MPI_LOR, comm));
2219: if (alg_scalable) {
2220: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2221: PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2222: PetscCall(PetscInfo(B, "Use scalable algorithm, BN %" PetscInt_FMT ", fill*nz_allocated %g\n", B->cmap->N, (double)(product->fill * nz_local)));
2223: }
2224: }
2226: /* Get runtime option */
2227: if (product->api_user) {
2228: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatTransposeMatMult", "Mat");
2229: PetscCall(PetscOptionsEList("-mattransposematmult_via", "Algorithmic approach", "MatTransposeMatMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2230: PetscOptionsEnd();
2231: } else {
2232: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_AtB", "Mat");
2233: PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatTransposeMatMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2234: PetscOptionsEnd();
2235: }
2236: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2238: C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ;
2239: PetscFunctionReturn(PETSC_SUCCESS);
2240: }
2242: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C)
2243: {
2244: Mat_Product *product = C->product;
2245: Mat A = product->A, P = product->B;
2246: MPI_Comm comm;
2247: PetscBool flg;
2248: PetscInt alg = 1; /* set default algorithm */
2249: #if !defined(PETSC_HAVE_HYPRE)
2250: const char *algTypes[5] = {"scalable", "nonscalable", "allatonce", "allatonce_merged", "backend"};
2251: PetscInt nalg = 5;
2252: #else
2253: const char *algTypes[6] = {"scalable", "nonscalable", "allatonce", "allatonce_merged", "backend", "hypre"};
2254: PetscInt nalg = 6;
2255: #endif
2256: PetscInt pN = P->cmap->N;
2258: PetscFunctionBegin;
2259: /* Check matrix local sizes */
2260: PetscCall(PetscObjectGetComm((PetscObject)C, &comm));
2261: PetscCheck(A->rmap->rstart == P->rmap->rstart && A->rmap->rend == P->rmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, Arow (%" PetscInt_FMT ", %" PetscInt_FMT ") != Prow (%" PetscInt_FMT ",%" PetscInt_FMT ")",
2262: A->rmap->rstart, A->rmap->rend, P->rmap->rstart, P->rmap->rend);
2263: PetscCheck(A->cmap->rstart == P->rmap->rstart && A->cmap->rend == P->rmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, Acol (%" PetscInt_FMT ", %" PetscInt_FMT ") != Prow (%" PetscInt_FMT ",%" PetscInt_FMT ")",
2264: A->cmap->rstart, A->cmap->rend, P->rmap->rstart, P->rmap->rend);
2266: /* Set "nonscalable" as default algorithm */
2267: PetscCall(PetscStrcmp(C->product->alg, "default", &flg));
2268: if (flg) {
2269: PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2271: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2272: if (pN > 100000) {
2273: MatInfo Ainfo, Pinfo;
2274: PetscInt nz_local;
2275: PetscBool alg_scalable_loc = PETSC_FALSE, alg_scalable;
2277: PetscCall(MatGetInfo(A, MAT_LOCAL, &Ainfo));
2278: PetscCall(MatGetInfo(P, MAT_LOCAL, &Pinfo));
2279: nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
2281: if (pN > product->fill * nz_local) alg_scalable_loc = PETSC_TRUE;
2282: PetscCall(MPIU_Allreduce(&alg_scalable_loc, &alg_scalable, 1, MPIU_BOOL, MPI_LOR, comm));
2284: if (alg_scalable) {
2285: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2286: PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2287: }
2288: }
2289: }
2291: /* Get runtime option */
2292: if (product->api_user) {
2293: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatPtAP", "Mat");
2294: PetscCall(PetscOptionsEList("-matptap_via", "Algorithmic approach", "MatPtAP", algTypes, nalg, algTypes[alg], &alg, &flg));
2295: PetscOptionsEnd();
2296: } else {
2297: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_PtAP", "Mat");
2298: PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatPtAP", algTypes, nalg, algTypes[alg], &alg, &flg));
2299: PetscOptionsEnd();
2300: }
2301: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2303: C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ;
2304: PetscFunctionReturn(PETSC_SUCCESS);
2305: }
2307: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C)
2308: {
2309: Mat_Product *product = C->product;
2310: Mat A = product->A, R = product->B;
2312: PetscFunctionBegin;
2313: /* Check matrix local sizes */
2314: PetscCheck(A->cmap->n == R->cmap->n && A->rmap->n == R->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrix local dimensions are incompatible, A local (%" PetscInt_FMT ", %" PetscInt_FMT "), R local (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->n,
2315: A->rmap->n, R->rmap->n, R->cmap->n);
2317: C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ;
2318: PetscFunctionReturn(PETSC_SUCCESS);
2319: }
2321: /*
2322: Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm
2323: */
2324: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C)
2325: {
2326: Mat_Product *product = C->product;
2327: PetscBool flg = PETSC_FALSE;
2328: PetscInt alg = 1; /* default algorithm */
2329: const char *algTypes[3] = {"scalable", "nonscalable", "seqmpi"};
2330: PetscInt nalg = 3;
2332: PetscFunctionBegin;
2333: /* Set default algorithm */
2334: PetscCall(PetscStrcmp(C->product->alg, "default", &flg));
2335: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2337: /* Get runtime option */
2338: if (product->api_user) {
2339: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatMatMatMult", "Mat");
2340: PetscCall(PetscOptionsEList("-matmatmatmult_via", "Algorithmic approach", "MatMatMatMult", algTypes, nalg, algTypes[alg], &alg, &flg));
2341: PetscOptionsEnd();
2342: } else {
2343: PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_ABC", "Mat");
2344: PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatProduct_ABC", algTypes, nalg, algTypes[alg], &alg, &flg));
2345: PetscOptionsEnd();
2346: }
2347: if (flg) PetscCall(MatProductSetAlgorithm(C, (MatProductAlgorithm)algTypes[alg]));
2349: C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ;
2350: C->ops->productsymbolic = MatProductSymbolic_ABC;
2351: PetscFunctionReturn(PETSC_SUCCESS);
2352: }
2354: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C)
2355: {
2356: Mat_Product *product = C->product;
2358: PetscFunctionBegin;
2359: switch (product->type) {
2360: case MATPRODUCT_AB:
2361: PetscCall(MatProductSetFromOptions_MPIAIJ_AB(C));
2362: break;
2363: case MATPRODUCT_ABt:
2364: PetscCall(MatProductSetFromOptions_MPIAIJ_ABt(C));
2365: break;
2366: case MATPRODUCT_AtB:
2367: PetscCall(MatProductSetFromOptions_MPIAIJ_AtB(C));
2368: break;
2369: case MATPRODUCT_PtAP:
2370: PetscCall(MatProductSetFromOptions_MPIAIJ_PtAP(C));
2371: break;
2372: case MATPRODUCT_RARt:
2373: PetscCall(MatProductSetFromOptions_MPIAIJ_RARt(C));
2374: break;
2375: case MATPRODUCT_ABC:
2376: PetscCall(MatProductSetFromOptions_MPIAIJ_ABC(C));
2377: break;
2378: default:
2379: break;
2380: }
2381: PetscFunctionReturn(PETSC_SUCCESS);
2382: }