xref: /petsc/src/mat/impls/composite/mcomposite.c (revision 2ef1f0ff6e3530e8731eb06ad663081f5844f49f)
1793850ffSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
3793850ffSBarry Smith 
4f4259b30SLisandro Dalcin const char *const MatCompositeMergeTypes[] = {"left", "right", "MatCompositeMergeType", "MAT_COMPOSITE_", NULL};
58c8613bfSJakub Kruzik 
6793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink;
7793850ffSBarry Smith struct _Mat_CompositeLink {
8793850ffSBarry Smith   Mat               mat;
96d7c1e57SBarry Smith   Vec               work;
106d7c1e57SBarry Smith   Mat_CompositeLink next, prev;
11793850ffSBarry Smith };
12793850ffSBarry Smith 
13793850ffSBarry Smith typedef struct {
146d7c1e57SBarry Smith   MatCompositeType      type;
156d7c1e57SBarry Smith   Mat_CompositeLink     head, tail;
16793850ffSBarry Smith   Vec                   work;
17e4fc5df0SBarry Smith   PetscScalar           scale;                                      /* scale factor supplied with MatScale() */
18e4fc5df0SBarry Smith   Vec                   left, right;                                /* left and right diagonal scaling provided with MatDiagonalScale() */
1903049c21SJunchao Zhang   Vec                   leftwork, rightwork, leftwork2, rightwork2; /* Two pairs of working vectors */
206a7ede75SJakub Kruzik   PetscInt              nmat;
214b2558d6SJakub Kruzik   PetscBool             merge;
228c8613bfSJakub Kruzik   MatCompositeMergeType mergetype;
233b35acafSJakub Kruzik   MatStructure          structure;
2403049c21SJunchao Zhang 
2503049c21SJunchao Zhang   PetscScalar *scalings;
2603049c21SJunchao Zhang   PetscBool    merge_mvctx; /* Whether need to merge mvctx of component matrices */
2703049c21SJunchao Zhang   Vec         *lvecs;       /* [nmat] Basically, they are Mvctx->lvec of each component matrix */
2803049c21SJunchao Zhang   PetscScalar *larray;      /* [len] Data arrays of lvecs[] are stored consecutively in larray */
2903049c21SJunchao Zhang   PetscInt     len;         /* Length of larray[] */
3003049c21SJunchao Zhang   Vec          gvec;        /* Union of lvecs[] without duplicated entries */
3103049c21SJunchao Zhang   PetscInt    *location;    /* A map that maps entries in garray[] to larray[] */
3203049c21SJunchao Zhang   VecScatter   Mvctx;
33793850ffSBarry Smith } Mat_Composite;
34793850ffSBarry Smith 
35d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_Composite(Mat mat)
36d71ae5a4SJacob Faibussowitsch {
372c33bbaeSSatish Balay   Mat_Composite    *shell = (Mat_Composite *)mat->data;
386d7c1e57SBarry Smith   Mat_CompositeLink next  = shell->head, oldnext;
3903049c21SJunchao Zhang   PetscInt          i;
40793850ffSBarry Smith 
41793850ffSBarry Smith   PetscFunctionBegin;
42793850ffSBarry Smith   while (next) {
439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&next->mat));
4448a46eb9SPierre Jolivet     if (next->work && (!next->next || next->work != next->next->work)) PetscCall(VecDestroy(&next->work));
456d7c1e57SBarry Smith     oldnext = next;
46793850ffSBarry Smith     next    = next->next;
479566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldnext));
48793850ffSBarry Smith   }
499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->work));
509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->left));
519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->right));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->leftwork));
539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->rightwork));
549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->leftwork2));
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->rightwork2));
5603049c21SJunchao Zhang 
5703049c21SJunchao Zhang   if (shell->Mvctx) {
589566063dSJacob Faibussowitsch     for (i = 0; i < shell->nmat; i++) PetscCall(VecDestroy(&shell->lvecs[i]));
599566063dSJacob Faibussowitsch     PetscCall(PetscFree3(shell->location, shell->larray, shell->lvecs));
609566063dSJacob Faibussowitsch     PetscCall(PetscFree(shell->larray));
619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&shell->gvec));
629566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&shell->Mvctx));
6303049c21SJunchao Zhang   }
6403049c21SJunchao Zhang 
659566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell->scalings));
662e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeAddMat_C", NULL));
672e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetType_C", NULL));
682e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetType_C", NULL));
692e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMergeType_C", NULL));
702e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMatStructure_C", NULL));
712e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMatStructure_C", NULL));
722e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeMerge_C", NULL));
732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetNumberMat_C", NULL));
742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMat_C", NULL));
752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetScalings_C", NULL));
769566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78793850ffSBarry Smith }
79793850ffSBarry Smith 
80d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Composite_Multiplicative(Mat A, Vec x, Vec y)
81d71ae5a4SJacob Faibussowitsch {
826d7c1e57SBarry Smith   Mat_Composite    *shell = (Mat_Composite *)A->data;
836d7c1e57SBarry Smith   Mat_CompositeLink next  = shell->head;
846d7c1e57SBarry Smith   Vec               in, out;
8503049c21SJunchao Zhang   PetscScalar       scale;
8603049c21SJunchao Zhang   PetscInt          i;
876d7c1e57SBarry Smith 
886d7c1e57SBarry Smith   PetscFunctionBegin;
8928b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
906d7c1e57SBarry Smith   in = x;
91e4fc5df0SBarry Smith   if (shell->right) {
9248a46eb9SPierre Jolivet     if (!shell->rightwork) PetscCall(VecDuplicate(shell->right, &shell->rightwork));
939566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->rightwork, shell->right, in));
94e4fc5df0SBarry Smith     in = shell->rightwork;
95e4fc5df0SBarry Smith   }
966d7c1e57SBarry Smith   while (next->next) {
976d7c1e57SBarry Smith     if (!next->work) { /* should reuse previous work if the same size */
989566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(next->mat, NULL, &next->work));
996d7c1e57SBarry Smith     }
1006d7c1e57SBarry Smith     out = next->work;
1019566063dSJacob Faibussowitsch     PetscCall(MatMult(next->mat, in, out));
1026d7c1e57SBarry Smith     in   = out;
1036d7c1e57SBarry Smith     next = next->next;
1046d7c1e57SBarry Smith   }
1059566063dSJacob Faibussowitsch   PetscCall(MatMult(next->mat, in, y));
1061baa6e33SBarry Smith   if (shell->left) PetscCall(VecPointwiseMult(y, shell->left, y));
10703049c21SJunchao Zhang   scale = shell->scale;
1089371c9d4SSatish Balay   if (shell->scalings) {
1099371c9d4SSatish Balay     for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
1109371c9d4SSatish Balay   }
1119566063dSJacob Faibussowitsch   PetscCall(VecScale(y, scale));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1136d7c1e57SBarry Smith }
1146d7c1e57SBarry Smith 
115d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A, Vec x, Vec y)
116d71ae5a4SJacob Faibussowitsch {
1176d7c1e57SBarry Smith   Mat_Composite    *shell = (Mat_Composite *)A->data;
1186d7c1e57SBarry Smith   Mat_CompositeLink tail  = shell->tail;
1196d7c1e57SBarry Smith   Vec               in, out;
12003049c21SJunchao Zhang   PetscScalar       scale;
12103049c21SJunchao Zhang   PetscInt          i;
1226d7c1e57SBarry Smith 
1236d7c1e57SBarry Smith   PetscFunctionBegin;
12428b400f6SJacob Faibussowitsch   PetscCheck(tail, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
1256d7c1e57SBarry Smith   in = x;
126e4fc5df0SBarry Smith   if (shell->left) {
12748a46eb9SPierre Jolivet     if (!shell->leftwork) PetscCall(VecDuplicate(shell->left, &shell->leftwork));
1289566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->leftwork, shell->left, in));
129e4fc5df0SBarry Smith     in = shell->leftwork;
130e4fc5df0SBarry Smith   }
1316d7c1e57SBarry Smith   while (tail->prev) {
1326d7c1e57SBarry Smith     if (!tail->prev->work) { /* should reuse previous work if the same size */
1339566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(tail->mat, NULL, &tail->prev->work));
1346d7c1e57SBarry Smith     }
1356d7c1e57SBarry Smith     out = tail->prev->work;
1369566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(tail->mat, in, out));
1376d7c1e57SBarry Smith     in   = out;
1386d7c1e57SBarry Smith     tail = tail->prev;
1396d7c1e57SBarry Smith   }
1409566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tail->mat, in, y));
1411baa6e33SBarry Smith   if (shell->right) PetscCall(VecPointwiseMult(y, shell->right, y));
14203049c21SJunchao Zhang 
14303049c21SJunchao Zhang   scale = shell->scale;
1449371c9d4SSatish Balay   if (shell->scalings) {
1459371c9d4SSatish Balay     for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
1469371c9d4SSatish Balay   }
1479566063dSJacob Faibussowitsch   PetscCall(VecScale(y, scale));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1496d7c1e57SBarry Smith }
1506d7c1e57SBarry Smith 
151d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Composite(Mat mat, Vec x, Vec y)
152d71ae5a4SJacob Faibussowitsch {
15303049c21SJunchao Zhang   Mat_Composite     *shell = (Mat_Composite *)mat->data;
15403049c21SJunchao Zhang   Mat_CompositeLink  cur   = shell->head;
15503049c21SJunchao Zhang   Vec                in, y2, xin;
15603049c21SJunchao Zhang   Mat                A, B;
15703049c21SJunchao Zhang   PetscInt           i, j, k, n, nuniq, lo, hi, mid, *gindices, *buf, *tmp, tot;
15803049c21SJunchao Zhang   const PetscScalar *vals;
15903049c21SJunchao Zhang   const PetscInt    *garray;
16003049c21SJunchao Zhang   IS                 ix, iy;
16103049c21SJunchao Zhang   PetscBool          match;
162793850ffSBarry Smith 
163793850ffSBarry Smith   PetscFunctionBegin;
16428b400f6SJacob Faibussowitsch   PetscCheck(cur, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
165e4fc5df0SBarry Smith   in = x;
166e4fc5df0SBarry Smith   if (shell->right) {
16748a46eb9SPierre Jolivet     if (!shell->rightwork) PetscCall(VecDuplicate(shell->right, &shell->rightwork));
1689566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->rightwork, shell->right, in));
169e4fc5df0SBarry Smith     in = shell->rightwork;
170e4fc5df0SBarry Smith   }
17103049c21SJunchao Zhang 
17203049c21SJunchao Zhang   /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time
17303049c21SJunchao Zhang      we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and
17403049c21SJunchao Zhang      it is legal to merge Mvctx, because all component matrices have the same size.
17503049c21SJunchao Zhang    */
17603049c21SJunchao Zhang   if (shell->merge_mvctx && !shell->Mvctx) {
17703049c21SJunchao Zhang     /* Currently only implemented for MATMPIAIJ */
17803049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
1799566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat, MATMPIAIJ, &match));
18003049c21SJunchao Zhang       if (!match) {
18103049c21SJunchao Zhang         shell->merge_mvctx = PETSC_FALSE;
18203049c21SJunchao Zhang         goto skip_merge_mvctx;
183e4fc5df0SBarry Smith       }
184e4fc5df0SBarry Smith     }
18503049c21SJunchao Zhang 
18603049c21SJunchao Zhang     /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */
18703049c21SJunchao Zhang     tot = 0;
18803049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
1899566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, NULL));
1909566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
19103049c21SJunchao Zhang       tot += n;
19203049c21SJunchao Zhang     }
1939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(tot, &shell->location, tot, &shell->larray, shell->nmat, &shell->lvecs));
19403049c21SJunchao Zhang     shell->len = tot;
19503049c21SJunchao Zhang 
19603049c21SJunchao Zhang     /* Go through matrices second time to sort off-diag columns and remove dups */
1979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot, &gindices)); /* No Malloc2() since we will give one to petsc and free the other */
1989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot, &buf));
19903049c21SJunchao Zhang     nuniq = 0; /* Number of unique nonzero columns */
20003049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
2019566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray));
2029566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
20303049c21SJunchao Zhang       /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */
20403049c21SJunchao Zhang       i = j = k = 0;
20503049c21SJunchao Zhang       while (i < n && j < nuniq) {
20603049c21SJunchao Zhang         if (garray[i] < gindices[j]) buf[k++] = garray[i++];
20703049c21SJunchao Zhang         else if (garray[i] > gindices[j]) buf[k++] = gindices[j++];
2089371c9d4SSatish Balay         else {
2099371c9d4SSatish Balay           buf[k++] = garray[i++];
2109371c9d4SSatish Balay           j++;
2119371c9d4SSatish Balay         }
21203049c21SJunchao Zhang       }
21303049c21SJunchao Zhang       /* Copy leftover in garray[] or gindices[] */
21403049c21SJunchao Zhang       if (i < n) {
2159566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf + k, garray + i, n - i));
21603049c21SJunchao Zhang         nuniq = k + n - i;
21703049c21SJunchao Zhang       } else if (j < nuniq) {
2189566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf + k, gindices + j, nuniq - j));
21903049c21SJunchao Zhang         nuniq = k + nuniq - j;
22003049c21SJunchao Zhang       } else nuniq = k;
22103049c21SJunchao Zhang       /* Swap gindices and buf to merge garray of the next matrix */
22203049c21SJunchao Zhang       tmp      = gindices;
22303049c21SJunchao Zhang       gindices = buf;
22403049c21SJunchao Zhang       buf      = tmp;
22503049c21SJunchao Zhang     }
2269566063dSJacob Faibussowitsch     PetscCall(PetscFree(buf));
22703049c21SJunchao Zhang 
22803049c21SJunchao Zhang     /* Go through matrices third time to build a map from gindices[] to garray[] */
22903049c21SJunchao Zhang     tot = 0;
23003049c21SJunchao Zhang     for (cur = shell->head, j = 0; cur; cur = cur->next, j++) { /* j-th matrix */
2319566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray));
2329566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
2339566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, n, NULL, &shell->lvecs[j]));
23403049c21SJunchao Zhang       /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */
23503049c21SJunchao Zhang       lo = 0;
23603049c21SJunchao Zhang       for (i = 0; i < n; i++) {
23703049c21SJunchao Zhang         hi = nuniq;
23803049c21SJunchao Zhang         while (hi - lo > 1) {
23903049c21SJunchao Zhang           mid = lo + (hi - lo) / 2;
24003049c21SJunchao Zhang           if (garray[i] < gindices[mid]) hi = mid;
24103049c21SJunchao Zhang           else lo = mid;
24203049c21SJunchao Zhang         }
24303049c21SJunchao Zhang         shell->location[tot + i] = lo; /* gindices[lo] = garray[i] */
24403049c21SJunchao Zhang         lo++;                          /* Since garray[i+1] > garray[i], we can safely advance lo */
24503049c21SJunchao Zhang       }
24603049c21SJunchao Zhang       tot += n;
24703049c21SJunchao Zhang     }
24803049c21SJunchao Zhang 
24903049c21SJunchao Zhang     /* Build merged Mvctx */
2509566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nuniq, gindices, PETSC_OWN_POINTER, &ix));
2519566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, nuniq, 0, 1, &iy));
2529566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &xin));
2539566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, nuniq, &shell->gvec));
2549566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(xin, ix, shell->gvec, iy, &shell->Mvctx));
2559566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xin));
2569566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ix));
2579566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iy));
25803049c21SJunchao Zhang   }
25903049c21SJunchao Zhang 
26003049c21SJunchao Zhang skip_merge_mvctx:
2619566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0));
2629566063dSJacob Faibussowitsch   if (!shell->leftwork2) PetscCall(VecDuplicate(y, &shell->leftwork2));
26303049c21SJunchao Zhang   y2 = shell->leftwork2;
26403049c21SJunchao Zhang 
26503049c21SJunchao Zhang   if (shell->Mvctx) { /* Have a merged Mvctx */
26603049c21SJunchao Zhang     /* Suppose we want to compute y = sMx, where s is the scaling factor and A, B are matrix M's diagonal/off-diagonal part. We could do
267da81f932SPierre Jolivet        in y = s(Ax1 + Bx2) or y = sAx1 + sBx2. The former incurs less FLOPS than the latter, but the latter provides an opportunity to
26803049c21SJunchao Zhang        overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach.
26903049c21SJunchao Zhang      */
2709566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(shell->Mvctx, in, shell->gvec, INSERT_VALUES, SCATTER_FORWARD));
2719566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(shell->Mvctx, in, shell->gvec, INSERT_VALUES, SCATTER_FORWARD));
27203049c21SJunchao Zhang 
2739566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(shell->gvec, &vals));
27403049c21SJunchao Zhang     for (i = 0; i < shell->len; i++) shell->larray[i] = vals[shell->location[i]];
2759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(shell->gvec, &vals));
27603049c21SJunchao Zhang 
27703049c21SJunchao Zhang     for (cur = shell->head, tot = i = 0; cur; cur = cur->next, i++) { /* i-th matrix */
2789566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, &A, &B, NULL));
279dbbe0bcdSBarry Smith       PetscUseTypeMethod(A, mult, in, y2);
2809566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
2819566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(shell->lvecs[i], &shell->larray[tot]));
2829566063dSJacob Faibussowitsch       PetscCall((*B->ops->multadd)(B, shell->lvecs[i], y2, y2));
2839566063dSJacob Faibussowitsch       PetscCall(VecResetArray(shell->lvecs[i]));
2849566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, (shell->scalings ? shell->scalings[i] : 1.0), y2));
28503049c21SJunchao Zhang       tot += n;
28603049c21SJunchao Zhang     }
28703049c21SJunchao Zhang   } else {
28803049c21SJunchao Zhang     if (shell->scalings) {
28903049c21SJunchao Zhang       for (cur = shell->head, i = 0; cur; cur = cur->next, i++) {
2909566063dSJacob Faibussowitsch         PetscCall(MatMult(cur->mat, in, y2));
2919566063dSJacob Faibussowitsch         PetscCall(VecAXPY(y, shell->scalings[i], y2));
29203049c21SJunchao Zhang       }
29303049c21SJunchao Zhang     } else {
2949566063dSJacob Faibussowitsch       for (cur = shell->head; cur; cur = cur->next) PetscCall(MatMultAdd(cur->mat, in, y, y));
29503049c21SJunchao Zhang     }
29603049c21SJunchao Zhang   }
29703049c21SJunchao Zhang 
2989566063dSJacob Faibussowitsch   if (shell->left) PetscCall(VecPointwiseMult(y, shell->left, y));
2999566063dSJacob Faibussowitsch   PetscCall(VecScale(y, shell->scale));
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301793850ffSBarry Smith }
302793850ffSBarry Smith 
303d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_Composite(Mat A, Vec x, Vec y)
304d71ae5a4SJacob Faibussowitsch {
305793850ffSBarry Smith   Mat_Composite    *shell = (Mat_Composite *)A->data;
306793850ffSBarry Smith   Mat_CompositeLink next  = shell->head;
30703049c21SJunchao Zhang   Vec               in, y2 = NULL;
30803049c21SJunchao Zhang   PetscInt          i;
309793850ffSBarry Smith 
310793850ffSBarry Smith   PetscFunctionBegin;
31128b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
312e4fc5df0SBarry Smith   in = x;
313e4fc5df0SBarry Smith   if (shell->left) {
31448a46eb9SPierre Jolivet     if (!shell->leftwork) PetscCall(VecDuplicate(shell->left, &shell->leftwork));
3159566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->leftwork, shell->left, in));
316e4fc5df0SBarry Smith     in = shell->leftwork;
317e4fc5df0SBarry Smith   }
31803049c21SJunchao Zhang 
3199566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(next->mat, in, y));
32003049c21SJunchao Zhang   if (shell->scalings) {
3219566063dSJacob Faibussowitsch     PetscCall(VecScale(y, shell->scalings[0]));
3229566063dSJacob Faibussowitsch     if (!shell->rightwork2) PetscCall(VecDuplicate(y, &shell->rightwork2));
32303049c21SJunchao Zhang     y2 = shell->rightwork2;
32403049c21SJunchao Zhang   }
32503049c21SJunchao Zhang   i = 1;
326e4fc5df0SBarry Smith   while ((next = next->next)) {
3279566063dSJacob Faibussowitsch     if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat, in, y, y));
32803049c21SJunchao Zhang     else {
3299566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(next->mat, in, y2));
3309566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, shell->scalings[i++], y2));
33103049c21SJunchao Zhang     }
332e4fc5df0SBarry Smith   }
3331baa6e33SBarry Smith   if (shell->right) PetscCall(VecPointwiseMult(y, shell->right, y));
3349566063dSJacob Faibussowitsch   PetscCall(VecScale(y, shell->scale));
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
336793850ffSBarry Smith }
337793850ffSBarry Smith 
338d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_Composite(Mat A, Vec x, Vec y, Vec z)
339d71ae5a4SJacob Faibussowitsch {
3407bf3a71bSJakub Kruzik   Mat_Composite *shell = (Mat_Composite *)A->data;
3417bf3a71bSJakub Kruzik 
3427bf3a71bSJakub Kruzik   PetscFunctionBegin;
3437bf3a71bSJakub Kruzik   if (y != z) {
3449566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, z));
3459566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, 1.0, y));
3467bf3a71bSJakub Kruzik   } else {
34748a46eb9SPierre Jolivet     if (!shell->leftwork) PetscCall(VecDuplicate(z, &shell->leftwork));
3489566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, shell->leftwork));
3499566063dSJacob Faibussowitsch     PetscCall(VecCopy(y, z));
3509566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, 1.0, shell->leftwork));
3517bf3a71bSJakub Kruzik   }
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3537bf3a71bSJakub Kruzik }
3547bf3a71bSJakub Kruzik 
355d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_Composite(Mat A, Vec x, Vec y, Vec z)
356d71ae5a4SJacob Faibussowitsch {
3577bf3a71bSJakub Kruzik   Mat_Composite *shell = (Mat_Composite *)A->data;
3587bf3a71bSJakub Kruzik 
3597bf3a71bSJakub Kruzik   PetscFunctionBegin;
3607bf3a71bSJakub Kruzik   if (y != z) {
3619566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, x, z));
3629566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, 1.0, y));
3637bf3a71bSJakub Kruzik   } else {
36448a46eb9SPierre Jolivet     if (!shell->rightwork) PetscCall(VecDuplicate(z, &shell->rightwork));
3659566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, x, shell->rightwork));
3669566063dSJacob Faibussowitsch     PetscCall(VecCopy(y, z));
3679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z, 1.0, shell->rightwork));
3687bf3a71bSJakub Kruzik   }
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3707bf3a71bSJakub Kruzik }
3717bf3a71bSJakub Kruzik 
372d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_Composite(Mat A, Vec v)
373d71ae5a4SJacob Faibussowitsch {
374793850ffSBarry Smith   Mat_Composite    *shell = (Mat_Composite *)A->data;
375793850ffSBarry Smith   Mat_CompositeLink next  = shell->head;
37603049c21SJunchao Zhang   PetscInt          i;
377793850ffSBarry Smith 
378793850ffSBarry Smith   PetscFunctionBegin;
37928b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
38008401ef6SPierre Jolivet   PetscCheck(!shell->right && !shell->left, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get diagonal if left or right scaling");
381e4fc5df0SBarry Smith 
3829566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(next->mat, v));
3839566063dSJacob Faibussowitsch   if (shell->scalings) PetscCall(VecScale(v, shell->scalings[0]));
38403049c21SJunchao Zhang 
38548a46eb9SPierre Jolivet   if (next->next && !shell->work) PetscCall(VecDuplicate(v, &shell->work));
38603049c21SJunchao Zhang   i = 1;
387793850ffSBarry Smith   while ((next = next->next)) {
3889566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(next->mat, shell->work));
3899566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v, (shell->scalings ? shell->scalings[i++] : 1.0), shell->work));
390793850ffSBarry Smith   }
3919566063dSJacob Faibussowitsch   PetscCall(VecScale(v, shell->scale));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393793850ffSBarry Smith }
394793850ffSBarry Smith 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_Composite(Mat Y, MatAssemblyType t)
396d71ae5a4SJacob Faibussowitsch {
3974b2558d6SJakub Kruzik   Mat_Composite *shell = (Mat_Composite *)Y->data;
398b52f573bSBarry Smith 
399793850ffSBarry Smith   PetscFunctionBegin;
4001baa6e33SBarry Smith   if (shell->merge) PetscCall(MatCompositeMerge(Y));
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
402793850ffSBarry Smith }
403793850ffSBarry Smith 
404d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_Composite(Mat inA, PetscScalar alpha)
405d71ae5a4SJacob Faibussowitsch {
406e4fc5df0SBarry Smith   Mat_Composite *a = (Mat_Composite *)inA->data;
407e4fc5df0SBarry Smith 
408e4fc5df0SBarry Smith   PetscFunctionBegin;
409321429dbSBarry Smith   a->scale *= alpha;
4103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
411e4fc5df0SBarry Smith }
412e4fc5df0SBarry Smith 
413d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_Composite(Mat inA, Vec left, Vec right)
414d71ae5a4SJacob Faibussowitsch {
415e4fc5df0SBarry Smith   Mat_Composite *a = (Mat_Composite *)inA->data;
416e4fc5df0SBarry Smith 
417e4fc5df0SBarry Smith   PetscFunctionBegin;
418e4fc5df0SBarry Smith   if (left) {
419321429dbSBarry Smith     if (!a->left) {
4209566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
4219566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
422321429dbSBarry Smith     } else {
4239566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
424321429dbSBarry Smith     }
425e4fc5df0SBarry Smith   }
426e4fc5df0SBarry Smith   if (right) {
427321429dbSBarry Smith     if (!a->right) {
4289566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
4299566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
430321429dbSBarry Smith     } else {
4319566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
432321429dbSBarry Smith     }
433e4fc5df0SBarry Smith   }
4343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
435e4fc5df0SBarry Smith }
436793850ffSBarry Smith 
437d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetFromOptions_Composite(Mat A, PetscOptionItems *PetscOptionsObject)
438d71ae5a4SJacob Faibussowitsch {
4394b2558d6SJakub Kruzik   Mat_Composite *a = (Mat_Composite *)A->data;
4404b2558d6SJakub Kruzik 
4414b2558d6SJakub Kruzik   PetscFunctionBegin;
442d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MATCOMPOSITE options");
4439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge", "Merge at MatAssemblyEnd", "MatCompositeMerge", a->merge, &a->merge, NULL));
4449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-mat_composite_merge_type", "Set composite merge direction", "MatCompositeSetMergeType", MatCompositeMergeTypes, (PetscEnum)a->mergetype, (PetscEnum *)&a->mergetype, NULL));
4459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx", "Merge MatMult() vecscat contexts", "MatCreateComposite", a->merge_mvctx, &a->merge_mvctx, NULL));
446d0609cedSBarry Smith   PetscOptionsHeadEnd();
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4484b2558d6SJakub Kruzik }
4494b2558d6SJakub Kruzik 
4502c0821f3SBarry Smith /*@
4518bd739bdSJakub Kruzik    MatCreateComposite - Creates a matrix as the sum or product of one or more matrices
452793850ffSBarry Smith 
453d083f849SBarry Smith   Collective
454793850ffSBarry Smith 
455793850ffSBarry Smith    Input Parameters:
456793850ffSBarry Smith +  comm - MPI communicator
457793850ffSBarry Smith .  nmat - number of matrices to put in
458793850ffSBarry Smith -  mats - the matrices
459793850ffSBarry Smith 
460793850ffSBarry Smith    Output Parameter:
461db36af27SMatthew G. Knepley .  mat - the matrix
462793850ffSBarry Smith 
4634b2558d6SJakub Kruzik    Options Database Keys:
46411a5261eSBarry Smith +  -mat_composite_merge         - merge in `MatAssemblyEnd()`
46511a5261eSBarry Smith .  -mat_composite_merge_mvctx   - merge Mvctx of component matrices to optimize communication in `MatMult()` for ADDITIVE matrices
466b28d0daaSJakub Kruzik -  -mat_composite_merge_type    - set merge direction
4674b2558d6SJakub Kruzik 
468793850ffSBarry Smith    Level: advanced
469793850ffSBarry Smith 
47011a5261eSBarry Smith    Note:
471793850ffSBarry Smith      Alternative construction
472*2ef1f0ffSBarry Smith .vb
473*2ef1f0ffSBarry Smith        MatCreate(comm,&mat);
474*2ef1f0ffSBarry Smith        MatSetSizes(mat,m,n,M,N);
475*2ef1f0ffSBarry Smith        MatSetType(mat,MATCOMPOSITE);
476*2ef1f0ffSBarry Smith        MatCompositeAddMat(mat,mats[0]);
477*2ef1f0ffSBarry Smith        ....
478*2ef1f0ffSBarry Smith        MatCompositeAddMat(mat,mats[nmat-1]);
479*2ef1f0ffSBarry Smith        MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
480*2ef1f0ffSBarry Smith        MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
481*2ef1f0ffSBarry Smith .ve
482793850ffSBarry Smith 
4836d7c1e57SBarry Smith      For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0]
4846d7c1e57SBarry Smith 
485*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCompositeGetMat()`, `MatCompositeMerge()`, `MatCompositeSetType()`, `MATCOMPOSITE`
486793850ffSBarry Smith @*/
487d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateComposite(MPI_Comm comm, PetscInt nmat, const Mat *mats, Mat *mat)
488d71ae5a4SJacob Faibussowitsch {
489793850ffSBarry Smith   PetscInt m, n, M, N, i;
490793850ffSBarry Smith 
491793850ffSBarry Smith   PetscFunctionBegin;
49208401ef6SPierre Jolivet   PetscCheck(nmat >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in at least one matrix");
493064a246eSJacob Faibussowitsch   PetscValidPointer(mat, 4);
494793850ffSBarry Smith 
4959566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[0], PETSC_IGNORE, &n));
4969566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[nmat - 1], &m, PETSC_IGNORE));
4979566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[0], PETSC_IGNORE, &N));
4989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[nmat - 1], &M, PETSC_IGNORE));
4999566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, mat));
5009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mat, m, n, M, N));
5019566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat, MATCOMPOSITE));
50248a46eb9SPierre Jolivet   for (i = 0; i < nmat; i++) PetscCall(MatCompositeAddMat(*mat, mats[i]));
5039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
5049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
506793850ffSBarry Smith }
507793850ffSBarry Smith 
508d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeAddMat_Composite(Mat mat, Mat smat)
509d71ae5a4SJacob Faibussowitsch {
510d7f81bf2SJakub Kruzik   Mat_Composite    *shell = (Mat_Composite *)mat->data;
511d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink, next = shell->head;
512d7f81bf2SJakub Kruzik 
513d7f81bf2SJakub Kruzik   PetscFunctionBegin;
5144dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
515f4259b30SLisandro Dalcin   ilink->next = NULL;
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)smat));
517d7f81bf2SJakub Kruzik   ilink->mat = smat;
518d7f81bf2SJakub Kruzik 
519d7f81bf2SJakub Kruzik   if (!next) shell->head = ilink;
520d7f81bf2SJakub Kruzik   else {
521ad540459SPierre Jolivet     while (next->next) next = next->next;
522d7f81bf2SJakub Kruzik     next->next  = ilink;
523d7f81bf2SJakub Kruzik     ilink->prev = next;
524d7f81bf2SJakub Kruzik   }
525d7f81bf2SJakub Kruzik   shell->tail = ilink;
526d7f81bf2SJakub Kruzik   shell->nmat += 1;
52703049c21SJunchao Zhang 
52803049c21SJunchao Zhang   /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */
52903049c21SJunchao Zhang   if (shell->scalings) {
5309566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(sizeof(PetscScalar) * shell->nmat, &shell->scalings));
53103049c21SJunchao Zhang     shell->scalings[shell->nmat - 1] = 1.0;
53203049c21SJunchao Zhang   }
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
534d7f81bf2SJakub Kruzik }
535d7f81bf2SJakub Kruzik 
536793850ffSBarry Smith /*@
5378bd739bdSJakub Kruzik     MatCompositeAddMat - Add another matrix to a composite matrix.
538793850ffSBarry Smith 
539c3339decSBarry Smith    Collective
540793850ffSBarry Smith 
541793850ffSBarry Smith     Input Parameters:
542793850ffSBarry Smith +   mat - the composite matrix
543793850ffSBarry Smith -   smat - the partial matrix
544793850ffSBarry Smith 
545793850ffSBarry Smith    Level: advanced
546793850ffSBarry Smith 
547*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
548793850ffSBarry Smith @*/
549d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeAddMat(Mat mat, Mat smat)
550d71ae5a4SJacob Faibussowitsch {
551793850ffSBarry Smith   PetscFunctionBegin;
5520700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
5530700a824SBarry Smith   PetscValidHeaderSpecific(smat, MAT_CLASSID, 2);
554cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeAddMat_C", (Mat, Mat), (mat, smat));
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
556d7f81bf2SJakub Kruzik }
557793850ffSBarry Smith 
558d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetType_Composite(Mat mat, MatCompositeType type)
559d71ae5a4SJacob Faibussowitsch {
560d7f81bf2SJakub Kruzik   Mat_Composite *b = (Mat_Composite *)mat->data;
561d7f81bf2SJakub Kruzik 
562d7f81bf2SJakub Kruzik   PetscFunctionBegin;
563ced1ac25SJakub Kruzik   b->type = type;
564d7f81bf2SJakub Kruzik   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
565f4259b30SLisandro Dalcin     mat->ops->getdiagonal   = NULL;
566d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite_Multiplicative;
567d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative;
56803049c21SJunchao Zhang     b->merge_mvctx          = PETSC_FALSE;
569d7f81bf2SJakub Kruzik   } else {
570d7f81bf2SJakub Kruzik     mat->ops->getdiagonal   = MatGetDiagonal_Composite;
571d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite;
572d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite;
573793850ffSBarry Smith   }
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5756d7c1e57SBarry Smith }
5766d7c1e57SBarry Smith 
5772c0821f3SBarry Smith /*@
5788bd739bdSJakub Kruzik    MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product.
5796d7c1e57SBarry Smith 
580c3339decSBarry Smith    Logically Collective
5816d7c1e57SBarry Smith 
5826d7c1e57SBarry Smith    Input Parameters:
5836d7c1e57SBarry Smith .  mat - the composite matrix
5846d7c1e57SBarry Smith 
5856d7c1e57SBarry Smith    Level: advanced
5866d7c1e57SBarry Smith 
587*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE`
5886d7c1e57SBarry Smith @*/
589d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetType(Mat mat, MatCompositeType type)
590d71ae5a4SJacob Faibussowitsch {
5916d7c1e57SBarry Smith   PetscFunctionBegin;
592d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
593b2b245abSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
594cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetType_C", (Mat, MatCompositeType), (mat, type));
5953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
596793850ffSBarry Smith }
597793850ffSBarry Smith 
598d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetType_Composite(Mat mat, MatCompositeType *type)
599d71ae5a4SJacob Faibussowitsch {
6006dbc55e5SJakub Kruzik   Mat_Composite *b = (Mat_Composite *)mat->data;
6016dbc55e5SJakub Kruzik 
6026dbc55e5SJakub Kruzik   PetscFunctionBegin;
6036dbc55e5SJakub Kruzik   *type = b->type;
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6056dbc55e5SJakub Kruzik }
6066dbc55e5SJakub Kruzik 
6076dbc55e5SJakub Kruzik /*@
6086dbc55e5SJakub Kruzik    MatCompositeGetType - Returns type of composite.
6096dbc55e5SJakub Kruzik 
6106dbc55e5SJakub Kruzik    Not Collective
6116dbc55e5SJakub Kruzik 
6126dbc55e5SJakub Kruzik    Input Parameter:
6136dbc55e5SJakub Kruzik .  mat - the composite matrix
6146dbc55e5SJakub Kruzik 
6156dbc55e5SJakub Kruzik    Output Parameter:
6166dbc55e5SJakub Kruzik .  type - type of composite
6176dbc55e5SJakub Kruzik 
6186dbc55e5SJakub Kruzik    Level: advanced
6196dbc55e5SJakub Kruzik 
620*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE`, `MatCompositeType`
6216dbc55e5SJakub Kruzik @*/
622d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetType(Mat mat, MatCompositeType *type)
623d71ae5a4SJacob Faibussowitsch {
6246dbc55e5SJakub Kruzik   PetscFunctionBegin;
6256dbc55e5SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6266dbc55e5SJakub Kruzik   PetscValidPointer(type, 2);
627cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetType_C", (Mat, MatCompositeType *), (mat, type));
6283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6296dbc55e5SJakub Kruzik }
6306dbc55e5SJakub Kruzik 
631d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat, MatStructure str)
632d71ae5a4SJacob Faibussowitsch {
6333b35acafSJakub Kruzik   Mat_Composite *b = (Mat_Composite *)mat->data;
6343b35acafSJakub Kruzik 
6353b35acafSJakub Kruzik   PetscFunctionBegin;
6363b35acafSJakub Kruzik   b->structure = str;
6373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6383b35acafSJakub Kruzik }
6393b35acafSJakub Kruzik 
6403b35acafSJakub Kruzik /*@
6413b35acafSJakub Kruzik    MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix.
6423b35acafSJakub Kruzik 
6433b35acafSJakub Kruzik    Not Collective
6443b35acafSJakub Kruzik 
6453b35acafSJakub Kruzik    Input Parameters:
6463b35acafSJakub Kruzik +  mat - the composite matrix
64711a5261eSBarry Smith -  str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN` (default) or `SUBSET_NONZERO_PATTERN`
6483b35acafSJakub Kruzik 
6493b35acafSJakub Kruzik    Level: advanced
6503b35acafSJakub Kruzik 
65111a5261eSBarry Smith    Note:
65211a5261eSBarry Smith     Information about the matrices structure is used in `MatCompositeMerge()` for additive composite matrix.
6533b35acafSJakub Kruzik 
654*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE`
6553b35acafSJakub Kruzik @*/
656d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMatStructure(Mat mat, MatStructure str)
657d71ae5a4SJacob Faibussowitsch {
6583b35acafSJakub Kruzik   PetscFunctionBegin;
6593b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
660cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMatStructure_C", (Mat, MatStructure), (mat, str));
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6623b35acafSJakub Kruzik }
6633b35acafSJakub Kruzik 
664d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat, MatStructure *str)
665d71ae5a4SJacob Faibussowitsch {
6663b35acafSJakub Kruzik   Mat_Composite *b = (Mat_Composite *)mat->data;
6673b35acafSJakub Kruzik 
6683b35acafSJakub Kruzik   PetscFunctionBegin;
6693b35acafSJakub Kruzik   *str = b->structure;
6703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6713b35acafSJakub Kruzik }
6723b35acafSJakub Kruzik 
6733b35acafSJakub Kruzik /*@
6743b35acafSJakub Kruzik    MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix.
6753b35acafSJakub Kruzik 
6763b35acafSJakub Kruzik    Not Collective
6773b35acafSJakub Kruzik 
6783b35acafSJakub Kruzik    Input Parameter:
6793b35acafSJakub Kruzik .  mat - the composite matrix
6803b35acafSJakub Kruzik 
6813b35acafSJakub Kruzik    Output Parameter:
6823b35acafSJakub Kruzik .  str - structure of the matrices
6833b35acafSJakub Kruzik 
6843b35acafSJakub Kruzik    Level: advanced
6853b35acafSJakub Kruzik 
686*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE`
6873b35acafSJakub Kruzik @*/
688d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMatStructure(Mat mat, MatStructure *str)
689d71ae5a4SJacob Faibussowitsch {
6903b35acafSJakub Kruzik   PetscFunctionBegin;
6913b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6923b35acafSJakub Kruzik   PetscValidPointer(str, 2);
693cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMatStructure_C", (Mat, MatStructure *), (mat, str));
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6953b35acafSJakub Kruzik }
6963b35acafSJakub Kruzik 
697d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat, MatCompositeMergeType type)
698d71ae5a4SJacob Faibussowitsch {
6998c8613bfSJakub Kruzik   Mat_Composite *shell = (Mat_Composite *)mat->data;
7008c8613bfSJakub Kruzik 
7018c8613bfSJakub Kruzik   PetscFunctionBegin;
7028c8613bfSJakub Kruzik   shell->mergetype = type;
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7048c8613bfSJakub Kruzik }
7058c8613bfSJakub Kruzik 
7068c8613bfSJakub Kruzik /*@
70711a5261eSBarry Smith    MatCompositeSetMergeType - Sets order of `MatCompositeMerge()`.
7088c8613bfSJakub Kruzik 
709c3339decSBarry Smith    Logically Collective
7108c8613bfSJakub Kruzik 
711d8d19677SJose E. Roman    Input Parameters:
7128c8613bfSJakub Kruzik +  mat - the composite matrix
71311a5261eSBarry Smith -  type - `MAT_COMPOSITE_MERGE RIGHT` (default) to start merge from right with the first added matrix (mat[0]),
71411a5261eSBarry Smith           `MAT_COMPOSITE_MERGE_LEFT` to start merge from left with the last added matrix (mat[nmat-1])
7158c8613bfSJakub Kruzik 
7168c8613bfSJakub Kruzik    Level: advanced
7178c8613bfSJakub Kruzik 
71811a5261eSBarry Smith    Note:
719da81f932SPierre Jolivet     The resulting matrix is the same regardless of the `MatCompositeMergeType`. Only the order of operation is changed.
72011a5261eSBarry Smith     If set to `MAT_COMPOSITE_MERGE_RIGHT` the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0])))
7218c8613bfSJakub Kruzik     otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0].
7228c8613bfSJakub Kruzik 
723*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE`
7248c8613bfSJakub Kruzik @*/
725d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMergeType(Mat mat, MatCompositeMergeType type)
726d71ae5a4SJacob Faibussowitsch {
7278c8613bfSJakub Kruzik   PetscFunctionBegin;
7288c8613bfSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
7298c8613bfSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
730cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMergeType_C", (Mat, MatCompositeMergeType), (mat, type));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7328c8613bfSJakub Kruzik }
7338c8613bfSJakub Kruzik 
734d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeMerge_Composite(Mat mat)
735d71ae5a4SJacob Faibussowitsch {
736b52f573bSBarry Smith   Mat_Composite    *shell = (Mat_Composite *)mat->data;
7376d7c1e57SBarry Smith   Mat_CompositeLink next = shell->head, prev = shell->tail;
7386d7c1e57SBarry Smith   Mat               tmat, newmat;
7391ba60692SJed Brown   Vec               left, right;
7401ba60692SJed Brown   PetscScalar       scale;
74103049c21SJunchao Zhang   PetscInt          i;
742b52f573bSBarry Smith 
743b52f573bSBarry Smith   PetscFunctionBegin;
74428b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
74503049c21SJunchao Zhang   scale = shell->scale;
7466d7c1e57SBarry Smith   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
7478c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
74803049c21SJunchao Zhang       i = 0;
7499566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
7509566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i++]));
75148a46eb9SPierre Jolivet       while ((next = next->next)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i++] : 1.0), next->mat, shell->structure));
7526d7c1e57SBarry Smith     } else {
75303049c21SJunchao Zhang       i = shell->nmat - 1;
7549566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
7559566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i--]));
75648a46eb9SPierre Jolivet       while ((prev = prev->prev)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i--] : 1.0), prev->mat, shell->structure));
7578c8613bfSJakub Kruzik     }
7588c8613bfSJakub Kruzik   } else {
7598c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
7609566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
761b6cfcaf5SJakub Kruzik       while ((next = next->next)) {
7629566063dSJacob Faibussowitsch         PetscCall(MatMatMult(next->mat, tmat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat));
7639566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
7646d7c1e57SBarry Smith         tmat = newmat;
7656d7c1e57SBarry Smith       }
76604d534ceSJakub Kruzik     } else {
7679566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
76804d534ceSJakub Kruzik       while ((prev = prev->prev)) {
7699566063dSJacob Faibussowitsch         PetscCall(MatMatMult(tmat, prev->mat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat));
7709566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
77104d534ceSJakub Kruzik         tmat = newmat;
77204d534ceSJakub Kruzik       }
77304d534ceSJakub Kruzik     }
7749371c9d4SSatish Balay     if (shell->scalings) {
7759371c9d4SSatish Balay       for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
7769371c9d4SSatish Balay     }
7776d7c1e57SBarry Smith   }
7781ba60692SJed Brown 
7799566063dSJacob Faibussowitsch   if ((left = shell->left)) PetscCall(PetscObjectReference((PetscObject)left));
7809566063dSJacob Faibussowitsch   if ((right = shell->right)) PetscCall(PetscObjectReference((PetscObject)right));
7811ba60692SJed Brown 
7829566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(mat, &tmat));
7831ba60692SJed Brown 
7849566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(mat, left, right));
7859566063dSJacob Faibussowitsch   PetscCall(MatScale(mat, scale));
7869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
7879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
7883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
789b52f573bSBarry Smith }
7906a7ede75SJakub Kruzik 
7916a7ede75SJakub Kruzik /*@
792d7f81bf2SJakub Kruzik    MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
7938bd739bdSJakub Kruzik      by summing or computing the product of all the matrices inside the composite matrix.
794d7f81bf2SJakub Kruzik 
795b2b245abSJakub Kruzik   Collective
796d7f81bf2SJakub Kruzik 
797f899ff85SJose E. Roman    Input Parameter:
798d7f81bf2SJakub Kruzik .  mat - the composite matrix
799d7f81bf2SJakub Kruzik 
8004b2558d6SJakub Kruzik    Options Database Keys:
80111a5261eSBarry Smith +  -mat_composite_merge - merge in `MatAssemblyEnd()`
802b28d0daaSJakub Kruzik -  -mat_composite_merge_type - set merge direction
803d7f81bf2SJakub Kruzik 
804d7f81bf2SJakub Kruzik    Level: advanced
805d7f81bf2SJakub Kruzik 
80611a5261eSBarry Smith    Note:
80711a5261eSBarry Smith       The `MatType` of the resulting matrix will be the same as the `MatType` of the FIRST matrix in the composite matrix.
808d7f81bf2SJakub Kruzik 
809*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE`
810d7f81bf2SJakub Kruzik @*/
811d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeMerge(Mat mat)
812d71ae5a4SJacob Faibussowitsch {
813d7f81bf2SJakub Kruzik   PetscFunctionBegin;
814d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
815cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeMerge_C", (Mat), (mat));
8163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
817d7f81bf2SJakub Kruzik }
818d7f81bf2SJakub Kruzik 
819d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat, PetscInt *nmat)
820d71ae5a4SJacob Faibussowitsch {
821d7f81bf2SJakub Kruzik   Mat_Composite *shell = (Mat_Composite *)mat->data;
822d7f81bf2SJakub Kruzik 
823d7f81bf2SJakub Kruzik   PetscFunctionBegin;
824d7f81bf2SJakub Kruzik   *nmat = shell->nmat;
8253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
826d7f81bf2SJakub Kruzik }
827d7f81bf2SJakub Kruzik 
828d7f81bf2SJakub Kruzik /*@
8296d0add67SJakub Kruzik    MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix.
8306a7ede75SJakub Kruzik 
8316a7ede75SJakub Kruzik    Not Collective
8326a7ede75SJakub Kruzik 
8336a7ede75SJakub Kruzik    Input Parameter:
834d7f81bf2SJakub Kruzik .  mat - the composite matrix
8356a7ede75SJakub Kruzik 
8366a7ede75SJakub Kruzik    Output Parameter:
8376d0add67SJakub Kruzik .  nmat - number of matrices in the composite matrix
8386a7ede75SJakub Kruzik 
8398b5c3584SJakub Kruzik    Level: advanced
8406a7ede75SJakub Kruzik 
841*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
8426a7ede75SJakub Kruzik @*/
843d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetNumberMat(Mat mat, PetscInt *nmat)
844d71ae5a4SJacob Faibussowitsch {
8456a7ede75SJakub Kruzik   PetscFunctionBegin;
846d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
847dadcf809SJacob Faibussowitsch   PetscValidIntPointer(nmat, 2);
848cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetNumberMat_C", (Mat, PetscInt *), (mat, nmat));
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
850d7f81bf2SJakub Kruzik }
851d7f81bf2SJakub Kruzik 
852d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMat_Composite(Mat mat, PetscInt i, Mat *Ai)
853d71ae5a4SJacob Faibussowitsch {
854d7f81bf2SJakub Kruzik   Mat_Composite    *shell = (Mat_Composite *)mat->data;
855d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink;
856d7f81bf2SJakub Kruzik   PetscInt          k;
857d7f81bf2SJakub Kruzik 
858d7f81bf2SJakub Kruzik   PetscFunctionBegin;
85908401ef6SPierre Jolivet   PetscCheck(i < shell->nmat, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT, i, shell->nmat);
860d7f81bf2SJakub Kruzik   ilink = shell->head;
861ad540459SPierre Jolivet   for (k = 0; k < i; k++) ilink = ilink->next;
862d7f81bf2SJakub Kruzik   *Ai = ilink->mat;
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8646a7ede75SJakub Kruzik }
8656a7ede75SJakub Kruzik 
8668b5c3584SJakub Kruzik /*@
8678bd739bdSJakub Kruzik    MatCompositeGetMat - Returns the ith matrix from the composite matrix.
8688b5c3584SJakub Kruzik 
869c3339decSBarry Smith    Logically Collective
8708b5c3584SJakub Kruzik 
871d8d19677SJose E. Roman    Input Parameters:
872d7f81bf2SJakub Kruzik +  mat - the composite matrix
8738b5c3584SJakub Kruzik -  i - the number of requested matrix
8748b5c3584SJakub Kruzik 
8758b5c3584SJakub Kruzik    Output Parameter:
8768b5c3584SJakub Kruzik .  Ai - ith matrix in composite
8778b5c3584SJakub Kruzik 
8788b5c3584SJakub Kruzik    Level: advanced
8798b5c3584SJakub Kruzik 
880*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE`
8818b5c3584SJakub Kruzik @*/
882d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMat(Mat mat, PetscInt i, Mat *Ai)
883d71ae5a4SJacob Faibussowitsch {
8848b5c3584SJakub Kruzik   PetscFunctionBegin;
885d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
886d7f81bf2SJakub Kruzik   PetscValidLogicalCollectiveInt(mat, i, 2);
8878b5c3584SJakub Kruzik   PetscValidPointer(Ai, 3);
888cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMat_C", (Mat, PetscInt, Mat *), (mat, i, Ai));
8893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8908b5c3584SJakub Kruzik }
8918b5c3584SJakub Kruzik 
892d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetScalings_Composite(Mat mat, const PetscScalar *scalings)
893d71ae5a4SJacob Faibussowitsch {
89403049c21SJunchao Zhang   Mat_Composite *shell = (Mat_Composite *)mat->data;
89503049c21SJunchao Zhang   PetscInt       nmat;
89603049c21SJunchao Zhang 
89703049c21SJunchao Zhang   PetscFunctionBegin;
8989566063dSJacob Faibussowitsch   PetscCall(MatCompositeGetNumberMat(mat, &nmat));
8999566063dSJacob Faibussowitsch   if (!shell->scalings) PetscCall(PetscMalloc1(nmat, &shell->scalings));
9009566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(shell->scalings, scalings, nmat));
9013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90203049c21SJunchao Zhang }
90303049c21SJunchao Zhang 
90403049c21SJunchao Zhang /*@
90503049c21SJunchao Zhang    MatCompositeSetScalings - Sets separate scaling factors for component matrices.
90603049c21SJunchao Zhang 
907c3339decSBarry Smith    Logically Collective
90803049c21SJunchao Zhang 
909d8d19677SJose E. Roman    Input Parameters:
91003049c21SJunchao Zhang +  mat      - the composite matrix
91103049c21SJunchao Zhang -  scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat)
91203049c21SJunchao Zhang 
91303049c21SJunchao Zhang    Level: advanced
91403049c21SJunchao Zhang 
915*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE`
91603049c21SJunchao Zhang @*/
917d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetScalings(Mat mat, const PetscScalar *scalings)
918d71ae5a4SJacob Faibussowitsch {
91903049c21SJunchao Zhang   PetscFunctionBegin;
92003049c21SJunchao Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
921dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(scalings, 2);
92203049c21SJunchao Zhang   PetscValidLogicalCollectiveScalar(mat, *scalings, 2);
923cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetScalings_C", (Mat, const PetscScalar *), (mat, scalings));
9243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92503049c21SJunchao Zhang }
92603049c21SJunchao Zhang 
927f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL,
928f4259b30SLisandro Dalcin                                        NULL,
929f4259b30SLisandro Dalcin                                        NULL,
93041cd0125SJakub Kruzik                                        MatMult_Composite,
9317bf3a71bSJakub Kruzik                                        MatMultAdd_Composite,
93241cd0125SJakub Kruzik                                        /*  5*/ MatMultTranspose_Composite,
9337bf3a71bSJakub Kruzik                                        MatMultTransposeAdd_Composite,
934f4259b30SLisandro Dalcin                                        NULL,
935f4259b30SLisandro Dalcin                                        NULL,
936f4259b30SLisandro Dalcin                                        NULL,
937f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
938f4259b30SLisandro Dalcin                                        NULL,
939f4259b30SLisandro Dalcin                                        NULL,
940f4259b30SLisandro Dalcin                                        NULL,
941f4259b30SLisandro Dalcin                                        NULL,
942f4259b30SLisandro Dalcin                                        /* 15*/ NULL,
943f4259b30SLisandro Dalcin                                        NULL,
94441cd0125SJakub Kruzik                                        MatGetDiagonal_Composite,
94541cd0125SJakub Kruzik                                        MatDiagonalScale_Composite,
946f4259b30SLisandro Dalcin                                        NULL,
947f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
94841cd0125SJakub Kruzik                                        MatAssemblyEnd_Composite,
949f4259b30SLisandro Dalcin                                        NULL,
950f4259b30SLisandro Dalcin                                        NULL,
951f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
952f4259b30SLisandro Dalcin                                        NULL,
953f4259b30SLisandro Dalcin                                        NULL,
954f4259b30SLisandro Dalcin                                        NULL,
955f4259b30SLisandro Dalcin                                        NULL,
956f4259b30SLisandro Dalcin                                        /* 29*/ NULL,
957f4259b30SLisandro Dalcin                                        NULL,
958f4259b30SLisandro Dalcin                                        NULL,
959f4259b30SLisandro Dalcin                                        NULL,
960f4259b30SLisandro Dalcin                                        NULL,
961f4259b30SLisandro Dalcin                                        /* 34*/ NULL,
962f4259b30SLisandro Dalcin                                        NULL,
963f4259b30SLisandro Dalcin                                        NULL,
964f4259b30SLisandro Dalcin                                        NULL,
965f4259b30SLisandro Dalcin                                        NULL,
966f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
967f4259b30SLisandro Dalcin                                        NULL,
968f4259b30SLisandro Dalcin                                        NULL,
969f4259b30SLisandro Dalcin                                        NULL,
970f4259b30SLisandro Dalcin                                        NULL,
971f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
97241cd0125SJakub Kruzik                                        MatScale_Composite,
97341cd0125SJakub Kruzik                                        MatShift_Basic,
974f4259b30SLisandro Dalcin                                        NULL,
975f4259b30SLisandro Dalcin                                        NULL,
976f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
977f4259b30SLisandro Dalcin                                        NULL,
978f4259b30SLisandro Dalcin                                        NULL,
979f4259b30SLisandro Dalcin                                        NULL,
980f4259b30SLisandro Dalcin                                        NULL,
981f4259b30SLisandro Dalcin                                        /* 54*/ NULL,
982f4259b30SLisandro Dalcin                                        NULL,
983f4259b30SLisandro Dalcin                                        NULL,
984f4259b30SLisandro Dalcin                                        NULL,
985f4259b30SLisandro Dalcin                                        NULL,
986f4259b30SLisandro Dalcin                                        /* 59*/ NULL,
98741cd0125SJakub Kruzik                                        MatDestroy_Composite,
988f4259b30SLisandro Dalcin                                        NULL,
989f4259b30SLisandro Dalcin                                        NULL,
990f4259b30SLisandro Dalcin                                        NULL,
991f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
992f4259b30SLisandro Dalcin                                        NULL,
993f4259b30SLisandro Dalcin                                        NULL,
994f4259b30SLisandro Dalcin                                        NULL,
995f4259b30SLisandro Dalcin                                        NULL,
996f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
997f4259b30SLisandro Dalcin                                        NULL,
998f4259b30SLisandro Dalcin                                        NULL,
999f4259b30SLisandro Dalcin                                        NULL,
1000f4259b30SLisandro Dalcin                                        NULL,
1001f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
1002f4259b30SLisandro Dalcin                                        NULL,
10034b2558d6SJakub Kruzik                                        MatSetFromOptions_Composite,
1004f4259b30SLisandro Dalcin                                        NULL,
1005f4259b30SLisandro Dalcin                                        NULL,
1006f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
1007f4259b30SLisandro Dalcin                                        NULL,
1008f4259b30SLisandro Dalcin                                        NULL,
1009f4259b30SLisandro Dalcin                                        NULL,
1010f4259b30SLisandro Dalcin                                        NULL,
1011f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
1012f4259b30SLisandro Dalcin                                        NULL,
1013f4259b30SLisandro Dalcin                                        NULL,
1014f4259b30SLisandro Dalcin                                        NULL,
1015f4259b30SLisandro Dalcin                                        NULL,
1016f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
1017f4259b30SLisandro Dalcin                                        NULL,
1018f4259b30SLisandro Dalcin                                        NULL,
1019f4259b30SLisandro Dalcin                                        NULL,
1020f4259b30SLisandro Dalcin                                        NULL,
1021f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
1022f4259b30SLisandro Dalcin                                        NULL,
1023f4259b30SLisandro Dalcin                                        NULL,
1024f4259b30SLisandro Dalcin                                        NULL,
1025f4259b30SLisandro Dalcin                                        NULL,
1026f4259b30SLisandro Dalcin                                        /*99*/ NULL,
1027f4259b30SLisandro Dalcin                                        NULL,
1028f4259b30SLisandro Dalcin                                        NULL,
1029f4259b30SLisandro Dalcin                                        NULL,
1030f4259b30SLisandro Dalcin                                        NULL,
1031f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1032f4259b30SLisandro Dalcin                                        NULL,
1033f4259b30SLisandro Dalcin                                        NULL,
1034f4259b30SLisandro Dalcin                                        NULL,
1035f4259b30SLisandro Dalcin                                        NULL,
1036f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1037f4259b30SLisandro Dalcin                                        NULL,
1038f4259b30SLisandro Dalcin                                        NULL,
1039f4259b30SLisandro Dalcin                                        NULL,
1040f4259b30SLisandro Dalcin                                        NULL,
1041f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1042f4259b30SLisandro Dalcin                                        NULL,
1043f4259b30SLisandro Dalcin                                        NULL,
1044f4259b30SLisandro Dalcin                                        NULL,
1045f4259b30SLisandro Dalcin                                        NULL,
1046f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1047f4259b30SLisandro Dalcin                                        NULL,
1048f4259b30SLisandro Dalcin                                        NULL,
1049f4259b30SLisandro Dalcin                                        NULL,
1050f4259b30SLisandro Dalcin                                        NULL,
1051f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1052f4259b30SLisandro Dalcin                                        NULL,
1053f4259b30SLisandro Dalcin                                        NULL,
1054f4259b30SLisandro Dalcin                                        NULL,
1055f4259b30SLisandro Dalcin                                        NULL,
1056f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1057f4259b30SLisandro Dalcin                                        NULL,
1058f4259b30SLisandro Dalcin                                        NULL,
1059f4259b30SLisandro Dalcin                                        NULL,
1060f4259b30SLisandro Dalcin                                        NULL,
1061f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1062f4259b30SLisandro Dalcin                                        NULL,
1063f4259b30SLisandro Dalcin                                        NULL,
1064f4259b30SLisandro Dalcin                                        NULL,
1065f4259b30SLisandro Dalcin                                        NULL,
1066f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1067f4259b30SLisandro Dalcin                                        NULL,
1068d70f29a3SPierre Jolivet                                        NULL,
1069d70f29a3SPierre Jolivet                                        NULL,
1070d70f29a3SPierre Jolivet                                        NULL,
1071d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1072d70f29a3SPierre Jolivet                                        NULL,
1073d70f29a3SPierre Jolivet                                        NULL,
107499a7f59eSMark Adams                                        NULL,
107599a7f59eSMark Adams                                        NULL,
10767fb60732SBarry Smith                                        NULL,
1077dec0b466SHong Zhang                                        /*150*/ NULL,
1078dec0b466SHong Zhang                                        NULL};
107941cd0125SJakub Kruzik 
108041cd0125SJakub Kruzik /*MC
108141cd0125SJakub Kruzik    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices.
108241cd0125SJakub Kruzik     The matrices need to have a correct size and parallel layout for the sum or product to be valid.
108341cd0125SJakub Kruzik 
1084*2ef1f0ffSBarry Smith   Level: advanced
1085*2ef1f0ffSBarry Smith 
108611a5261eSBarry Smith    Note:
108711a5261eSBarry Smith    To use the product of the matrices call `MatCompositeSetType`(mat,`MAT_COMPOSITE_MULTIPLICATIVE`);
108841cd0125SJakub Kruzik 
1089*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`,
109011a5261eSBarry Smith           `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()`
109141cd0125SJakub Kruzik M*/
109241cd0125SJakub Kruzik 
1093d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A)
1094d71ae5a4SJacob Faibussowitsch {
109541cd0125SJakub Kruzik   Mat_Composite *b;
109641cd0125SJakub Kruzik 
109741cd0125SJakub Kruzik   PetscFunctionBegin;
10984dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
109941cd0125SJakub Kruzik   A->data = (void *)b;
11009566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A->ops, &MatOps_Values, sizeof(struct _MatOps)));
110141cd0125SJakub Kruzik 
11029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
11039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
110441cd0125SJakub Kruzik 
110541cd0125SJakub Kruzik   A->assembled    = PETSC_TRUE;
110641cd0125SJakub Kruzik   A->preallocated = PETSC_TRUE;
110741cd0125SJakub Kruzik   b->type         = MAT_COMPOSITE_ADDITIVE;
110841cd0125SJakub Kruzik   b->scale        = 1.0;
110941cd0125SJakub Kruzik   b->nmat         = 0;
11104b2558d6SJakub Kruzik   b->merge        = PETSC_FALSE;
11118c8613bfSJakub Kruzik   b->mergetype    = MAT_COMPOSITE_MERGE_RIGHT;
11123b35acafSJakub Kruzik   b->structure    = DIFFERENT_NONZERO_PATTERN;
111303049c21SJunchao Zhang   b->merge_mvctx  = PETSC_TRUE;
111403049c21SJunchao Zhang 
11159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeAddMat_C", MatCompositeAddMat_Composite));
11169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetType_C", MatCompositeSetType_Composite));
11179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetType_C", MatCompositeGetType_Composite));
11189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMergeType_C", MatCompositeSetMergeType_Composite));
11199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMatStructure_C", MatCompositeSetMatStructure_Composite));
11209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMatStructure_C", MatCompositeGetMatStructure_Composite));
11219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeMerge_C", MatCompositeMerge_Composite));
11229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetNumberMat_C", MatCompositeGetNumberMat_Composite));
11239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMat_C", MatCompositeGetMat_Composite));
11249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetScalings_C", MatCompositeSetScalings_Composite));
112541cd0125SJakub Kruzik 
11269566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCOMPOSITE));
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112841cd0125SJakub Kruzik }
1129