xref: /petsc/src/mat/impls/composite/mcomposite.c (revision b9c875b8ccb0466cfbdbcc08c79b487b3a7d395f)
1326b7573SPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2793850ffSBarry Smith 
3f4259b30SLisandro Dalcin const char *const MatCompositeMergeTypes[] = {"left", "right", "MatCompositeMergeType", "MAT_COMPOSITE_", NULL};
48c8613bfSJakub Kruzik 
5793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink;
6793850ffSBarry Smith struct _Mat_CompositeLink {
7793850ffSBarry Smith   Mat               mat;
86d7c1e57SBarry Smith   Vec               work;
96d7c1e57SBarry Smith   Mat_CompositeLink next, prev;
10793850ffSBarry Smith };
11793850ffSBarry Smith 
12793850ffSBarry Smith typedef struct {
136d7c1e57SBarry Smith   MatCompositeType      type;
146d7c1e57SBarry Smith   Mat_CompositeLink     head, tail;
15793850ffSBarry Smith   Vec                   work;
166a7ede75SJakub Kruzik   PetscInt              nmat;
174b2558d6SJakub Kruzik   PetscBool             merge;
188c8613bfSJakub Kruzik   MatCompositeMergeType mergetype;
193b35acafSJakub Kruzik   MatStructure          structure;
2003049c21SJunchao Zhang 
2103049c21SJunchao Zhang   PetscScalar *scalings;
2203049c21SJunchao Zhang   PetscBool    merge_mvctx; /* Whether need to merge mvctx of component matrices */
2303049c21SJunchao Zhang   Vec         *lvecs;       /* [nmat] Basically, they are Mvctx->lvec of each component matrix */
2403049c21SJunchao Zhang   PetscScalar *larray;      /* [len] Data arrays of lvecs[] are stored consecutively in larray */
2503049c21SJunchao Zhang   PetscInt     len;         /* Length of larray[] */
2603049c21SJunchao Zhang   Vec          gvec;        /* Union of lvecs[] without duplicated entries */
2703049c21SJunchao Zhang   PetscInt    *location;    /* A map that maps entries in garray[] to larray[] */
2803049c21SJunchao Zhang   VecScatter   Mvctx;
29793850ffSBarry Smith } Mat_Composite;
30793850ffSBarry Smith 
3166976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_Composite(Mat mat)
32d71ae5a4SJacob Faibussowitsch {
33326b7573SPierre Jolivet   Mat_Composite    *shell;
34326b7573SPierre Jolivet   Mat_CompositeLink next, oldnext;
3503049c21SJunchao Zhang   PetscInt          i;
36793850ffSBarry Smith 
37793850ffSBarry Smith   PetscFunctionBegin;
38326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
39326b7573SPierre Jolivet   next = shell->head;
40793850ffSBarry Smith   while (next) {
419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&next->mat));
4248a46eb9SPierre Jolivet     if (next->work && (!next->next || next->work != next->next->work)) PetscCall(VecDestroy(&next->work));
436d7c1e57SBarry Smith     oldnext = next;
44793850ffSBarry Smith     next    = next->next;
459566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldnext));
46793850ffSBarry Smith   }
479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->work));
4803049c21SJunchao Zhang 
4903049c21SJunchao Zhang   if (shell->Mvctx) {
509566063dSJacob Faibussowitsch     for (i = 0; i < shell->nmat; i++) PetscCall(VecDestroy(&shell->lvecs[i]));
519566063dSJacob Faibussowitsch     PetscCall(PetscFree3(shell->location, shell->larray, shell->lvecs));
529566063dSJacob Faibussowitsch     PetscCall(PetscFree(shell->larray));
539566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&shell->gvec));
549566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&shell->Mvctx));
5503049c21SJunchao Zhang   }
5603049c21SJunchao Zhang 
579566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell->scalings));
582e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeAddMat_C", NULL));
592e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetType_C", NULL));
602e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetType_C", NULL));
612e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMergeType_C", NULL));
622e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMatStructure_C", NULL));
632e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMatStructure_C", NULL));
642e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeMerge_C", NULL));
652e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetNumberMat_C", NULL));
662e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMat_C", NULL));
672e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetScalings_C", NULL));
68326b7573SPierre Jolivet   PetscCall(PetscFree(shell));
69326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL)); // needed to avoid a call to MatShellSetContext_Immutable()
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71793850ffSBarry Smith }
72793850ffSBarry Smith 
7366976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_Composite_Multiplicative(Mat A, Vec x, Vec y)
74d71ae5a4SJacob Faibussowitsch {
75326b7573SPierre Jolivet   Mat_Composite    *shell;
76326b7573SPierre Jolivet   Mat_CompositeLink next;
77326b7573SPierre Jolivet   Vec               out;
786d7c1e57SBarry Smith 
796d7c1e57SBarry Smith   PetscFunctionBegin;
80326b7573SPierre Jolivet   PetscCall(MatShellGetContext(A, &shell));
81326b7573SPierre Jolivet   next = shell->head;
8228b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
836d7c1e57SBarry Smith   while (next->next) {
846d7c1e57SBarry Smith     if (!next->work) { /* should reuse previous work if the same size */
859566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(next->mat, NULL, &next->work));
866d7c1e57SBarry Smith     }
876d7c1e57SBarry Smith     out = next->work;
88326b7573SPierre Jolivet     PetscCall(MatMult(next->mat, x, out));
89326b7573SPierre Jolivet     x    = out;
906d7c1e57SBarry Smith     next = next->next;
916d7c1e57SBarry Smith   }
92326b7573SPierre Jolivet   PetscCall(MatMult(next->mat, x, y));
939371c9d4SSatish Balay   if (shell->scalings) {
94326b7573SPierre Jolivet     PetscScalar scale = 1.0;
95326b7573SPierre Jolivet     for (PetscInt i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
969566063dSJacob Faibussowitsch     PetscCall(VecScale(y, scale));
97326b7573SPierre Jolivet   }
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
996d7c1e57SBarry Smith }
1006d7c1e57SBarry Smith 
10166976f2fSJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A, Vec x, Vec y)
102d71ae5a4SJacob Faibussowitsch {
103326b7573SPierre Jolivet   Mat_Composite    *shell;
104326b7573SPierre Jolivet   Mat_CompositeLink tail;
105326b7573SPierre Jolivet   Vec               out;
1066d7c1e57SBarry Smith 
1076d7c1e57SBarry Smith   PetscFunctionBegin;
108326b7573SPierre Jolivet   PetscCall(MatShellGetContext(A, &shell));
109326b7573SPierre Jolivet   tail = shell->tail;
11028b400f6SJacob Faibussowitsch   PetscCheck(tail, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
1116d7c1e57SBarry Smith   while (tail->prev) {
1126d7c1e57SBarry Smith     if (!tail->prev->work) { /* should reuse previous work if the same size */
1139566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(tail->mat, NULL, &tail->prev->work));
1146d7c1e57SBarry Smith     }
1156d7c1e57SBarry Smith     out = tail->prev->work;
116326b7573SPierre Jolivet     PetscCall(MatMultTranspose(tail->mat, x, out));
117326b7573SPierre Jolivet     x    = out;
1186d7c1e57SBarry Smith     tail = tail->prev;
1196d7c1e57SBarry Smith   }
120326b7573SPierre Jolivet   PetscCall(MatMultTranspose(tail->mat, x, y));
1219371c9d4SSatish Balay   if (shell->scalings) {
122326b7573SPierre Jolivet     PetscScalar scale = 1.0;
123326b7573SPierre Jolivet     for (PetscInt i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
1249566063dSJacob Faibussowitsch     PetscCall(VecScale(y, scale));
125326b7573SPierre Jolivet   }
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1276d7c1e57SBarry Smith }
1286d7c1e57SBarry Smith 
12966976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_Composite(Mat mat, Vec x, Vec y)
130d71ae5a4SJacob Faibussowitsch {
131326b7573SPierre Jolivet   Mat_Composite     *shell;
132326b7573SPierre Jolivet   Mat_CompositeLink  cur;
133326b7573SPierre Jolivet   Vec                y2, xin;
13403049c21SJunchao Zhang   Mat                A, B;
13503049c21SJunchao Zhang   PetscInt           i, j, k, n, nuniq, lo, hi, mid, *gindices, *buf, *tmp, tot;
13603049c21SJunchao Zhang   const PetscScalar *vals;
13703049c21SJunchao Zhang   const PetscInt    *garray;
13803049c21SJunchao Zhang   IS                 ix, iy;
13903049c21SJunchao Zhang   PetscBool          match;
140793850ffSBarry Smith 
141793850ffSBarry Smith   PetscFunctionBegin;
142326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
143326b7573SPierre Jolivet   cur = shell->head;
14428b400f6SJacob Faibussowitsch   PetscCheck(cur, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
14503049c21SJunchao Zhang 
14603049c21SJunchao Zhang   /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time
14703049c21SJunchao Zhang      we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and
14803049c21SJunchao Zhang      it is legal to merge Mvctx, because all component matrices have the same size.
14903049c21SJunchao Zhang    */
15003049c21SJunchao Zhang   if (shell->merge_mvctx && !shell->Mvctx) {
15103049c21SJunchao Zhang     /* Currently only implemented for MATMPIAIJ */
15203049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
1539566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat, MATMPIAIJ, &match));
15403049c21SJunchao Zhang       if (!match) {
15503049c21SJunchao Zhang         shell->merge_mvctx = PETSC_FALSE;
15603049c21SJunchao Zhang         goto skip_merge_mvctx;
157e4fc5df0SBarry Smith       }
158e4fc5df0SBarry Smith     }
15903049c21SJunchao Zhang 
16003049c21SJunchao Zhang     /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */
16103049c21SJunchao Zhang     tot = 0;
16203049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
1639566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, NULL));
1649566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
16503049c21SJunchao Zhang       tot += n;
16603049c21SJunchao Zhang     }
1679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(tot, &shell->location, tot, &shell->larray, shell->nmat, &shell->lvecs));
16803049c21SJunchao Zhang     shell->len = tot;
16903049c21SJunchao Zhang 
17003049c21SJunchao Zhang     /* Go through matrices second time to sort off-diag columns and remove dups */
1719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot, &gindices)); /* No Malloc2() since we will give one to petsc and free the other */
1729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot, &buf));
17303049c21SJunchao Zhang     nuniq = 0; /* Number of unique nonzero columns */
17403049c21SJunchao Zhang     for (cur = shell->head; cur; cur = cur->next) {
1759566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray));
1769566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
17703049c21SJunchao Zhang       /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */
17803049c21SJunchao Zhang       i = j = k = 0;
17903049c21SJunchao Zhang       while (i < n && j < nuniq) {
18003049c21SJunchao Zhang         if (garray[i] < gindices[j]) buf[k++] = garray[i++];
18103049c21SJunchao Zhang         else if (garray[i] > gindices[j]) buf[k++] = gindices[j++];
1829371c9d4SSatish Balay         else {
1839371c9d4SSatish Balay           buf[k++] = garray[i++];
1849371c9d4SSatish Balay           j++;
1859371c9d4SSatish Balay         }
18603049c21SJunchao Zhang       }
18703049c21SJunchao Zhang       /* Copy leftover in garray[] or gindices[] */
18803049c21SJunchao Zhang       if (i < n) {
1899566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf + k, garray + i, n - i));
19003049c21SJunchao Zhang         nuniq = k + n - i;
19103049c21SJunchao Zhang       } else if (j < nuniq) {
1929566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf + k, gindices + j, nuniq - j));
19303049c21SJunchao Zhang         nuniq = k + nuniq - j;
19403049c21SJunchao Zhang       } else nuniq = k;
19503049c21SJunchao Zhang       /* Swap gindices and buf to merge garray of the next matrix */
19603049c21SJunchao Zhang       tmp      = gindices;
19703049c21SJunchao Zhang       gindices = buf;
19803049c21SJunchao Zhang       buf      = tmp;
19903049c21SJunchao Zhang     }
2009566063dSJacob Faibussowitsch     PetscCall(PetscFree(buf));
20103049c21SJunchao Zhang 
20203049c21SJunchao Zhang     /* Go through matrices third time to build a map from gindices[] to garray[] */
20303049c21SJunchao Zhang     tot = 0;
20403049c21SJunchao Zhang     for (cur = shell->head, j = 0; cur; cur = cur->next, j++) { /* j-th matrix */
2059566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray));
2069566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
2079566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, n, NULL, &shell->lvecs[j]));
20803049c21SJunchao Zhang       /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */
20903049c21SJunchao Zhang       lo = 0;
21003049c21SJunchao Zhang       for (i = 0; i < n; i++) {
21103049c21SJunchao Zhang         hi = nuniq;
21203049c21SJunchao Zhang         while (hi - lo > 1) {
21303049c21SJunchao Zhang           mid = lo + (hi - lo) / 2;
21403049c21SJunchao Zhang           if (garray[i] < gindices[mid]) hi = mid;
21503049c21SJunchao Zhang           else lo = mid;
21603049c21SJunchao Zhang         }
21703049c21SJunchao Zhang         shell->location[tot + i] = lo; /* gindices[lo] = garray[i] */
21803049c21SJunchao Zhang         lo++;                          /* Since garray[i+1] > garray[i], we can safely advance lo */
21903049c21SJunchao Zhang       }
22003049c21SJunchao Zhang       tot += n;
22103049c21SJunchao Zhang     }
22203049c21SJunchao Zhang 
22303049c21SJunchao Zhang     /* Build merged Mvctx */
2249566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nuniq, gindices, PETSC_OWN_POINTER, &ix));
2259566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, nuniq, 0, 1, &iy));
2269566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &xin));
2279566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, nuniq, &shell->gvec));
2289566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(xin, ix, shell->gvec, iy, &shell->Mvctx));
2299566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xin));
2309566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ix));
2319566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iy));
23203049c21SJunchao Zhang   }
23303049c21SJunchao Zhang 
23403049c21SJunchao Zhang skip_merge_mvctx:
2359566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0));
236326b7573SPierre Jolivet   if (!((Mat_Shell *)mat->data)->left_work) PetscCall(VecDuplicate(y, &(((Mat_Shell *)mat->data)->left_work)));
237326b7573SPierre Jolivet   y2 = ((Mat_Shell *)mat->data)->left_work;
23803049c21SJunchao Zhang 
23903049c21SJunchao Zhang   if (shell->Mvctx) { /* Have a merged Mvctx */
24003049c21SJunchao 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
241da81f932SPierre 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
24203049c21SJunchao Zhang        overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach.
24303049c21SJunchao Zhang      */
244326b7573SPierre Jolivet     PetscCall(VecScatterBegin(shell->Mvctx, x, shell->gvec, INSERT_VALUES, SCATTER_FORWARD));
245326b7573SPierre Jolivet     PetscCall(VecScatterEnd(shell->Mvctx, x, shell->gvec, INSERT_VALUES, SCATTER_FORWARD));
24603049c21SJunchao Zhang 
2479566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(shell->gvec, &vals));
24803049c21SJunchao Zhang     for (i = 0; i < shell->len; i++) shell->larray[i] = vals[shell->location[i]];
2499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(shell->gvec, &vals));
25003049c21SJunchao Zhang 
25103049c21SJunchao Zhang     for (cur = shell->head, tot = i = 0; cur; cur = cur->next, i++) { /* i-th matrix */
2529566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, &A, &B, NULL));
253326b7573SPierre Jolivet       PetscUseTypeMethod(A, mult, x, y2);
2549566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B, NULL, &n));
2559566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(shell->lvecs[i], &shell->larray[tot]));
2569927e4dfSBarry Smith       PetscUseTypeMethod(B, multadd, shell->lvecs[i], y2, y2);
2579566063dSJacob Faibussowitsch       PetscCall(VecResetArray(shell->lvecs[i]));
2589566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, (shell->scalings ? shell->scalings[i] : 1.0), y2));
25903049c21SJunchao Zhang       tot += n;
26003049c21SJunchao Zhang     }
26103049c21SJunchao Zhang   } else {
26203049c21SJunchao Zhang     if (shell->scalings) {
26303049c21SJunchao Zhang       for (cur = shell->head, i = 0; cur; cur = cur->next, i++) {
264326b7573SPierre Jolivet         PetscCall(MatMult(cur->mat, x, y2));
2659566063dSJacob Faibussowitsch         PetscCall(VecAXPY(y, shell->scalings[i], y2));
26603049c21SJunchao Zhang       }
26703049c21SJunchao Zhang     } else {
268326b7573SPierre Jolivet       for (cur = shell->head; cur; cur = cur->next) PetscCall(MatMultAdd(cur->mat, x, y, y));
26903049c21SJunchao Zhang     }
27003049c21SJunchao Zhang   }
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272793850ffSBarry Smith }
273793850ffSBarry Smith 
27466976f2fSJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Composite(Mat A, Vec x, Vec y)
275d71ae5a4SJacob Faibussowitsch {
276326b7573SPierre Jolivet   Mat_Composite    *shell;
277326b7573SPierre Jolivet   Mat_CompositeLink next;
278326b7573SPierre Jolivet   Vec               y2 = NULL;
27903049c21SJunchao Zhang   PetscInt          i;
280793850ffSBarry Smith 
281793850ffSBarry Smith   PetscFunctionBegin;
282326b7573SPierre Jolivet   PetscCall(MatShellGetContext(A, &shell));
283326b7573SPierre Jolivet   next = shell->head;
28428b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
28503049c21SJunchao Zhang 
286326b7573SPierre Jolivet   PetscCall(MatMultTranspose(next->mat, x, y));
28703049c21SJunchao Zhang   if (shell->scalings) {
2889566063dSJacob Faibussowitsch     PetscCall(VecScale(y, shell->scalings[0]));
289326b7573SPierre Jolivet     if (!((Mat_Shell *)A->data)->right_work) PetscCall(VecDuplicate(y, &(((Mat_Shell *)A->data)->right_work)));
290326b7573SPierre Jolivet     y2 = ((Mat_Shell *)A->data)->right_work;
29103049c21SJunchao Zhang   }
29203049c21SJunchao Zhang   i = 1;
293e4fc5df0SBarry Smith   while ((next = next->next)) {
294326b7573SPierre Jolivet     if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat, x, y, y));
29503049c21SJunchao Zhang     else {
296326b7573SPierre Jolivet       PetscCall(MatMultTranspose(next->mat, x, y2));
2979566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, shell->scalings[i++], y2));
29803049c21SJunchao Zhang     }
299e4fc5df0SBarry Smith   }
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3017bf3a71bSJakub Kruzik }
3027bf3a71bSJakub Kruzik 
30366976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Composite(Mat A, Vec v)
304d71ae5a4SJacob Faibussowitsch {
305326b7573SPierre Jolivet   Mat_Composite    *shell;
306326b7573SPierre Jolivet   Mat_CompositeLink next;
30703049c21SJunchao Zhang   PetscInt          i;
308793850ffSBarry Smith 
309793850ffSBarry Smith   PetscFunctionBegin;
310326b7573SPierre Jolivet   PetscCall(MatShellGetContext(A, &shell));
311326b7573SPierre Jolivet   next = shell->head;
31228b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
3139566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(next->mat, v));
3149566063dSJacob Faibussowitsch   if (shell->scalings) PetscCall(VecScale(v, shell->scalings[0]));
31503049c21SJunchao Zhang 
31648a46eb9SPierre Jolivet   if (next->next && !shell->work) PetscCall(VecDuplicate(v, &shell->work));
31703049c21SJunchao Zhang   i = 1;
318793850ffSBarry Smith   while ((next = next->next)) {
3199566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(next->mat, shell->work));
3209566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v, (shell->scalings ? shell->scalings[i++] : 1.0), shell->work));
321793850ffSBarry Smith   }
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323793850ffSBarry Smith }
324793850ffSBarry Smith 
32566976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_Composite(Mat Y, MatAssemblyType t)
326d71ae5a4SJacob Faibussowitsch {
327326b7573SPierre Jolivet   Mat_Composite *shell;
328b52f573bSBarry Smith 
329793850ffSBarry Smith   PetscFunctionBegin;
330326b7573SPierre Jolivet   PetscCall(MatShellGetContext(Y, &shell));
3311baa6e33SBarry Smith   if (shell->merge) PetscCall(MatCompositeMerge(Y));
332326b7573SPierre Jolivet   else PetscCall(MatAssemblyEnd_Shell(Y, t));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
334e4fc5df0SBarry Smith }
335793850ffSBarry Smith 
33666976f2fSJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Composite(Mat A, PetscOptionItems *PetscOptionsObject)
337d71ae5a4SJacob Faibussowitsch {
338326b7573SPierre Jolivet   Mat_Composite *a;
3394b2558d6SJakub Kruzik 
3404b2558d6SJakub Kruzik   PetscFunctionBegin;
341326b7573SPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
342d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MATCOMPOSITE options");
3439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge", "Merge at MatAssemblyEnd", "MatCompositeMerge", a->merge, &a->merge, NULL));
3449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-mat_composite_merge_type", "Set composite merge direction", "MatCompositeSetMergeType", MatCompositeMergeTypes, (PetscEnum)a->mergetype, (PetscEnum *)&a->mergetype, NULL));
3459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx", "Merge MatMult() vecscat contexts", "MatCreateComposite", a->merge_mvctx, &a->merge_mvctx, NULL));
346d0609cedSBarry Smith   PetscOptionsHeadEnd();
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3484b2558d6SJakub Kruzik }
3494b2558d6SJakub Kruzik 
3502c0821f3SBarry Smith /*@
3518bd739bdSJakub Kruzik   MatCreateComposite - Creates a matrix as the sum or product of one or more matrices
352793850ffSBarry Smith 
353d083f849SBarry Smith   Collective
354793850ffSBarry Smith 
355793850ffSBarry Smith   Input Parameters:
356793850ffSBarry Smith + comm - MPI communicator
357793850ffSBarry Smith . nmat - number of matrices to put in
358793850ffSBarry Smith - mats - the matrices
359793850ffSBarry Smith 
360793850ffSBarry Smith   Output Parameter:
361db36af27SMatthew G. Knepley . mat - the matrix
362793850ffSBarry Smith 
3634b2558d6SJakub Kruzik   Options Database Keys:
36411a5261eSBarry Smith + -mat_composite_merge       - merge in `MatAssemblyEnd()`
36511a5261eSBarry Smith . -mat_composite_merge_mvctx - merge Mvctx of component matrices to optimize communication in `MatMult()` for ADDITIVE matrices
366b28d0daaSJakub Kruzik - -mat_composite_merge_type  - set merge direction
3674b2558d6SJakub Kruzik 
368793850ffSBarry Smith   Level: advanced
369793850ffSBarry Smith 
37011a5261eSBarry Smith   Note:
371793850ffSBarry Smith   Alternative construction
3722ef1f0ffSBarry Smith .vb
3732ef1f0ffSBarry Smith        MatCreate(comm,&mat);
3742ef1f0ffSBarry Smith        MatSetSizes(mat,m,n,M,N);
3752ef1f0ffSBarry Smith        MatSetType(mat,MATCOMPOSITE);
3762ef1f0ffSBarry Smith        MatCompositeAddMat(mat,mats[0]);
3772ef1f0ffSBarry Smith        ....
3782ef1f0ffSBarry Smith        MatCompositeAddMat(mat,mats[nmat-1]);
3792ef1f0ffSBarry Smith        MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
3802ef1f0ffSBarry Smith        MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
3812ef1f0ffSBarry Smith .ve
382793850ffSBarry Smith 
3836d7c1e57SBarry Smith   For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0]
3846d7c1e57SBarry Smith 
3851cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCompositeGetMat()`, `MatCompositeMerge()`, `MatCompositeSetType()`,
38620f4b53cSBarry Smith           `MATCOMPOSITE`, `MatCompositeType`
387793850ffSBarry Smith @*/
388d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateComposite(MPI_Comm comm, PetscInt nmat, const Mat *mats, Mat *mat)
389d71ae5a4SJacob Faibussowitsch {
390793850ffSBarry Smith   PetscInt m, n, M, N, i;
391793850ffSBarry Smith 
392793850ffSBarry Smith   PetscFunctionBegin;
39308401ef6SPierre Jolivet   PetscCheck(nmat >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in at least one matrix");
3944f572ea9SToby Isaac   PetscAssertPointer(mat, 4);
395793850ffSBarry Smith 
3969566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[0], PETSC_IGNORE, &n));
3979566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[nmat - 1], &m, PETSC_IGNORE));
3989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[0], PETSC_IGNORE, &N));
3999566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[nmat - 1], &M, PETSC_IGNORE));
4009566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, mat));
4019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mat, m, n, M, N));
4029566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat, MATCOMPOSITE));
40348a46eb9SPierre Jolivet   for (i = 0; i < nmat; i++) PetscCall(MatCompositeAddMat(*mat, mats[i]));
4049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
4059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
407793850ffSBarry Smith }
408793850ffSBarry Smith 
409d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeAddMat_Composite(Mat mat, Mat smat)
410d71ae5a4SJacob Faibussowitsch {
411326b7573SPierre Jolivet   Mat_Composite    *shell;
412326b7573SPierre Jolivet   Mat_CompositeLink ilink, next;
4135c8dc79eSRichard Tran Mills   VecType           vtype_mat, vtype_smat;
4145c8dc79eSRichard Tran Mills   PetscBool         match;
415d7f81bf2SJakub Kruzik 
416d7f81bf2SJakub Kruzik   PetscFunctionBegin;
417326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
418326b7573SPierre Jolivet   next = shell->head;
4194dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
420f4259b30SLisandro Dalcin   ilink->next = NULL;
4219566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)smat));
422d7f81bf2SJakub Kruzik   ilink->mat = smat;
423d7f81bf2SJakub Kruzik 
424d7f81bf2SJakub Kruzik   if (!next) shell->head = ilink;
425d7f81bf2SJakub Kruzik   else {
426ad540459SPierre Jolivet     while (next->next) next = next->next;
427d7f81bf2SJakub Kruzik     next->next  = ilink;
428d7f81bf2SJakub Kruzik     ilink->prev = next;
429d7f81bf2SJakub Kruzik   }
430d7f81bf2SJakub Kruzik   shell->tail = ilink;
431d7f81bf2SJakub Kruzik   shell->nmat += 1;
43203049c21SJunchao Zhang 
4335c8dc79eSRichard Tran Mills   /* If all of the partial matrices have the same default vector type, then the composite matrix should also have this default type.
4345c8dc79eSRichard Tran Mills      Otherwise, the default type should be "standard". */
4355c8dc79eSRichard Tran Mills   PetscCall(MatGetVecType(smat, &vtype_smat));
4365c8dc79eSRichard Tran Mills   if (shell->nmat == 1) PetscCall(MatSetVecType(mat, vtype_smat));
4375c8dc79eSRichard Tran Mills   else {
4385c8dc79eSRichard Tran Mills     PetscCall(MatGetVecType(mat, &vtype_mat));
4395c8dc79eSRichard Tran Mills     PetscCall(PetscStrcmp(vtype_smat, vtype_mat, &match));
4405c8dc79eSRichard Tran Mills     if (!match) PetscCall(MatSetVecType(mat, VECSTANDARD));
4415c8dc79eSRichard Tran Mills   }
4425c8dc79eSRichard Tran Mills 
44303049c21SJunchao Zhang   /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */
44403049c21SJunchao Zhang   if (shell->scalings) {
4459566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(sizeof(PetscScalar) * shell->nmat, &shell->scalings));
44603049c21SJunchao Zhang     shell->scalings[shell->nmat - 1] = 1.0;
44703049c21SJunchao Zhang   }
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
449d7f81bf2SJakub Kruzik }
450d7f81bf2SJakub Kruzik 
451793850ffSBarry Smith /*@
4528bd739bdSJakub Kruzik   MatCompositeAddMat - Add another matrix to a composite matrix.
453793850ffSBarry Smith 
454c3339decSBarry Smith   Collective
455793850ffSBarry Smith 
456793850ffSBarry Smith   Input Parameters:
457793850ffSBarry Smith + mat  - the composite matrix
458793850ffSBarry Smith - smat - the partial matrix
459793850ffSBarry Smith 
460793850ffSBarry Smith   Level: advanced
461793850ffSBarry Smith 
4621cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
463793850ffSBarry Smith @*/
464d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeAddMat(Mat mat, Mat smat)
465d71ae5a4SJacob Faibussowitsch {
466793850ffSBarry Smith   PetscFunctionBegin;
4670700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4680700a824SBarry Smith   PetscValidHeaderSpecific(smat, MAT_CLASSID, 2);
469cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeAddMat_C", (Mat, Mat), (mat, smat));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471d7f81bf2SJakub Kruzik }
472793850ffSBarry Smith 
473d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetType_Composite(Mat mat, MatCompositeType type)
474d71ae5a4SJacob Faibussowitsch {
475326b7573SPierre Jolivet   Mat_Composite *b;
476d7f81bf2SJakub Kruzik 
477d7f81bf2SJakub Kruzik   PetscFunctionBegin;
478326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &b));
479ced1ac25SJakub Kruzik   b->type = type;
480d7f81bf2SJakub Kruzik   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
481326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, NULL));
482326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Composite_Multiplicative));
483326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Composite_Multiplicative));
48403049c21SJunchao Zhang     b->merge_mvctx = PETSC_FALSE;
485d7f81bf2SJakub Kruzik   } else {
486326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Composite));
487326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (void (*)(void))MatMult_Composite));
488326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Composite));
489793850ffSBarry Smith   }
4903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4916d7c1e57SBarry Smith }
4926d7c1e57SBarry Smith 
4932c0821f3SBarry Smith /*@
4948bd739bdSJakub Kruzik   MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product.
4956d7c1e57SBarry Smith 
496c3339decSBarry Smith   Logically Collective
4976d7c1e57SBarry Smith 
4986d7c1e57SBarry Smith   Input Parameters:
49920f4b53cSBarry Smith + mat  - the composite matrix
50020f4b53cSBarry Smith - type - the `MatCompositeType` to use for the matrix
5016d7c1e57SBarry Smith 
5026d7c1e57SBarry Smith   Level: advanced
5036d7c1e57SBarry Smith 
5041cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE`,
50520f4b53cSBarry Smith           `MatCompositeType`
5066d7c1e57SBarry Smith @*/
507d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetType(Mat mat, MatCompositeType type)
508d71ae5a4SJacob Faibussowitsch {
5096d7c1e57SBarry Smith   PetscFunctionBegin;
510d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
511b2b245abSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
512cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetType_C", (Mat, MatCompositeType), (mat, type));
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
514793850ffSBarry Smith }
515793850ffSBarry Smith 
516d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetType_Composite(Mat mat, MatCompositeType *type)
517d71ae5a4SJacob Faibussowitsch {
518326b7573SPierre Jolivet   Mat_Composite *shell;
5196dbc55e5SJakub Kruzik 
5206dbc55e5SJakub Kruzik   PetscFunctionBegin;
521326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
522326b7573SPierre Jolivet   *type = shell->type;
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5246dbc55e5SJakub Kruzik }
5256dbc55e5SJakub Kruzik 
5266dbc55e5SJakub Kruzik /*@
5276dbc55e5SJakub Kruzik   MatCompositeGetType - Returns type of composite.
5286dbc55e5SJakub Kruzik 
5296dbc55e5SJakub Kruzik   Not Collective
5306dbc55e5SJakub Kruzik 
5316dbc55e5SJakub Kruzik   Input Parameter:
5326dbc55e5SJakub Kruzik . mat - the composite matrix
5336dbc55e5SJakub Kruzik 
5346dbc55e5SJakub Kruzik   Output Parameter:
5356dbc55e5SJakub Kruzik . type - type of composite
5366dbc55e5SJakub Kruzik 
5376dbc55e5SJakub Kruzik   Level: advanced
5386dbc55e5SJakub Kruzik 
5391cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE`, `MatCompositeType`
5406dbc55e5SJakub Kruzik @*/
541d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetType(Mat mat, MatCompositeType *type)
542d71ae5a4SJacob Faibussowitsch {
5436dbc55e5SJakub Kruzik   PetscFunctionBegin;
5446dbc55e5SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
5454f572ea9SToby Isaac   PetscAssertPointer(type, 2);
546cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetType_C", (Mat, MatCompositeType *), (mat, type));
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5486dbc55e5SJakub Kruzik }
5496dbc55e5SJakub Kruzik 
550d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat, MatStructure str)
551d71ae5a4SJacob Faibussowitsch {
552326b7573SPierre Jolivet   Mat_Composite *shell;
5533b35acafSJakub Kruzik 
5543b35acafSJakub Kruzik   PetscFunctionBegin;
555326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
556326b7573SPierre Jolivet   shell->structure = str;
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5583b35acafSJakub Kruzik }
5593b35acafSJakub Kruzik 
5603b35acafSJakub Kruzik /*@
5613b35acafSJakub Kruzik   MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix.
5623b35acafSJakub Kruzik 
5633b35acafSJakub Kruzik   Not Collective
5643b35acafSJakub Kruzik 
5653b35acafSJakub Kruzik   Input Parameters:
5663b35acafSJakub Kruzik + mat - the composite matrix
56711a5261eSBarry Smith - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN` (default) or `SUBSET_NONZERO_PATTERN`
5683b35acafSJakub Kruzik 
5693b35acafSJakub Kruzik   Level: advanced
5703b35acafSJakub Kruzik 
57111a5261eSBarry Smith   Note:
57211a5261eSBarry Smith   Information about the matrices structure is used in `MatCompositeMerge()` for additive composite matrix.
5733b35acafSJakub Kruzik 
5741cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE`
5753b35acafSJakub Kruzik @*/
576d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMatStructure(Mat mat, MatStructure str)
577d71ae5a4SJacob Faibussowitsch {
5783b35acafSJakub Kruzik   PetscFunctionBegin;
5793b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
580cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMatStructure_C", (Mat, MatStructure), (mat, str));
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5823b35acafSJakub Kruzik }
5833b35acafSJakub Kruzik 
584d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat, MatStructure *str)
585d71ae5a4SJacob Faibussowitsch {
586326b7573SPierre Jolivet   Mat_Composite *shell;
5873b35acafSJakub Kruzik 
5883b35acafSJakub Kruzik   PetscFunctionBegin;
589326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
590326b7573SPierre Jolivet   *str = shell->structure;
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5923b35acafSJakub Kruzik }
5933b35acafSJakub Kruzik 
5943b35acafSJakub Kruzik /*@
5953b35acafSJakub Kruzik   MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix.
5963b35acafSJakub Kruzik 
5973b35acafSJakub Kruzik   Not Collective
5983b35acafSJakub Kruzik 
5993b35acafSJakub Kruzik   Input Parameter:
6003b35acafSJakub Kruzik . mat - the composite matrix
6013b35acafSJakub Kruzik 
6023b35acafSJakub Kruzik   Output Parameter:
6033b35acafSJakub Kruzik . str - structure of the matrices
6043b35acafSJakub Kruzik 
6053b35acafSJakub Kruzik   Level: advanced
6063b35acafSJakub Kruzik 
6071cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE`
6083b35acafSJakub Kruzik @*/
609d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMatStructure(Mat mat, MatStructure *str)
610d71ae5a4SJacob Faibussowitsch {
6113b35acafSJakub Kruzik   PetscFunctionBegin;
6123b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6134f572ea9SToby Isaac   PetscAssertPointer(str, 2);
614cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMatStructure_C", (Mat, MatStructure *), (mat, str));
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6163b35acafSJakub Kruzik }
6173b35acafSJakub Kruzik 
618d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat, MatCompositeMergeType type)
619d71ae5a4SJacob Faibussowitsch {
620326b7573SPierre Jolivet   Mat_Composite *shell;
6218c8613bfSJakub Kruzik 
6228c8613bfSJakub Kruzik   PetscFunctionBegin;
623326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
6248c8613bfSJakub Kruzik   shell->mergetype = type;
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6268c8613bfSJakub Kruzik }
6278c8613bfSJakub Kruzik 
6288c8613bfSJakub Kruzik /*@
62911a5261eSBarry Smith   MatCompositeSetMergeType - Sets order of `MatCompositeMerge()`.
6308c8613bfSJakub Kruzik 
631c3339decSBarry Smith   Logically Collective
6328c8613bfSJakub Kruzik 
633d8d19677SJose E. Roman   Input Parameters:
6348c8613bfSJakub Kruzik + mat  - the composite matrix
63511a5261eSBarry Smith - type - `MAT_COMPOSITE_MERGE RIGHT` (default) to start merge from right with the first added matrix (mat[0]),
63611a5261eSBarry Smith           `MAT_COMPOSITE_MERGE_LEFT` to start merge from left with the last added matrix (mat[nmat-1])
6378c8613bfSJakub Kruzik 
6388c8613bfSJakub Kruzik   Level: advanced
6398c8613bfSJakub Kruzik 
64011a5261eSBarry Smith   Note:
641da81f932SPierre Jolivet   The resulting matrix is the same regardless of the `MatCompositeMergeType`. Only the order of operation is changed.
64211a5261eSBarry Smith   If set to `MAT_COMPOSITE_MERGE_RIGHT` the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0])))
6438c8613bfSJakub Kruzik   otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0].
6448c8613bfSJakub Kruzik 
6451cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE`
6468c8613bfSJakub Kruzik @*/
647d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMergeType(Mat mat, MatCompositeMergeType type)
648d71ae5a4SJacob Faibussowitsch {
6498c8613bfSJakub Kruzik   PetscFunctionBegin;
6508c8613bfSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6518c8613bfSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
652cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMergeType_C", (Mat, MatCompositeMergeType), (mat, type));
6533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6548c8613bfSJakub Kruzik }
6558c8613bfSJakub Kruzik 
656d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeMerge_Composite(Mat mat)
657d71ae5a4SJacob Faibussowitsch {
658326b7573SPierre Jolivet   Mat_Composite    *shell;
659326b7573SPierre Jolivet   Mat_CompositeLink next, prev;
6606d7c1e57SBarry Smith   Mat               tmat, newmat;
661326b7573SPierre Jolivet   Vec               left, right, dshift;
662326b7573SPierre Jolivet   PetscScalar       scale, shift;
66303049c21SJunchao Zhang   PetscInt          i;
664b52f573bSBarry Smith 
665b52f573bSBarry Smith   PetscFunctionBegin;
666326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
667326b7573SPierre Jolivet   next = shell->head;
668326b7573SPierre Jolivet   prev = shell->tail;
66928b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
670*b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(mat, &shift, &scale, &dshift, &left, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
6716d7c1e57SBarry Smith   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
6728c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
67303049c21SJunchao Zhang       i = 0;
6749566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
6759566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i++]));
67648a46eb9SPierre Jolivet       while ((next = next->next)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i++] : 1.0), next->mat, shell->structure));
6776d7c1e57SBarry Smith     } else {
67803049c21SJunchao Zhang       i = shell->nmat - 1;
6799566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
6809566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i--]));
68148a46eb9SPierre Jolivet       while ((prev = prev->prev)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i--] : 1.0), prev->mat, shell->structure));
6828c8613bfSJakub Kruzik     }
6838c8613bfSJakub Kruzik   } else {
6848c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
6859566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
686b6cfcaf5SJakub Kruzik       while ((next = next->next)) {
6879566063dSJacob Faibussowitsch         PetscCall(MatMatMult(next->mat, tmat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat));
6889566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
6896d7c1e57SBarry Smith         tmat = newmat;
6906d7c1e57SBarry Smith       }
69104d534ceSJakub Kruzik     } else {
6929566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
69304d534ceSJakub Kruzik       while ((prev = prev->prev)) {
6949566063dSJacob Faibussowitsch         PetscCall(MatMatMult(tmat, prev->mat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat));
6959566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
69604d534ceSJakub Kruzik         tmat = newmat;
69704d534ceSJakub Kruzik       }
69804d534ceSJakub Kruzik     }
6999371c9d4SSatish Balay     if (shell->scalings) {
7009371c9d4SSatish Balay       for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
7019371c9d4SSatish Balay     }
7026d7c1e57SBarry Smith   }
7031ba60692SJed Brown 
704*b9c875b8SPierre Jolivet   if (left) PetscCall(PetscObjectReference((PetscObject)left));
705*b9c875b8SPierre Jolivet   if (right) PetscCall(PetscObjectReference((PetscObject)right));
706*b9c875b8SPierre Jolivet   if (dshift) PetscCall(PetscObjectReference((PetscObject)dshift));
7071ba60692SJed Brown 
7089566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(mat, &tmat));
7091ba60692SJed Brown 
7109566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(mat, left, right));
7119566063dSJacob Faibussowitsch   PetscCall(MatScale(mat, scale));
712326b7573SPierre Jolivet   PetscCall(MatShift(mat, shift));
7139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
7149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
715326b7573SPierre Jolivet   if (dshift) {
716326b7573SPierre Jolivet     PetscCall(MatDiagonalSet(mat, dshift, ADD_VALUES));
717326b7573SPierre Jolivet     PetscCall(VecDestroy(&dshift));
718326b7573SPierre Jolivet   }
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
720b52f573bSBarry Smith }
7216a7ede75SJakub Kruzik 
7226a7ede75SJakub Kruzik /*@
723d7f81bf2SJakub Kruzik   MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
7248bd739bdSJakub Kruzik   by summing or computing the product of all the matrices inside the composite matrix.
725d7f81bf2SJakub Kruzik 
726b2b245abSJakub Kruzik   Collective
727d7f81bf2SJakub Kruzik 
728f899ff85SJose E. Roman   Input Parameter:
729d7f81bf2SJakub Kruzik . mat - the composite matrix
730d7f81bf2SJakub Kruzik 
7314b2558d6SJakub Kruzik   Options Database Keys:
73211a5261eSBarry Smith + -mat_composite_merge      - merge in `MatAssemblyEnd()`
733b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction
734d7f81bf2SJakub Kruzik 
735d7f81bf2SJakub Kruzik   Level: advanced
736d7f81bf2SJakub Kruzik 
73711a5261eSBarry Smith   Note:
73811a5261eSBarry Smith   The `MatType` of the resulting matrix will be the same as the `MatType` of the FIRST matrix in the composite matrix.
739d7f81bf2SJakub Kruzik 
7401cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE`
741d7f81bf2SJakub Kruzik @*/
742d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeMerge(Mat mat)
743d71ae5a4SJacob Faibussowitsch {
744d7f81bf2SJakub Kruzik   PetscFunctionBegin;
745d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
746cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeMerge_C", (Mat), (mat));
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
748d7f81bf2SJakub Kruzik }
749d7f81bf2SJakub Kruzik 
750d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat, PetscInt *nmat)
751d71ae5a4SJacob Faibussowitsch {
752326b7573SPierre Jolivet   Mat_Composite *shell;
753d7f81bf2SJakub Kruzik 
754d7f81bf2SJakub Kruzik   PetscFunctionBegin;
755326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
756d7f81bf2SJakub Kruzik   *nmat = shell->nmat;
7573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
758d7f81bf2SJakub Kruzik }
759d7f81bf2SJakub Kruzik 
760d7f81bf2SJakub Kruzik /*@
7616d0add67SJakub Kruzik   MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix.
7626a7ede75SJakub Kruzik 
7636a7ede75SJakub Kruzik   Not Collective
7646a7ede75SJakub Kruzik 
7656a7ede75SJakub Kruzik   Input Parameter:
766d7f81bf2SJakub Kruzik . mat - the composite matrix
7676a7ede75SJakub Kruzik 
7686a7ede75SJakub Kruzik   Output Parameter:
7696d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix
7706a7ede75SJakub Kruzik 
7718b5c3584SJakub Kruzik   Level: advanced
7726a7ede75SJakub Kruzik 
7731cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
7746a7ede75SJakub Kruzik @*/
775d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetNumberMat(Mat mat, PetscInt *nmat)
776d71ae5a4SJacob Faibussowitsch {
7776a7ede75SJakub Kruzik   PetscFunctionBegin;
778d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
7794f572ea9SToby Isaac   PetscAssertPointer(nmat, 2);
780cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetNumberMat_C", (Mat, PetscInt *), (mat, nmat));
7813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
782d7f81bf2SJakub Kruzik }
783d7f81bf2SJakub Kruzik 
784d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMat_Composite(Mat mat, PetscInt i, Mat *Ai)
785d71ae5a4SJacob Faibussowitsch {
786326b7573SPierre Jolivet   Mat_Composite    *shell;
787d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink;
788d7f81bf2SJakub Kruzik   PetscInt          k;
789d7f81bf2SJakub Kruzik 
790d7f81bf2SJakub Kruzik   PetscFunctionBegin;
791326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
79208401ef6SPierre Jolivet   PetscCheck(i < shell->nmat, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT, i, shell->nmat);
793d7f81bf2SJakub Kruzik   ilink = shell->head;
794ad540459SPierre Jolivet   for (k = 0; k < i; k++) ilink = ilink->next;
795d7f81bf2SJakub Kruzik   *Ai = ilink->mat;
7963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7976a7ede75SJakub Kruzik }
7986a7ede75SJakub Kruzik 
7998b5c3584SJakub Kruzik /*@
8008bd739bdSJakub Kruzik   MatCompositeGetMat - Returns the ith matrix from the composite matrix.
8018b5c3584SJakub Kruzik 
802c3339decSBarry Smith   Logically Collective
8038b5c3584SJakub Kruzik 
804d8d19677SJose E. Roman   Input Parameters:
805d7f81bf2SJakub Kruzik + mat - the composite matrix
8068b5c3584SJakub Kruzik - i   - the number of requested matrix
8078b5c3584SJakub Kruzik 
8088b5c3584SJakub Kruzik   Output Parameter:
8098b5c3584SJakub Kruzik . Ai - ith matrix in composite
8108b5c3584SJakub Kruzik 
8118b5c3584SJakub Kruzik   Level: advanced
8128b5c3584SJakub Kruzik 
8131cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE`
8148b5c3584SJakub Kruzik @*/
815d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMat(Mat mat, PetscInt i, Mat *Ai)
816d71ae5a4SJacob Faibussowitsch {
8178b5c3584SJakub Kruzik   PetscFunctionBegin;
818d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
819d7f81bf2SJakub Kruzik   PetscValidLogicalCollectiveInt(mat, i, 2);
8204f572ea9SToby Isaac   PetscAssertPointer(Ai, 3);
821cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMat_C", (Mat, PetscInt, Mat *), (mat, i, Ai));
8223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8238b5c3584SJakub Kruzik }
8248b5c3584SJakub Kruzik 
82566976f2fSJacob Faibussowitsch static PetscErrorCode MatCompositeSetScalings_Composite(Mat mat, const PetscScalar *scalings)
826d71ae5a4SJacob Faibussowitsch {
827326b7573SPierre Jolivet   Mat_Composite *shell;
82803049c21SJunchao Zhang   PetscInt       nmat;
82903049c21SJunchao Zhang 
83003049c21SJunchao Zhang   PetscFunctionBegin;
831326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
8329566063dSJacob Faibussowitsch   PetscCall(MatCompositeGetNumberMat(mat, &nmat));
8339566063dSJacob Faibussowitsch   if (!shell->scalings) PetscCall(PetscMalloc1(nmat, &shell->scalings));
8349566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(shell->scalings, scalings, nmat));
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83603049c21SJunchao Zhang }
83703049c21SJunchao Zhang 
83803049c21SJunchao Zhang /*@
83903049c21SJunchao Zhang   MatCompositeSetScalings - Sets separate scaling factors for component matrices.
84003049c21SJunchao Zhang 
841c3339decSBarry Smith   Logically Collective
84203049c21SJunchao Zhang 
843d8d19677SJose E. Roman   Input Parameters:
84403049c21SJunchao Zhang + mat      - the composite matrix
84503049c21SJunchao Zhang - scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat)
84603049c21SJunchao Zhang 
84703049c21SJunchao Zhang   Level: advanced
84803049c21SJunchao Zhang 
8491cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE`
85003049c21SJunchao Zhang @*/
851d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetScalings(Mat mat, const PetscScalar *scalings)
852d71ae5a4SJacob Faibussowitsch {
85303049c21SJunchao Zhang   PetscFunctionBegin;
85403049c21SJunchao Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
8554f572ea9SToby Isaac   PetscAssertPointer(scalings, 2);
85603049c21SJunchao Zhang   PetscValidLogicalCollectiveScalar(mat, *scalings, 2);
857cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetScalings_C", (Mat, const PetscScalar *), (mat, scalings));
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85903049c21SJunchao Zhang }
86003049c21SJunchao Zhang 
86141cd0125SJakub Kruzik /*MC
86241cd0125SJakub Kruzik    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices.
86341cd0125SJakub Kruzik     The matrices need to have a correct size and parallel layout for the sum or product to be valid.
86441cd0125SJakub Kruzik 
8652ef1f0ffSBarry Smith   Level: advanced
8662ef1f0ffSBarry Smith 
86711a5261eSBarry Smith    Note:
86811a5261eSBarry Smith    To use the product of the matrices call `MatCompositeSetType`(mat,`MAT_COMPOSITE_MULTIPLICATIVE`);
86941cd0125SJakub Kruzik 
870326b7573SPierre Jolivet   Developer Notes:
871326b7573SPierre Jolivet   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
872326b7573SPierre Jolivet 
873326b7573SPierre Jolivet   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
874326b7573SPierre Jolivet 
8751cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`,
87611a5261eSBarry Smith           `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()`
87741cd0125SJakub Kruzik M*/
87841cd0125SJakub Kruzik 
879d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A)
880d71ae5a4SJacob Faibussowitsch {
88141cd0125SJakub Kruzik   Mat_Composite *b;
88241cd0125SJakub Kruzik 
88341cd0125SJakub Kruzik   PetscFunctionBegin;
8844dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
88541cd0125SJakub Kruzik 
88641cd0125SJakub Kruzik   b->type        = MAT_COMPOSITE_ADDITIVE;
88741cd0125SJakub Kruzik   b->nmat        = 0;
8884b2558d6SJakub Kruzik   b->merge       = PETSC_FALSE;
8898c8613bfSJakub Kruzik   b->mergetype   = MAT_COMPOSITE_MERGE_RIGHT;
8903b35acafSJakub Kruzik   b->structure   = DIFFERENT_NONZERO_PATTERN;
89103049c21SJunchao Zhang   b->merge_mvctx = PETSC_TRUE;
89203049c21SJunchao Zhang 
893326b7573SPierre Jolivet   PetscCall(MatSetType(A, MATSHELL));
894326b7573SPierre Jolivet   PetscCall(MatShellSetContext(A, b));
895326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_Composite));
896326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Composite));
897326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Composite));
898326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Composite));
899326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_Composite));
900326b7573SPierre Jolivet   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_Composite));
9019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeAddMat_C", MatCompositeAddMat_Composite));
9029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetType_C", MatCompositeSetType_Composite));
9039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetType_C", MatCompositeGetType_Composite));
9049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMergeType_C", MatCompositeSetMergeType_Composite));
9059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMatStructure_C", MatCompositeSetMatStructure_Composite));
9069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMatStructure_C", MatCompositeGetMatStructure_Composite));
9079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeMerge_C", MatCompositeMerge_Composite));
9089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetNumberMat_C", MatCompositeGetNumberMat_Composite));
9099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMat_C", MatCompositeGetMat_Composite));
9109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetScalings_C", MatCompositeSetScalings_Composite));
911326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
912326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
913326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
9149566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCOMPOSITE));
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91641cd0125SJakub Kruzik }
917