xref: /petsc/src/mat/impls/composite/mcomposite.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
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 
MatDestroy_Composite(Mat mat)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 
MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y)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 
MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y)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 
MatMult_Composite(Mat mat,Vec x,Vec y)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 */
171f0b74427SPierre Jolivet     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]));
25857508eceSPierre Jolivet       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 
MatMultTranspose_Composite(Mat A,Vec x,Vec y)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 
MatGetDiagonal_Composite(Mat A,Vec v)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));
32057508eceSPierre Jolivet     PetscCall(VecAXPY(v, shell->scalings ? shell->scalings[i++] : 1.0, shell->work));
321793850ffSBarry Smith   }
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323793850ffSBarry Smith }
324793850ffSBarry Smith 
MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)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 
MatSetFromOptions_Composite(Mat A,PetscOptionItems PetscOptionsObject)336ce78bad3SBarry Smith 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 @*/
MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat * mats,Mat * mat)388d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateComposite(MPI_Comm comm, PetscInt nmat, const Mat *mats, Mat *mat)
389d71ae5a4SJacob Faibussowitsch {
390793850ffSBarry Smith   PetscFunctionBegin;
39108401ef6SPierre Jolivet   PetscCheck(nmat >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in at least one matrix");
3924f572ea9SToby Isaac   PetscAssertPointer(mat, 4);
3939566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, mat));
3949566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat, MATCOMPOSITE));
3957aeb3b90SRichard Tran Mills   for (PetscInt i = 0; i < nmat; i++) PetscCall(MatCompositeAddMat(*mat, mats[i]));
3969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
3979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
399793850ffSBarry Smith }
400793850ffSBarry Smith 
MatCompositeAddMat_Composite(Mat mat,Mat smat)401d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeAddMat_Composite(Mat mat, Mat smat)
402d71ae5a4SJacob Faibussowitsch {
403326b7573SPierre Jolivet   Mat_Composite    *shell;
404326b7573SPierre Jolivet   Mat_CompositeLink ilink, next;
4055c8dc79eSRichard Tran Mills   VecType           vtype_mat, vtype_smat;
4065c8dc79eSRichard Tran Mills   PetscBool         match;
407d7f81bf2SJakub Kruzik 
408d7f81bf2SJakub Kruzik   PetscFunctionBegin;
409326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
410326b7573SPierre Jolivet   next = shell->head;
4114dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
412f4259b30SLisandro Dalcin   ilink->next = NULL;
4139566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)smat));
414d7f81bf2SJakub Kruzik   ilink->mat = smat;
415d7f81bf2SJakub Kruzik 
416d7f81bf2SJakub Kruzik   if (!next) shell->head = ilink;
417d7f81bf2SJakub Kruzik   else {
418ad540459SPierre Jolivet     while (next->next) next = next->next;
419d7f81bf2SJakub Kruzik     next->next  = ilink;
420d7f81bf2SJakub Kruzik     ilink->prev = next;
421d7f81bf2SJakub Kruzik   }
422d7f81bf2SJakub Kruzik   shell->tail = ilink;
423d7f81bf2SJakub Kruzik   shell->nmat += 1;
42403049c21SJunchao Zhang 
4255c8dc79eSRichard Tran Mills   /* If all of the partial matrices have the same default vector type, then the composite matrix should also have this default type.
4265c8dc79eSRichard Tran Mills      Otherwise, the default type should be "standard". */
4275c8dc79eSRichard Tran Mills   PetscCall(MatGetVecType(smat, &vtype_smat));
4285c8dc79eSRichard Tran Mills   if (shell->nmat == 1) PetscCall(MatSetVecType(mat, vtype_smat));
4295c8dc79eSRichard Tran Mills   else {
4305c8dc79eSRichard Tran Mills     PetscCall(MatGetVecType(mat, &vtype_mat));
4315c8dc79eSRichard Tran Mills     PetscCall(PetscStrcmp(vtype_smat, vtype_mat, &match));
4325c8dc79eSRichard Tran Mills     if (!match) PetscCall(MatSetVecType(mat, VECSTANDARD));
4335c8dc79eSRichard Tran Mills   }
4345c8dc79eSRichard Tran Mills 
43503049c21SJunchao Zhang   /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */
43603049c21SJunchao Zhang   if (shell->scalings) {
4379566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(sizeof(PetscScalar) * shell->nmat, &shell->scalings));
43803049c21SJunchao Zhang     shell->scalings[shell->nmat - 1] = 1.0;
43903049c21SJunchao Zhang   }
4407aeb3b90SRichard Tran Mills 
4417aeb3b90SRichard Tran Mills   /* The composite matrix requires PetscLayouts for its rows and columns; we copy these from the constituent partial matrices. */
4427aeb3b90SRichard Tran Mills   if (shell->nmat == 1) PetscCall(PetscLayoutReference(smat->cmap, &mat->cmap));
4437aeb3b90SRichard Tran Mills   PetscCall(PetscLayoutReference(smat->rmap, &mat->rmap));
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
445d7f81bf2SJakub Kruzik }
446d7f81bf2SJakub Kruzik 
447793850ffSBarry Smith /*@
4488bd739bdSJakub Kruzik   MatCompositeAddMat - Add another matrix to a composite matrix.
449793850ffSBarry Smith 
450c3339decSBarry Smith   Collective
451793850ffSBarry Smith 
452793850ffSBarry Smith   Input Parameters:
453793850ffSBarry Smith + mat  - the composite matrix
454793850ffSBarry Smith - smat - the partial matrix
455793850ffSBarry Smith 
456793850ffSBarry Smith   Level: advanced
457793850ffSBarry Smith 
4581cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
459793850ffSBarry Smith @*/
MatCompositeAddMat(Mat mat,Mat smat)460d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeAddMat(Mat mat, Mat smat)
461d71ae5a4SJacob Faibussowitsch {
462793850ffSBarry Smith   PetscFunctionBegin;
4630700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
4640700a824SBarry Smith   PetscValidHeaderSpecific(smat, MAT_CLASSID, 2);
465cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeAddMat_C", (Mat, Mat), (mat, smat));
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467d7f81bf2SJakub Kruzik }
468793850ffSBarry Smith 
MatCompositeSetType_Composite(Mat mat,MatCompositeType type)469d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetType_Composite(Mat mat, MatCompositeType type)
470d71ae5a4SJacob Faibussowitsch {
471326b7573SPierre Jolivet   Mat_Composite *b;
472d7f81bf2SJakub Kruzik 
473d7f81bf2SJakub Kruzik   PetscFunctionBegin;
474326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &b));
475ced1ac25SJakub Kruzik   b->type = type;
476d7f81bf2SJakub Kruzik   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
477326b7573SPierre Jolivet     PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, NULL));
478*57d50842SBarry Smith     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Composite_Multiplicative));
479*57d50842SBarry Smith     PetscCall(MatShellSetOperation(mat, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Composite_Multiplicative));
48003049c21SJunchao Zhang     b->merge_mvctx = PETSC_FALSE;
481d7f81bf2SJakub Kruzik   } else {
482*57d50842SBarry Smith     PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_Composite));
483*57d50842SBarry Smith     PetscCall(MatShellSetOperation(mat, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Composite));
484*57d50842SBarry Smith     PetscCall(MatShellSetOperation(mat, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Composite));
485793850ffSBarry Smith   }
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4876d7c1e57SBarry Smith }
4886d7c1e57SBarry Smith 
4892c0821f3SBarry Smith /*@
4908bd739bdSJakub Kruzik   MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product.
4916d7c1e57SBarry Smith 
492c3339decSBarry Smith   Logically Collective
4936d7c1e57SBarry Smith 
4946d7c1e57SBarry Smith   Input Parameters:
49520f4b53cSBarry Smith + mat  - the composite matrix
49620f4b53cSBarry Smith - type - the `MatCompositeType` to use for the matrix
4976d7c1e57SBarry Smith 
4986d7c1e57SBarry Smith   Level: advanced
4996d7c1e57SBarry Smith 
5001cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE`,
50120f4b53cSBarry Smith           `MatCompositeType`
5026d7c1e57SBarry Smith @*/
MatCompositeSetType(Mat mat,MatCompositeType type)503d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetType(Mat mat, MatCompositeType type)
504d71ae5a4SJacob Faibussowitsch {
5056d7c1e57SBarry Smith   PetscFunctionBegin;
506d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
507b2b245abSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
508cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetType_C", (Mat, MatCompositeType), (mat, type));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
510793850ffSBarry Smith }
511793850ffSBarry Smith 
MatCompositeGetType_Composite(Mat mat,MatCompositeType * type)512d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetType_Composite(Mat mat, MatCompositeType *type)
513d71ae5a4SJacob Faibussowitsch {
514326b7573SPierre Jolivet   Mat_Composite *shell;
5156dbc55e5SJakub Kruzik 
5166dbc55e5SJakub Kruzik   PetscFunctionBegin;
517326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
518326b7573SPierre Jolivet   *type = shell->type;
5193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5206dbc55e5SJakub Kruzik }
5216dbc55e5SJakub Kruzik 
5226dbc55e5SJakub Kruzik /*@
5236dbc55e5SJakub Kruzik   MatCompositeGetType - Returns type of composite.
5246dbc55e5SJakub Kruzik 
5256dbc55e5SJakub Kruzik   Not Collective
5266dbc55e5SJakub Kruzik 
5276dbc55e5SJakub Kruzik   Input Parameter:
5286dbc55e5SJakub Kruzik . mat - the composite matrix
5296dbc55e5SJakub Kruzik 
5306dbc55e5SJakub Kruzik   Output Parameter:
5316dbc55e5SJakub Kruzik . type - type of composite
5326dbc55e5SJakub Kruzik 
5336dbc55e5SJakub Kruzik   Level: advanced
5346dbc55e5SJakub Kruzik 
5351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE`, `MatCompositeType`
5366dbc55e5SJakub Kruzik @*/
MatCompositeGetType(Mat mat,MatCompositeType * type)537d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetType(Mat mat, MatCompositeType *type)
538d71ae5a4SJacob Faibussowitsch {
5396dbc55e5SJakub Kruzik   PetscFunctionBegin;
5406dbc55e5SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
5414f572ea9SToby Isaac   PetscAssertPointer(type, 2);
542cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetType_C", (Mat, MatCompositeType *), (mat, type));
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5446dbc55e5SJakub Kruzik }
5456dbc55e5SJakub Kruzik 
MatCompositeSetMatStructure_Composite(Mat mat,MatStructure str)546d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat, MatStructure str)
547d71ae5a4SJacob Faibussowitsch {
548326b7573SPierre Jolivet   Mat_Composite *shell;
5493b35acafSJakub Kruzik 
5503b35acafSJakub Kruzik   PetscFunctionBegin;
551326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
552326b7573SPierre Jolivet   shell->structure = str;
5533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5543b35acafSJakub Kruzik }
5553b35acafSJakub Kruzik 
5563b35acafSJakub Kruzik /*@
5573b35acafSJakub Kruzik   MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix.
5583b35acafSJakub Kruzik 
5593b35acafSJakub Kruzik   Not Collective
5603b35acafSJakub Kruzik 
5613b35acafSJakub Kruzik   Input Parameters:
5623b35acafSJakub Kruzik + mat - the composite matrix
56311a5261eSBarry Smith - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN` (default) or `SUBSET_NONZERO_PATTERN`
5643b35acafSJakub Kruzik 
5653b35acafSJakub Kruzik   Level: advanced
5663b35acafSJakub Kruzik 
56711a5261eSBarry Smith   Note:
56811a5261eSBarry Smith   Information about the matrices structure is used in `MatCompositeMerge()` for additive composite matrix.
5693b35acafSJakub Kruzik 
5701cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE`
5713b35acafSJakub Kruzik @*/
MatCompositeSetMatStructure(Mat mat,MatStructure str)572d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMatStructure(Mat mat, MatStructure str)
573d71ae5a4SJacob Faibussowitsch {
5743b35acafSJakub Kruzik   PetscFunctionBegin;
5753b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
576cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMatStructure_C", (Mat, MatStructure), (mat, str));
5773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5783b35acafSJakub Kruzik }
5793b35acafSJakub Kruzik 
MatCompositeGetMatStructure_Composite(Mat mat,MatStructure * str)580d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat, MatStructure *str)
581d71ae5a4SJacob Faibussowitsch {
582326b7573SPierre Jolivet   Mat_Composite *shell;
5833b35acafSJakub Kruzik 
5843b35acafSJakub Kruzik   PetscFunctionBegin;
585326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
586326b7573SPierre Jolivet   *str = shell->structure;
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5883b35acafSJakub Kruzik }
5893b35acafSJakub Kruzik 
5903b35acafSJakub Kruzik /*@
5913b35acafSJakub Kruzik   MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix.
5923b35acafSJakub Kruzik 
5933b35acafSJakub Kruzik   Not Collective
5943b35acafSJakub Kruzik 
5953b35acafSJakub Kruzik   Input Parameter:
5963b35acafSJakub Kruzik . mat - the composite matrix
5973b35acafSJakub Kruzik 
5983b35acafSJakub Kruzik   Output Parameter:
5993b35acafSJakub Kruzik . str - structure of the matrices
6003b35acafSJakub Kruzik 
6013b35acafSJakub Kruzik   Level: advanced
6023b35acafSJakub Kruzik 
6031cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE`
6043b35acafSJakub Kruzik @*/
MatCompositeGetMatStructure(Mat mat,MatStructure * str)605d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMatStructure(Mat mat, MatStructure *str)
606d71ae5a4SJacob Faibussowitsch {
6073b35acafSJakub Kruzik   PetscFunctionBegin;
6083b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6094f572ea9SToby Isaac   PetscAssertPointer(str, 2);
610cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMatStructure_C", (Mat, MatStructure *), (mat, str));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6123b35acafSJakub Kruzik }
6133b35acafSJakub Kruzik 
MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type)614d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat, MatCompositeMergeType type)
615d71ae5a4SJacob Faibussowitsch {
616326b7573SPierre Jolivet   Mat_Composite *shell;
6178c8613bfSJakub Kruzik 
6188c8613bfSJakub Kruzik   PetscFunctionBegin;
619326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
6208c8613bfSJakub Kruzik   shell->mergetype = type;
6213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6228c8613bfSJakub Kruzik }
6238c8613bfSJakub Kruzik 
6248c8613bfSJakub Kruzik /*@
62511a5261eSBarry Smith   MatCompositeSetMergeType - Sets order of `MatCompositeMerge()`.
6268c8613bfSJakub Kruzik 
627c3339decSBarry Smith   Logically Collective
6288c8613bfSJakub Kruzik 
629d8d19677SJose E. Roman   Input Parameters:
6308c8613bfSJakub Kruzik + mat  - the composite matrix
63111a5261eSBarry Smith - type - `MAT_COMPOSITE_MERGE RIGHT` (default) to start merge from right with the first added matrix (mat[0]),
63211a5261eSBarry Smith           `MAT_COMPOSITE_MERGE_LEFT` to start merge from left with the last added matrix (mat[nmat-1])
6338c8613bfSJakub Kruzik 
6348c8613bfSJakub Kruzik   Level: advanced
6358c8613bfSJakub Kruzik 
63611a5261eSBarry Smith   Note:
637da81f932SPierre Jolivet   The resulting matrix is the same regardless of the `MatCompositeMergeType`. Only the order of operation is changed.
63811a5261eSBarry Smith   If set to `MAT_COMPOSITE_MERGE_RIGHT` the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0])))
6398c8613bfSJakub Kruzik   otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0].
6408c8613bfSJakub Kruzik 
6411cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE`
6428c8613bfSJakub Kruzik @*/
MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type)643d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetMergeType(Mat mat, MatCompositeMergeType type)
644d71ae5a4SJacob Faibussowitsch {
6458c8613bfSJakub Kruzik   PetscFunctionBegin;
6468c8613bfSJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
6478c8613bfSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat, type, 2);
648cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetMergeType_C", (Mat, MatCompositeMergeType), (mat, type));
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6508c8613bfSJakub Kruzik }
6518c8613bfSJakub Kruzik 
MatCompositeMerge_Composite(Mat mat)652d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeMerge_Composite(Mat mat)
653d71ae5a4SJacob Faibussowitsch {
654326b7573SPierre Jolivet   Mat_Composite    *shell;
655326b7573SPierre Jolivet   Mat_CompositeLink next, prev;
6566d7c1e57SBarry Smith   Mat               tmat, newmat;
657326b7573SPierre Jolivet   Vec               left, right, dshift;
658326b7573SPierre Jolivet   PetscScalar       scale, shift;
65903049c21SJunchao Zhang   PetscInt          i;
660b52f573bSBarry Smith 
661b52f573bSBarry Smith   PetscFunctionBegin;
662326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
663326b7573SPierre Jolivet   next = shell->head;
664326b7573SPierre Jolivet   prev = shell->tail;
66528b400f6SJacob Faibussowitsch   PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()");
666b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(mat, &shift, &scale, &dshift, &left, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
6676d7c1e57SBarry Smith   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
6688c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
66903049c21SJunchao Zhang       i = 0;
6709566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
6719566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i++]));
67257508eceSPierre Jolivet       while ((next = next->next)) PetscCall(MatAXPY(tmat, shell->scalings ? shell->scalings[i++] : 1.0, next->mat, shell->structure));
6736d7c1e57SBarry Smith     } else {
67403049c21SJunchao Zhang       i = shell->nmat - 1;
6759566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
6769566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i--]));
67757508eceSPierre Jolivet       while ((prev = prev->prev)) PetscCall(MatAXPY(tmat, shell->scalings ? shell->scalings[i--] : 1.0, prev->mat, shell->structure));
6788c8613bfSJakub Kruzik     }
6798c8613bfSJakub Kruzik   } else {
6808c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
6819566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat));
682b6cfcaf5SJakub Kruzik       while ((next = next->next)) {
683fb842aefSJose E. Roman         PetscCall(MatMatMult(next->mat, tmat, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &newmat));
6849566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
6856d7c1e57SBarry Smith         tmat = newmat;
6866d7c1e57SBarry Smith       }
68704d534ceSJakub Kruzik     } else {
6889566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat));
68904d534ceSJakub Kruzik       while ((prev = prev->prev)) {
690fb842aefSJose E. Roman         PetscCall(MatMatMult(tmat, prev->mat, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &newmat));
6919566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
69204d534ceSJakub Kruzik         tmat = newmat;
69304d534ceSJakub Kruzik       }
69404d534ceSJakub Kruzik     }
6959371c9d4SSatish Balay     if (shell->scalings) {
6969371c9d4SSatish Balay       for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i];
6979371c9d4SSatish Balay     }
6986d7c1e57SBarry Smith   }
6991ba60692SJed Brown 
700b9c875b8SPierre Jolivet   if (left) PetscCall(PetscObjectReference((PetscObject)left));
701b9c875b8SPierre Jolivet   if (right) PetscCall(PetscObjectReference((PetscObject)right));
702b9c875b8SPierre Jolivet   if (dshift) PetscCall(PetscObjectReference((PetscObject)dshift));
7031ba60692SJed Brown 
7049566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(mat, &tmat));
7051ba60692SJed Brown 
7069566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(mat, left, right));
7079566063dSJacob Faibussowitsch   PetscCall(MatScale(mat, scale));
708326b7573SPierre Jolivet   PetscCall(MatShift(mat, shift));
7099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
7109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
711326b7573SPierre Jolivet   if (dshift) {
712326b7573SPierre Jolivet     PetscCall(MatDiagonalSet(mat, dshift, ADD_VALUES));
713326b7573SPierre Jolivet     PetscCall(VecDestroy(&dshift));
714326b7573SPierre Jolivet   }
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
716b52f573bSBarry Smith }
7176a7ede75SJakub Kruzik 
7186a7ede75SJakub Kruzik /*@
719d7f81bf2SJakub Kruzik   MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
7208bd739bdSJakub Kruzik   by summing or computing the product of all the matrices inside the composite matrix.
721d7f81bf2SJakub Kruzik 
722b2b245abSJakub Kruzik   Collective
723d7f81bf2SJakub Kruzik 
724f899ff85SJose E. Roman   Input Parameter:
725d7f81bf2SJakub Kruzik . mat - the composite matrix
726d7f81bf2SJakub Kruzik 
7274b2558d6SJakub Kruzik   Options Database Keys:
72811a5261eSBarry Smith + -mat_composite_merge      - merge in `MatAssemblyEnd()`
729b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction
730d7f81bf2SJakub Kruzik 
731d7f81bf2SJakub Kruzik   Level: advanced
732d7f81bf2SJakub Kruzik 
73311a5261eSBarry Smith   Note:
73411a5261eSBarry Smith   The `MatType` of the resulting matrix will be the same as the `MatType` of the FIRST matrix in the composite matrix.
735d7f81bf2SJakub Kruzik 
7361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE`
737d7f81bf2SJakub Kruzik @*/
MatCompositeMerge(Mat mat)738d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeMerge(Mat mat)
739d71ae5a4SJacob Faibussowitsch {
740d7f81bf2SJakub Kruzik   PetscFunctionBegin;
741d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
742cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeMerge_C", (Mat), (mat));
7433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
744d7f81bf2SJakub Kruzik }
745d7f81bf2SJakub Kruzik 
MatCompositeGetNumberMat_Composite(Mat mat,PetscInt * nmat)746d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat, PetscInt *nmat)
747d71ae5a4SJacob Faibussowitsch {
748326b7573SPierre Jolivet   Mat_Composite *shell;
749d7f81bf2SJakub Kruzik 
750d7f81bf2SJakub Kruzik   PetscFunctionBegin;
751326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
752d7f81bf2SJakub Kruzik   *nmat = shell->nmat;
7533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
754d7f81bf2SJakub Kruzik }
755d7f81bf2SJakub Kruzik 
756d7f81bf2SJakub Kruzik /*@
7576d0add67SJakub Kruzik   MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix.
7586a7ede75SJakub Kruzik 
7596a7ede75SJakub Kruzik   Not Collective
7606a7ede75SJakub Kruzik 
7616a7ede75SJakub Kruzik   Input Parameter:
762d7f81bf2SJakub Kruzik . mat - the composite matrix
7636a7ede75SJakub Kruzik 
7646a7ede75SJakub Kruzik   Output Parameter:
7656d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix
7666a7ede75SJakub Kruzik 
7678b5c3584SJakub Kruzik   Level: advanced
7686a7ede75SJakub Kruzik 
7691cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
7706a7ede75SJakub Kruzik @*/
MatCompositeGetNumberMat(Mat mat,PetscInt * nmat)771d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetNumberMat(Mat mat, PetscInt *nmat)
772d71ae5a4SJacob Faibussowitsch {
7736a7ede75SJakub Kruzik   PetscFunctionBegin;
774d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
7754f572ea9SToby Isaac   PetscAssertPointer(nmat, 2);
776cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetNumberMat_C", (Mat, PetscInt *), (mat, nmat));
7773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
778d7f81bf2SJakub Kruzik }
779d7f81bf2SJakub Kruzik 
MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat * Ai)780d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCompositeGetMat_Composite(Mat mat, PetscInt i, Mat *Ai)
781d71ae5a4SJacob Faibussowitsch {
782326b7573SPierre Jolivet   Mat_Composite    *shell;
783d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink;
784d7f81bf2SJakub Kruzik   PetscInt          k;
785d7f81bf2SJakub Kruzik 
786d7f81bf2SJakub Kruzik   PetscFunctionBegin;
787326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
78808401ef6SPierre Jolivet   PetscCheck(i < shell->nmat, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT, i, shell->nmat);
789d7f81bf2SJakub Kruzik   ilink = shell->head;
790ad540459SPierre Jolivet   for (k = 0; k < i; k++) ilink = ilink->next;
791d7f81bf2SJakub Kruzik   *Ai = ilink->mat;
7923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7936a7ede75SJakub Kruzik }
7946a7ede75SJakub Kruzik 
7958b5c3584SJakub Kruzik /*@
7968bd739bdSJakub Kruzik   MatCompositeGetMat - Returns the ith matrix from the composite matrix.
7978b5c3584SJakub Kruzik 
798c3339decSBarry Smith   Logically Collective
7998b5c3584SJakub Kruzik 
800d8d19677SJose E. Roman   Input Parameters:
801d7f81bf2SJakub Kruzik + mat - the composite matrix
8028b5c3584SJakub Kruzik - i   - the number of requested matrix
8038b5c3584SJakub Kruzik 
8048b5c3584SJakub Kruzik   Output Parameter:
8058b5c3584SJakub Kruzik . Ai - ith matrix in composite
8068b5c3584SJakub Kruzik 
8078b5c3584SJakub Kruzik   Level: advanced
8088b5c3584SJakub Kruzik 
8091cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE`
8108b5c3584SJakub Kruzik @*/
MatCompositeGetMat(Mat mat,PetscInt i,Mat * Ai)811d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeGetMat(Mat mat, PetscInt i, Mat *Ai)
812d71ae5a4SJacob Faibussowitsch {
8138b5c3584SJakub Kruzik   PetscFunctionBegin;
814d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
815d7f81bf2SJakub Kruzik   PetscValidLogicalCollectiveInt(mat, i, 2);
8164f572ea9SToby Isaac   PetscAssertPointer(Ai, 3);
817cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeGetMat_C", (Mat, PetscInt, Mat *), (mat, i, Ai));
8183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8198b5c3584SJakub Kruzik }
8208b5c3584SJakub Kruzik 
MatCompositeSetScalings_Composite(Mat mat,const PetscScalar * scalings)82166976f2fSJacob Faibussowitsch static PetscErrorCode MatCompositeSetScalings_Composite(Mat mat, const PetscScalar *scalings)
822d71ae5a4SJacob Faibussowitsch {
823326b7573SPierre Jolivet   Mat_Composite *shell;
82403049c21SJunchao Zhang   PetscInt       nmat;
82503049c21SJunchao Zhang 
82603049c21SJunchao Zhang   PetscFunctionBegin;
827326b7573SPierre Jolivet   PetscCall(MatShellGetContext(mat, &shell));
8289566063dSJacob Faibussowitsch   PetscCall(MatCompositeGetNumberMat(mat, &nmat));
8299566063dSJacob Faibussowitsch   if (!shell->scalings) PetscCall(PetscMalloc1(nmat, &shell->scalings));
8309566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(shell->scalings, scalings, nmat));
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83203049c21SJunchao Zhang }
83303049c21SJunchao Zhang 
83403049c21SJunchao Zhang /*@
83503049c21SJunchao Zhang   MatCompositeSetScalings - Sets separate scaling factors for component matrices.
83603049c21SJunchao Zhang 
837c3339decSBarry Smith   Logically Collective
83803049c21SJunchao Zhang 
839d8d19677SJose E. Roman   Input Parameters:
84003049c21SJunchao Zhang + mat      - the composite matrix
84103049c21SJunchao Zhang - scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat)
84203049c21SJunchao Zhang 
84303049c21SJunchao Zhang   Level: advanced
84403049c21SJunchao Zhang 
8451cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE`
84603049c21SJunchao Zhang @*/
MatCompositeSetScalings(Mat mat,const PetscScalar * scalings)847d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCompositeSetScalings(Mat mat, const PetscScalar *scalings)
848d71ae5a4SJacob Faibussowitsch {
84903049c21SJunchao Zhang   PetscFunctionBegin;
85003049c21SJunchao Zhang   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
8514f572ea9SToby Isaac   PetscAssertPointer(scalings, 2);
85203049c21SJunchao Zhang   PetscValidLogicalCollectiveScalar(mat, *scalings, 2);
853cac4c232SBarry Smith   PetscUseMethod(mat, "MatCompositeSetScalings_C", (Mat, const PetscScalar *), (mat, scalings));
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85503049c21SJunchao Zhang }
85603049c21SJunchao Zhang 
85741cd0125SJakub Kruzik /*MC
85841cd0125SJakub Kruzik    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices.
85941cd0125SJakub Kruzik     The matrices need to have a correct size and parallel layout for the sum or product to be valid.
86041cd0125SJakub Kruzik 
8612ef1f0ffSBarry Smith   Level: advanced
8622ef1f0ffSBarry Smith 
86311a5261eSBarry Smith    Note:
86411a5261eSBarry Smith    To use the product of the matrices call `MatCompositeSetType`(mat,`MAT_COMPOSITE_MULTIPLICATIVE`);
86541cd0125SJakub Kruzik 
866326b7573SPierre Jolivet   Developer Notes:
867326b7573SPierre Jolivet   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
868326b7573SPierre Jolivet 
869326b7573SPierre Jolivet   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
870326b7573SPierre Jolivet 
8711cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`,
87211a5261eSBarry Smith           `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()`
87341cd0125SJakub Kruzik M*/
87441cd0125SJakub Kruzik 
MatCreate_Composite(Mat A)875d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A)
876d71ae5a4SJacob Faibussowitsch {
87741cd0125SJakub Kruzik   Mat_Composite *b;
87841cd0125SJakub Kruzik 
87941cd0125SJakub Kruzik   PetscFunctionBegin;
8804dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
88141cd0125SJakub Kruzik 
88241cd0125SJakub Kruzik   b->type        = MAT_COMPOSITE_ADDITIVE;
88341cd0125SJakub Kruzik   b->nmat        = 0;
8844b2558d6SJakub Kruzik   b->merge       = PETSC_FALSE;
8858c8613bfSJakub Kruzik   b->mergetype   = MAT_COMPOSITE_MERGE_RIGHT;
8863b35acafSJakub Kruzik   b->structure   = DIFFERENT_NONZERO_PATTERN;
88703049c21SJunchao Zhang   b->merge_mvctx = PETSC_TRUE;
88803049c21SJunchao Zhang 
889326b7573SPierre Jolivet   PetscCall(MatSetType(A, MATSHELL));
890326b7573SPierre Jolivet   PetscCall(MatShellSetContext(A, b));
891*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_Composite));
892*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Composite));
893*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Composite));
894*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_Composite));
895*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (PetscErrorCodeFn *)MatAssemblyEnd_Composite));
896*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (PetscErrorCodeFn *)MatSetFromOptions_Composite));
8979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeAddMat_C", MatCompositeAddMat_Composite));
8989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetType_C", MatCompositeSetType_Composite));
8999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetType_C", MatCompositeGetType_Composite));
9009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMergeType_C", MatCompositeSetMergeType_Composite));
9019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMatStructure_C", MatCompositeSetMatStructure_Composite));
9029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMatStructure_C", MatCompositeGetMatStructure_Composite));
9039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeMerge_C", MatCompositeMerge_Composite));
9049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetNumberMat_C", MatCompositeGetNumberMat_Composite));
9059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMat_C", MatCompositeGetMat_Composite));
9069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetScalings_C", MatCompositeSetScalings_Composite));
907326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
908326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
909326b7573SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
9109566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCOMPOSITE));
9113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91241cd0125SJakub Kruzik }
913