xref: /petsc/src/mat/impls/composite/mcomposite.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
1793850ffSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
3793850ffSBarry Smith 
4f4259b30SLisandro Dalcin const char *const MatCompositeMergeTypes[] = {"left","right","MatCompositeMergeType","MAT_COMPOSITE_",NULL};
58c8613bfSJakub Kruzik 
6793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink;
7793850ffSBarry Smith struct _Mat_CompositeLink {
8793850ffSBarry Smith   Mat               mat;
96d7c1e57SBarry Smith   Vec               work;
106d7c1e57SBarry Smith   Mat_CompositeLink next,prev;
11793850ffSBarry Smith };
12793850ffSBarry Smith 
13793850ffSBarry Smith typedef struct {
146d7c1e57SBarry Smith   MatCompositeType      type;
156d7c1e57SBarry Smith   Mat_CompositeLink     head,tail;
16793850ffSBarry Smith   Vec                   work;
17e4fc5df0SBarry Smith   PetscScalar           scale;        /* scale factor supplied with MatScale() */
18e4fc5df0SBarry Smith   Vec                   left,right;   /* left and right diagonal scaling provided with MatDiagonalScale() */
1903049c21SJunchao Zhang   Vec                   leftwork,rightwork,leftwork2,rightwork2; /* Two pairs of working vectors */
206a7ede75SJakub Kruzik   PetscInt              nmat;
214b2558d6SJakub Kruzik   PetscBool             merge;
228c8613bfSJakub Kruzik   MatCompositeMergeType mergetype;
233b35acafSJakub Kruzik   MatStructure          structure;
2403049c21SJunchao Zhang 
2503049c21SJunchao Zhang   PetscScalar           *scalings;
2603049c21SJunchao Zhang   PetscBool             merge_mvctx;  /* Whether need to merge mvctx of component matrices */
2703049c21SJunchao Zhang   Vec                   *lvecs;       /* [nmat] Basically, they are Mvctx->lvec of each component matrix */
2803049c21SJunchao Zhang   PetscScalar           *larray;      /* [len] Data arrays of lvecs[] are stored consecutively in larray */
2903049c21SJunchao Zhang   PetscInt              len;          /* Length of larray[] */
3003049c21SJunchao Zhang   Vec                   gvec;         /* Union of lvecs[] without duplicated entries */
3103049c21SJunchao Zhang   PetscInt              *location;    /* A map that maps entries in garray[] to larray[] */
3203049c21SJunchao Zhang   VecScatter            Mvctx;
33793850ffSBarry Smith } Mat_Composite;
34793850ffSBarry Smith 
35793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat)
36793850ffSBarry Smith {
372c33bbaeSSatish Balay   Mat_Composite     *shell = (Mat_Composite*)mat->data;
386d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head,oldnext;
3903049c21SJunchao Zhang   PetscInt          i;
40793850ffSBarry Smith 
41793850ffSBarry Smith   PetscFunctionBegin;
42793850ffSBarry Smith   while (next) {
439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&next->mat));
446d7c1e57SBarry Smith     if (next->work && (!next->next || next->work != next->next->work)) {
459566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&next->work));
466d7c1e57SBarry Smith     }
476d7c1e57SBarry Smith     oldnext = next;
48793850ffSBarry Smith     next    = next->next;
499566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldnext));
50793850ffSBarry Smith   }
519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->work));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->left));
539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->right));
549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->leftwork));
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->rightwork));
569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->leftwork2));
579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&shell->rightwork2));
5803049c21SJunchao Zhang 
5903049c21SJunchao Zhang   if (shell->Mvctx) {
609566063dSJacob Faibussowitsch     for (i=0; i<shell->nmat; i++) PetscCall(VecDestroy(&shell->lvecs[i]));
619566063dSJacob Faibussowitsch     PetscCall(PetscFree3(shell->location,shell->larray,shell->lvecs));
629566063dSJacob Faibussowitsch     PetscCall(PetscFree(shell->larray));
639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&shell->gvec));
649566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&shell->Mvctx));
6503049c21SJunchao Zhang   }
6603049c21SJunchao Zhang 
679566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell->scalings));
68*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeAddMat_C",NULL));
69*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetType_C",NULL));
70*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetType_C",NULL));
71*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetMergeType_C",NULL));
72*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetMatStructure_C",NULL));
73*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetMatStructure_C",NULL));
74*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeMerge_C",NULL));
75*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetNumberMat_C",NULL));
76*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetMat_C",NULL));
77*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetScalings_C",NULL));
789566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
79793850ffSBarry Smith   PetscFunctionReturn(0);
80793850ffSBarry Smith }
81793850ffSBarry Smith 
826d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y)
836d7c1e57SBarry Smith {
846d7c1e57SBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
856d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head;
866d7c1e57SBarry Smith   Vec               in,out;
8703049c21SJunchao Zhang   PetscScalar       scale;
8803049c21SJunchao Zhang   PetscInt          i;
896d7c1e57SBarry Smith 
906d7c1e57SBarry Smith   PetscFunctionBegin;
9128b400f6SJacob Faibussowitsch   PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
926d7c1e57SBarry Smith   in = x;
93e4fc5df0SBarry Smith   if (shell->right) {
94e4fc5df0SBarry Smith     if (!shell->rightwork) {
959566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(shell->right,&shell->rightwork));
96e4fc5df0SBarry Smith     }
979566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->rightwork,shell->right,in));
98e4fc5df0SBarry Smith     in   = shell->rightwork;
99e4fc5df0SBarry Smith   }
1006d7c1e57SBarry Smith   while (next->next) {
1016d7c1e57SBarry Smith     if (!next->work) { /* should reuse previous work if the same size */
1029566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(next->mat,NULL,&next->work));
1036d7c1e57SBarry Smith     }
1046d7c1e57SBarry Smith     out  = next->work;
1059566063dSJacob Faibussowitsch     PetscCall(MatMult(next->mat,in,out));
1066d7c1e57SBarry Smith     in   = out;
1076d7c1e57SBarry Smith     next = next->next;
1086d7c1e57SBarry Smith   }
1099566063dSJacob Faibussowitsch   PetscCall(MatMult(next->mat,in,y));
110e4fc5df0SBarry Smith   if (shell->left) {
1119566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y,shell->left,y));
112e4fc5df0SBarry Smith   }
11303049c21SJunchao Zhang   scale = shell->scale;
11403049c21SJunchao Zhang   if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
1159566063dSJacob Faibussowitsch   PetscCall(VecScale(y,scale));
1166d7c1e57SBarry Smith   PetscFunctionReturn(0);
1176d7c1e57SBarry Smith }
1186d7c1e57SBarry Smith 
1196d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y)
1206d7c1e57SBarry Smith {
1216d7c1e57SBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
1226d7c1e57SBarry Smith   Mat_CompositeLink tail   = shell->tail;
1236d7c1e57SBarry Smith   Vec               in,out;
12403049c21SJunchao Zhang   PetscScalar       scale;
12503049c21SJunchao Zhang   PetscInt          i;
1266d7c1e57SBarry Smith 
1276d7c1e57SBarry Smith   PetscFunctionBegin;
12828b400f6SJacob Faibussowitsch   PetscCheck(tail,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
1296d7c1e57SBarry Smith   in = x;
130e4fc5df0SBarry Smith   if (shell->left) {
131e4fc5df0SBarry Smith     if (!shell->leftwork) {
1329566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(shell->left,&shell->leftwork));
133e4fc5df0SBarry Smith     }
1349566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in));
135e4fc5df0SBarry Smith     in   = shell->leftwork;
136e4fc5df0SBarry Smith   }
1376d7c1e57SBarry Smith   while (tail->prev) {
1386d7c1e57SBarry Smith     if (!tail->prev->work) { /* should reuse previous work if the same size */
1399566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(tail->mat,NULL,&tail->prev->work));
1406d7c1e57SBarry Smith     }
1416d7c1e57SBarry Smith     out  = tail->prev->work;
1429566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(tail->mat,in,out));
1436d7c1e57SBarry Smith     in   = out;
1446d7c1e57SBarry Smith     tail = tail->prev;
1456d7c1e57SBarry Smith   }
1469566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(tail->mat,in,y));
147e4fc5df0SBarry Smith   if (shell->right) {
1489566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y,shell->right,y));
149e4fc5df0SBarry Smith   }
15003049c21SJunchao Zhang 
15103049c21SJunchao Zhang   scale = shell->scale;
15203049c21SJunchao Zhang   if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
1539566063dSJacob Faibussowitsch   PetscCall(VecScale(y,scale));
1546d7c1e57SBarry Smith   PetscFunctionReturn(0);
1556d7c1e57SBarry Smith }
1566d7c1e57SBarry Smith 
15703049c21SJunchao Zhang PetscErrorCode MatMult_Composite(Mat mat,Vec x,Vec y)
158793850ffSBarry Smith {
15903049c21SJunchao Zhang   Mat_Composite     *shell = (Mat_Composite*)mat->data;
16003049c21SJunchao Zhang   Mat_CompositeLink cur = shell->head;
16103049c21SJunchao Zhang   Vec               in,y2,xin;
16203049c21SJunchao Zhang   Mat               A,B;
16303049c21SJunchao Zhang   PetscInt          i,j,k,n,nuniq,lo,hi,mid,*gindices,*buf,*tmp,tot;
16403049c21SJunchao Zhang   const PetscScalar *vals;
16503049c21SJunchao Zhang   const PetscInt    *garray;
16603049c21SJunchao Zhang   IS                ix,iy;
16703049c21SJunchao Zhang   PetscBool         match;
168793850ffSBarry Smith 
169793850ffSBarry Smith   PetscFunctionBegin;
17028b400f6SJacob Faibussowitsch   PetscCheck(cur,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
171e4fc5df0SBarry Smith   in = x;
172e4fc5df0SBarry Smith   if (shell->right) {
173e4fc5df0SBarry Smith     if (!shell->rightwork) {
1749566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(shell->right,&shell->rightwork));
175793850ffSBarry Smith     }
1769566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->rightwork,shell->right,in));
177e4fc5df0SBarry Smith     in   = shell->rightwork;
178e4fc5df0SBarry Smith   }
17903049c21SJunchao Zhang 
18003049c21SJunchao Zhang   /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time
18103049c21SJunchao Zhang      we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and
18203049c21SJunchao Zhang      it is legal to merge Mvctx, because all component matrices have the same size.
18303049c21SJunchao Zhang    */
18403049c21SJunchao Zhang   if (shell->merge_mvctx && !shell->Mvctx) {
18503049c21SJunchao Zhang     /* Currently only implemented for MATMPIAIJ */
18603049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
1879566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat,MATMPIAIJ,&match));
18803049c21SJunchao Zhang       if (!match) {
18903049c21SJunchao Zhang         shell->merge_mvctx = PETSC_FALSE;
19003049c21SJunchao Zhang         goto skip_merge_mvctx;
191e4fc5df0SBarry Smith       }
192e4fc5df0SBarry Smith     }
19303049c21SJunchao Zhang 
19403049c21SJunchao Zhang     /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */
19503049c21SJunchao Zhang     tot = 0;
19603049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
1979566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,NULL));
1989566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,NULL,&n));
19903049c21SJunchao Zhang       tot += n;
20003049c21SJunchao Zhang     }
2019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(tot,&shell->location,tot,&shell->larray,shell->nmat,&shell->lvecs));
20203049c21SJunchao Zhang     shell->len = tot;
20303049c21SJunchao Zhang 
20403049c21SJunchao Zhang     /* Go through matrices second time to sort off-diag columns and remove dups */
2059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot,&gindices)); /* No Malloc2() since we will give one to petsc and free the other */
2069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tot,&buf));
20703049c21SJunchao Zhang     nuniq = 0; /* Number of unique nonzero columns */
20803049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
2099566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray));
2109566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,NULL,&n));
21103049c21SJunchao Zhang       /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */
21203049c21SJunchao Zhang       i = j = k = 0;
21303049c21SJunchao Zhang       while (i < n && j < nuniq) {
21403049c21SJunchao Zhang         if (garray[i] < gindices[j]) buf[k++] = garray[i++];
21503049c21SJunchao Zhang         else if (garray[i] > gindices[j]) buf[k++] = gindices[j++];
21603049c21SJunchao Zhang         else {buf[k++] = garray[i++]; j++;}
21703049c21SJunchao Zhang       }
21803049c21SJunchao Zhang       /* Copy leftover in garray[] or gindices[] */
21903049c21SJunchao Zhang       if (i < n) {
2209566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf+k,garray+i,n-i));
22103049c21SJunchao Zhang         nuniq = k + n-i;
22203049c21SJunchao Zhang       } else if (j < nuniq) {
2239566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(buf+k,gindices+j,nuniq-j));
22403049c21SJunchao Zhang         nuniq = k + nuniq-j;
22503049c21SJunchao Zhang       } else nuniq = k;
22603049c21SJunchao Zhang       /* Swap gindices and buf to merge garray of the next matrix */
22703049c21SJunchao Zhang       tmp      = gindices;
22803049c21SJunchao Zhang       gindices = buf;
22903049c21SJunchao Zhang       buf      = tmp;
23003049c21SJunchao Zhang     }
2319566063dSJacob Faibussowitsch     PetscCall(PetscFree(buf));
23203049c21SJunchao Zhang 
23303049c21SJunchao Zhang     /* Go through matrices third time to build a map from gindices[] to garray[] */
23403049c21SJunchao Zhang     tot = 0;
23503049c21SJunchao Zhang     for (cur=shell->head,j=0; cur; cur=cur->next,j++) { /* j-th matrix */
2369566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray));
2379566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,NULL,&n));
2389566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&shell->lvecs[j]));
23903049c21SJunchao Zhang       /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */
24003049c21SJunchao Zhang       lo   = 0;
24103049c21SJunchao Zhang       for (i=0; i<n; i++) {
24203049c21SJunchao Zhang         hi = nuniq;
24303049c21SJunchao Zhang         while (hi - lo > 1) {
24403049c21SJunchao Zhang           mid = lo + (hi - lo)/2;
24503049c21SJunchao Zhang           if (garray[i] < gindices[mid]) hi = mid;
24603049c21SJunchao Zhang           else lo = mid;
24703049c21SJunchao Zhang         }
24803049c21SJunchao Zhang         shell->location[tot+i] = lo; /* gindices[lo] = garray[i] */
24903049c21SJunchao Zhang         lo++; /* Since garray[i+1] > garray[i], we can safely advance lo */
25003049c21SJunchao Zhang       }
25103049c21SJunchao Zhang       tot += n;
25203049c21SJunchao Zhang     }
25303049c21SJunchao Zhang 
25403049c21SJunchao Zhang     /* Build merged Mvctx */
2559566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nuniq,gindices,PETSC_OWN_POINTER,&ix));
2569566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,nuniq,0,1,&iy));
2579566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&xin));
2589566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,nuniq,&shell->gvec));
2599566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(xin,ix,shell->gvec,iy,&shell->Mvctx));
2609566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xin));
2619566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ix));
2629566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iy));
26303049c21SJunchao Zhang   }
26403049c21SJunchao Zhang 
26503049c21SJunchao Zhang skip_merge_mvctx:
2669566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0));
2679566063dSJacob Faibussowitsch   if (!shell->leftwork2) PetscCall(VecDuplicate(y,&shell->leftwork2));
26803049c21SJunchao Zhang   y2 = shell->leftwork2;
26903049c21SJunchao Zhang 
27003049c21SJunchao Zhang   if (shell->Mvctx) { /* Have a merged Mvctx */
27103049c21SJunchao 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
27203049c21SJunchao Zhang        in y = s(Ax1 + Bx2) or y = sAx1 + sBx2. The former incurs less FLOPS than the latter, but the latter provides an oppertunity to
27303049c21SJunchao Zhang        overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach.
27403049c21SJunchao Zhang      */
2759566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD));
2769566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD));
27703049c21SJunchao Zhang 
2789566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(shell->gvec,&vals));
27903049c21SJunchao Zhang     for (i=0; i<shell->len; i++) shell->larray[i] = vals[shell->location[i]];
2809566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(shell->gvec,&vals));
28103049c21SJunchao Zhang 
28203049c21SJunchao Zhang     for (cur=shell->head,tot=i=0; cur; cur=cur->next,i++) { /* i-th matrix */
2839566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,&A,&B,NULL));
2849566063dSJacob Faibussowitsch       PetscCall((*A->ops->mult)(A,in,y2));
2859566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(B,NULL,&n));
2869566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(shell->lvecs[i],&shell->larray[tot]));
2879566063dSJacob Faibussowitsch       PetscCall((*B->ops->multadd)(B,shell->lvecs[i],y2,y2));
2889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(shell->lvecs[i]));
2899566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,(shell->scalings ? shell->scalings[i] : 1.0),y2));
29003049c21SJunchao Zhang       tot += n;
29103049c21SJunchao Zhang     }
29203049c21SJunchao Zhang   } else {
29303049c21SJunchao Zhang     if (shell->scalings) {
29403049c21SJunchao Zhang       for (cur=shell->head,i=0; cur; cur=cur->next,i++) {
2959566063dSJacob Faibussowitsch         PetscCall(MatMult(cur->mat,in,y2));
2969566063dSJacob Faibussowitsch         PetscCall(VecAXPY(y,shell->scalings[i],y2));
29703049c21SJunchao Zhang       }
29803049c21SJunchao Zhang     } else {
2999566063dSJacob Faibussowitsch       for (cur=shell->head; cur; cur=cur->next) PetscCall(MatMultAdd(cur->mat,in,y,y));
30003049c21SJunchao Zhang     }
30103049c21SJunchao Zhang   }
30203049c21SJunchao Zhang 
3039566063dSJacob Faibussowitsch   if (shell->left) PetscCall(VecPointwiseMult(y,shell->left,y));
3049566063dSJacob Faibussowitsch   PetscCall(VecScale(y,shell->scale));
305793850ffSBarry Smith   PetscFunctionReturn(0);
306793850ffSBarry Smith }
307793850ffSBarry Smith 
308793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y)
309793850ffSBarry Smith {
310793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
311793850ffSBarry Smith   Mat_CompositeLink next   = shell->head;
31203049c21SJunchao Zhang   Vec               in,y2 = NULL;
31303049c21SJunchao Zhang   PetscInt          i;
314793850ffSBarry Smith 
315793850ffSBarry Smith   PetscFunctionBegin;
31628b400f6SJacob Faibussowitsch   PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
317e4fc5df0SBarry Smith   in = x;
318e4fc5df0SBarry Smith   if (shell->left) {
319e4fc5df0SBarry Smith     if (!shell->leftwork) {
3209566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(shell->left,&shell->leftwork));
321793850ffSBarry Smith     }
3229566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in));
323e4fc5df0SBarry Smith     in   = shell->leftwork;
324e4fc5df0SBarry Smith   }
32503049c21SJunchao Zhang 
3269566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(next->mat,in,y));
32703049c21SJunchao Zhang   if (shell->scalings) {
3289566063dSJacob Faibussowitsch     PetscCall(VecScale(y,shell->scalings[0]));
3299566063dSJacob Faibussowitsch     if (!shell->rightwork2) PetscCall(VecDuplicate(y,&shell->rightwork2));
33003049c21SJunchao Zhang     y2 = shell->rightwork2;
33103049c21SJunchao Zhang   }
33203049c21SJunchao Zhang   i = 1;
333e4fc5df0SBarry Smith   while ((next = next->next)) {
3349566063dSJacob Faibussowitsch     if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat,in,y,y));
33503049c21SJunchao Zhang     else {
3369566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(next->mat,in,y2));
3379566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,shell->scalings[i++],y2));
33803049c21SJunchao Zhang     }
339e4fc5df0SBarry Smith   }
340e4fc5df0SBarry Smith   if (shell->right) {
3419566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y,shell->right,y));
342e4fc5df0SBarry Smith   }
3439566063dSJacob Faibussowitsch   PetscCall(VecScale(y,shell->scale));
344793850ffSBarry Smith   PetscFunctionReturn(0);
345793850ffSBarry Smith }
346793850ffSBarry Smith 
3477bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z)
3487bf3a71bSJakub Kruzik {
3497bf3a71bSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)A->data;
3507bf3a71bSJakub Kruzik 
3517bf3a71bSJakub Kruzik   PetscFunctionBegin;
3527bf3a71bSJakub Kruzik   if (y != z) {
3539566063dSJacob Faibussowitsch     PetscCall(MatMult(A,x,z));
3549566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,1.0,y));
3557bf3a71bSJakub Kruzik   } else {
3567bf3a71bSJakub Kruzik     if (!shell->leftwork) {
3579566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(z,&shell->leftwork));
3587bf3a71bSJakub Kruzik     }
3599566063dSJacob Faibussowitsch     PetscCall(MatMult(A,x,shell->leftwork));
3609566063dSJacob Faibussowitsch     PetscCall(VecCopy(y,z));
3619566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,1.0,shell->leftwork));
3627bf3a71bSJakub Kruzik   }
3637bf3a71bSJakub Kruzik   PetscFunctionReturn(0);
3647bf3a71bSJakub Kruzik }
3657bf3a71bSJakub Kruzik 
3667bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z)
3677bf3a71bSJakub Kruzik {
3687bf3a71bSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)A->data;
3697bf3a71bSJakub Kruzik 
3707bf3a71bSJakub Kruzik   PetscFunctionBegin;
3717bf3a71bSJakub Kruzik   if (y != z) {
3729566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,x,z));
3739566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,1.0,y));
3747bf3a71bSJakub Kruzik   } else {
3757bf3a71bSJakub Kruzik     if (!shell->rightwork) {
3769566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(z,&shell->rightwork));
3777bf3a71bSJakub Kruzik     }
3789566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,x,shell->rightwork));
3799566063dSJacob Faibussowitsch     PetscCall(VecCopy(y,z));
3809566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,1.0,shell->rightwork));
3817bf3a71bSJakub Kruzik   }
3827bf3a71bSJakub Kruzik   PetscFunctionReturn(0);
3837bf3a71bSJakub Kruzik }
3847bf3a71bSJakub Kruzik 
385793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v)
386793850ffSBarry Smith {
387793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
388793850ffSBarry Smith   Mat_CompositeLink next   = shell->head;
38903049c21SJunchao Zhang   PetscInt          i;
390793850ffSBarry Smith 
391793850ffSBarry Smith   PetscFunctionBegin;
39228b400f6SJacob Faibussowitsch   PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
39308401ef6SPierre Jolivet   PetscCheck(!shell->right && !shell->left,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling");
394e4fc5df0SBarry Smith 
3959566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(next->mat,v));
3969566063dSJacob Faibussowitsch   if (shell->scalings) PetscCall(VecScale(v,shell->scalings[0]));
39703049c21SJunchao Zhang 
398793850ffSBarry Smith   if (next->next && !shell->work) {
3999566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v,&shell->work));
400793850ffSBarry Smith   }
40103049c21SJunchao Zhang   i = 1;
402793850ffSBarry Smith   while ((next = next->next)) {
4039566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(next->mat,shell->work));
4049566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v,(shell->scalings ? shell->scalings[i++] : 1.0),shell->work));
405793850ffSBarry Smith   }
4069566063dSJacob Faibussowitsch   PetscCall(VecScale(v,shell->scale));
407793850ffSBarry Smith   PetscFunctionReturn(0);
408793850ffSBarry Smith }
409793850ffSBarry Smith 
410793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)
411793850ffSBarry Smith {
4124b2558d6SJakub Kruzik   Mat_Composite  *shell = (Mat_Composite*)Y->data;
413b52f573bSBarry Smith 
414793850ffSBarry Smith   PetscFunctionBegin;
4154b2558d6SJakub Kruzik   if (shell->merge) {
4169566063dSJacob Faibussowitsch     PetscCall(MatCompositeMerge(Y));
417b52f573bSBarry Smith   }
418793850ffSBarry Smith   PetscFunctionReturn(0);
419793850ffSBarry Smith }
420793850ffSBarry Smith 
421e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha)
422e4fc5df0SBarry Smith {
423e4fc5df0SBarry Smith   Mat_Composite *a = (Mat_Composite*)inA->data;
424e4fc5df0SBarry Smith 
425e4fc5df0SBarry Smith   PetscFunctionBegin;
426321429dbSBarry Smith   a->scale *= alpha;
427e4fc5df0SBarry Smith   PetscFunctionReturn(0);
428e4fc5df0SBarry Smith }
429e4fc5df0SBarry Smith 
430e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right)
431e4fc5df0SBarry Smith {
432e4fc5df0SBarry Smith   Mat_Composite  *a = (Mat_Composite*)inA->data;
433e4fc5df0SBarry Smith 
434e4fc5df0SBarry Smith   PetscFunctionBegin;
435e4fc5df0SBarry Smith   if (left) {
436321429dbSBarry Smith     if (!a->left) {
4379566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left,&a->left));
4389566063dSJacob Faibussowitsch       PetscCall(VecCopy(left,a->left));
439321429dbSBarry Smith     } else {
4409566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left,left,a->left));
441321429dbSBarry Smith     }
442e4fc5df0SBarry Smith   }
443e4fc5df0SBarry Smith   if (right) {
444321429dbSBarry Smith     if (!a->right) {
4459566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right,&a->right));
4469566063dSJacob Faibussowitsch       PetscCall(VecCopy(right,a->right));
447321429dbSBarry Smith     } else {
4489566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right,right,a->right));
449321429dbSBarry Smith     }
450e4fc5df0SBarry Smith   }
451e4fc5df0SBarry Smith   PetscFunctionReturn(0);
452e4fc5df0SBarry Smith }
453793850ffSBarry Smith 
4544b2558d6SJakub Kruzik PetscErrorCode MatSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,Mat A)
4554b2558d6SJakub Kruzik {
4564b2558d6SJakub Kruzik   Mat_Composite  *a = (Mat_Composite*)A->data;
4574b2558d6SJakub Kruzik 
4584b2558d6SJakub Kruzik   PetscFunctionBegin;
459d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"MATCOMPOSITE options");
4609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge","Merge at MatAssemblyEnd","MatCompositeMerge",a->merge,&a->merge,NULL));
4619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-mat_composite_merge_type","Set composite merge direction","MatCompositeSetMergeType",MatCompositeMergeTypes,(PetscEnum)a->mergetype,(PetscEnum*)&a->mergetype,NULL));
4629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx","Merge MatMult() vecscat contexts","MatCreateComposite",a->merge_mvctx,&a->merge_mvctx,NULL));
463d0609cedSBarry Smith   PetscOptionsHeadEnd();
4644b2558d6SJakub Kruzik   PetscFunctionReturn(0);
4654b2558d6SJakub Kruzik }
4664b2558d6SJakub Kruzik 
4672c0821f3SBarry Smith /*@
4688bd739bdSJakub Kruzik    MatCreateComposite - Creates a matrix as the sum or product of one or more matrices
469793850ffSBarry Smith 
470d083f849SBarry Smith   Collective
471793850ffSBarry Smith 
472793850ffSBarry Smith    Input Parameters:
473793850ffSBarry Smith +  comm - MPI communicator
474793850ffSBarry Smith .  nmat - number of matrices to put in
475793850ffSBarry Smith -  mats - the matrices
476793850ffSBarry Smith 
477793850ffSBarry Smith    Output Parameter:
478db36af27SMatthew G. Knepley .  mat - the matrix
479793850ffSBarry Smith 
4804b2558d6SJakub Kruzik    Options Database Keys:
4814b2558d6SJakub Kruzik +  -mat_composite_merge         - merge in MatAssemblyEnd()
48203049c21SJunchao Zhang .  -mat_composite_merge_mvctx   - merge Mvctx of component matrices to optimize communication in MatMult() for ADDITIVE matrices
483b28d0daaSJakub Kruzik -  -mat_composite_merge_type    - set merge direction
4844b2558d6SJakub Kruzik 
485793850ffSBarry Smith    Level: advanced
486793850ffSBarry Smith 
487793850ffSBarry Smith    Notes:
488793850ffSBarry Smith      Alternative construction
489793850ffSBarry Smith $       MatCreate(comm,&mat);
490793850ffSBarry Smith $       MatSetSizes(mat,m,n,M,N);
491793850ffSBarry Smith $       MatSetType(mat,MATCOMPOSITE);
492793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[0]);
493793850ffSBarry Smith $       ....
494793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[nmat-1]);
495b52f573bSBarry Smith $       MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
496b52f573bSBarry Smith $       MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
497793850ffSBarry Smith 
4986d7c1e57SBarry Smith      For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0]
4996d7c1e57SBarry Smith 
500db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCompositeGetMat()`, `MatCompositeMerge()`, `MatCompositeSetType()`, `MATCOMPOSITE`
501793850ffSBarry Smith 
502793850ffSBarry Smith @*/
5037087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat)
504793850ffSBarry Smith {
505793850ffSBarry Smith   PetscInt       m,n,M,N,i;
506793850ffSBarry Smith 
507793850ffSBarry Smith   PetscFunctionBegin;
50808401ef6SPierre Jolivet   PetscCheck(nmat >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix");
509064a246eSJacob Faibussowitsch   PetscValidPointer(mat,4);
510793850ffSBarry Smith 
5119566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[0],PETSC_IGNORE,&n));
5129566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE));
5139566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[0],PETSC_IGNORE,&N));
5149566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mats[nmat-1],&M,PETSC_IGNORE));
5159566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,mat));
5169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mat,m,n,M,N));
5179566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat,MATCOMPOSITE));
518793850ffSBarry Smith   for (i=0; i<nmat; i++) {
5199566063dSJacob Faibussowitsch     PetscCall(MatCompositeAddMat(*mat,mats[i]));
520793850ffSBarry Smith   }
5219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY));
5229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY));
523793850ffSBarry Smith   PetscFunctionReturn(0);
524793850ffSBarry Smith }
525793850ffSBarry Smith 
526d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat)
527d7f81bf2SJakub Kruzik {
528d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
529d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink,next = shell->head;
530d7f81bf2SJakub Kruzik 
531d7f81bf2SJakub Kruzik   PetscFunctionBegin;
5329566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(mat,&ilink));
533f4259b30SLisandro Dalcin   ilink->next = NULL;
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)smat));
535d7f81bf2SJakub Kruzik   ilink->mat  = smat;
536d7f81bf2SJakub Kruzik 
537d7f81bf2SJakub Kruzik   if (!next) shell->head = ilink;
538d7f81bf2SJakub Kruzik   else {
539d7f81bf2SJakub Kruzik     while (next->next) {
540d7f81bf2SJakub Kruzik       next = next->next;
541d7f81bf2SJakub Kruzik     }
542d7f81bf2SJakub Kruzik     next->next  = ilink;
543d7f81bf2SJakub Kruzik     ilink->prev = next;
544d7f81bf2SJakub Kruzik   }
545d7f81bf2SJakub Kruzik   shell->tail =  ilink;
546d7f81bf2SJakub Kruzik   shell->nmat += 1;
54703049c21SJunchao Zhang 
54803049c21SJunchao Zhang   /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */
54903049c21SJunchao Zhang   if (shell->scalings) {
5509566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(sizeof(PetscScalar)*shell->nmat,&shell->scalings));
55103049c21SJunchao Zhang     shell->scalings[shell->nmat-1] = 1.0;
55203049c21SJunchao Zhang   }
553d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
554d7f81bf2SJakub Kruzik }
555d7f81bf2SJakub Kruzik 
556793850ffSBarry Smith /*@
5578bd739bdSJakub Kruzik     MatCompositeAddMat - Add another matrix to a composite matrix.
558793850ffSBarry Smith 
559793850ffSBarry Smith    Collective on Mat
560793850ffSBarry Smith 
561793850ffSBarry Smith     Input Parameters:
562793850ffSBarry Smith +   mat - the composite matrix
563793850ffSBarry Smith -   smat - the partial matrix
564793850ffSBarry Smith 
565793850ffSBarry Smith    Level: advanced
566793850ffSBarry Smith 
567db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
568793850ffSBarry Smith @*/
5697087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat)
570793850ffSBarry Smith {
571793850ffSBarry Smith   PetscFunctionBegin;
5720700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5730700a824SBarry Smith   PetscValidHeaderSpecific(smat,MAT_CLASSID,2);
574cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));
575d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
576d7f81bf2SJakub Kruzik }
577793850ffSBarry Smith 
578d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type)
579d7f81bf2SJakub Kruzik {
580d7f81bf2SJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
581d7f81bf2SJakub Kruzik 
582d7f81bf2SJakub Kruzik   PetscFunctionBegin;
583ced1ac25SJakub Kruzik   b->type = type;
584d7f81bf2SJakub Kruzik   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
585f4259b30SLisandro Dalcin     mat->ops->getdiagonal   = NULL;
586d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite_Multiplicative;
587d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative;
58803049c21SJunchao Zhang     b->merge_mvctx          = PETSC_FALSE;
589d7f81bf2SJakub Kruzik   } else {
590d7f81bf2SJakub Kruzik     mat->ops->getdiagonal   = MatGetDiagonal_Composite;
591d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite;
592d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite;
593793850ffSBarry Smith   }
5946d7c1e57SBarry Smith   PetscFunctionReturn(0);
5956d7c1e57SBarry Smith }
5966d7c1e57SBarry Smith 
5972c0821f3SBarry Smith /*@
5988bd739bdSJakub Kruzik    MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product.
5996d7c1e57SBarry Smith 
600b2b245abSJakub Kruzik    Logically Collective on Mat
6016d7c1e57SBarry Smith 
6026d7c1e57SBarry Smith    Input Parameters:
6036d7c1e57SBarry Smith .  mat - the composite matrix
6046d7c1e57SBarry Smith 
6056d7c1e57SBarry Smith    Level: advanced
6066d7c1e57SBarry Smith 
607db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE`
6086d7c1e57SBarry Smith 
6096d7c1e57SBarry Smith @*/
6107087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type)
6116d7c1e57SBarry Smith {
6126d7c1e57SBarry Smith   PetscFunctionBegin;
613d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
614b2b245abSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat,type,2);
615cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));
616793850ffSBarry Smith   PetscFunctionReturn(0);
617793850ffSBarry Smith }
618793850ffSBarry Smith 
6196dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type)
6206dbc55e5SJakub Kruzik {
6216dbc55e5SJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
6226dbc55e5SJakub Kruzik 
6236dbc55e5SJakub Kruzik   PetscFunctionBegin;
6246dbc55e5SJakub Kruzik   *type = b->type;
6256dbc55e5SJakub Kruzik   PetscFunctionReturn(0);
6266dbc55e5SJakub Kruzik }
6276dbc55e5SJakub Kruzik 
6286dbc55e5SJakub Kruzik /*@
6296dbc55e5SJakub Kruzik    MatCompositeGetType - Returns type of composite.
6306dbc55e5SJakub Kruzik 
6316dbc55e5SJakub Kruzik    Not Collective
6326dbc55e5SJakub Kruzik 
6336dbc55e5SJakub Kruzik    Input Parameter:
6346dbc55e5SJakub Kruzik .  mat - the composite matrix
6356dbc55e5SJakub Kruzik 
6366dbc55e5SJakub Kruzik    Output Parameter:
6376dbc55e5SJakub Kruzik .  type - type of composite
6386dbc55e5SJakub Kruzik 
6396dbc55e5SJakub Kruzik    Level: advanced
6406dbc55e5SJakub Kruzik 
641db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE`
6426dbc55e5SJakub Kruzik 
6436dbc55e5SJakub Kruzik @*/
6446dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type)
6456dbc55e5SJakub Kruzik {
6466dbc55e5SJakub Kruzik   PetscFunctionBegin;
6476dbc55e5SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6486dbc55e5SJakub Kruzik   PetscValidPointer(type,2);
649cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type));
6506dbc55e5SJakub Kruzik   PetscFunctionReturn(0);
6516dbc55e5SJakub Kruzik }
6526dbc55e5SJakub Kruzik 
6533b35acafSJakub Kruzik static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat,MatStructure str)
6543b35acafSJakub Kruzik {
6553b35acafSJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
6563b35acafSJakub Kruzik 
6573b35acafSJakub Kruzik   PetscFunctionBegin;
6583b35acafSJakub Kruzik   b->structure = str;
6593b35acafSJakub Kruzik   PetscFunctionReturn(0);
6603b35acafSJakub Kruzik }
6613b35acafSJakub Kruzik 
6623b35acafSJakub Kruzik /*@
6633b35acafSJakub Kruzik    MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix.
6643b35acafSJakub Kruzik 
6653b35acafSJakub Kruzik    Not Collective
6663b35acafSJakub Kruzik 
6673b35acafSJakub Kruzik    Input Parameters:
6683b35acafSJakub Kruzik +  mat - the composite matrix
6693b35acafSJakub Kruzik -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN (default) or SUBSET_NONZERO_PATTERN
6703b35acafSJakub Kruzik 
6713b35acafSJakub Kruzik    Level: advanced
6723b35acafSJakub Kruzik 
6733b35acafSJakub Kruzik    Notes:
6743b35acafSJakub Kruzik     Information about the matrices structure is used in MatCompositeMerge() for additive composite matrix.
6753b35acafSJakub Kruzik 
676db781477SPatrick Sanan .seealso: `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE`
6773b35acafSJakub Kruzik 
6783b35acafSJakub Kruzik @*/
6793b35acafSJakub Kruzik PetscErrorCode MatCompositeSetMatStructure(Mat mat,MatStructure str)
6803b35acafSJakub Kruzik {
6813b35acafSJakub Kruzik   PetscFunctionBegin;
6823b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
683cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeSetMatStructure_C",(Mat,MatStructure),(mat,str));
6843b35acafSJakub Kruzik   PetscFunctionReturn(0);
6853b35acafSJakub Kruzik }
6863b35acafSJakub Kruzik 
6873b35acafSJakub Kruzik static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat,MatStructure *str)
6883b35acafSJakub Kruzik {
6893b35acafSJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
6903b35acafSJakub Kruzik 
6913b35acafSJakub Kruzik   PetscFunctionBegin;
6923b35acafSJakub Kruzik   *str = b->structure;
6933b35acafSJakub Kruzik   PetscFunctionReturn(0);
6943b35acafSJakub Kruzik }
6953b35acafSJakub Kruzik 
6963b35acafSJakub Kruzik /*@
6973b35acafSJakub Kruzik    MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix.
6983b35acafSJakub Kruzik 
6993b35acafSJakub Kruzik    Not Collective
7003b35acafSJakub Kruzik 
7013b35acafSJakub Kruzik    Input Parameter:
7023b35acafSJakub Kruzik .  mat - the composite matrix
7033b35acafSJakub Kruzik 
7043b35acafSJakub Kruzik    Output Parameter:
7053b35acafSJakub Kruzik .  str - structure of the matrices
7063b35acafSJakub Kruzik 
7073b35acafSJakub Kruzik    Level: advanced
7083b35acafSJakub Kruzik 
709db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE`
7103b35acafSJakub Kruzik 
7113b35acafSJakub Kruzik @*/
7123b35acafSJakub Kruzik PetscErrorCode MatCompositeGetMatStructure(Mat mat,MatStructure *str)
7133b35acafSJakub Kruzik {
7143b35acafSJakub Kruzik   PetscFunctionBegin;
7153b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7163b35acafSJakub Kruzik   PetscValidPointer(str,2);
717cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeGetMatStructure_C",(Mat,MatStructure*),(mat,str));
7183b35acafSJakub Kruzik   PetscFunctionReturn(0);
7193b35acafSJakub Kruzik }
7203b35acafSJakub Kruzik 
7218c8613bfSJakub Kruzik static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type)
7228c8613bfSJakub Kruzik {
7238c8613bfSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
7248c8613bfSJakub Kruzik 
7258c8613bfSJakub Kruzik   PetscFunctionBegin;
7268c8613bfSJakub Kruzik   shell->mergetype = type;
7278c8613bfSJakub Kruzik   PetscFunctionReturn(0);
7288c8613bfSJakub Kruzik }
7298c8613bfSJakub Kruzik 
7308c8613bfSJakub Kruzik /*@
7318c8613bfSJakub Kruzik    MatCompositeSetMergeType - Sets order of MatCompositeMerge().
7328c8613bfSJakub Kruzik 
7338c8613bfSJakub Kruzik    Logically Collective on Mat
7348c8613bfSJakub Kruzik 
735d8d19677SJose E. Roman    Input Parameters:
7368c8613bfSJakub Kruzik +  mat - the composite matrix
7378c8613bfSJakub Kruzik -  type - MAT_COMPOSITE_MERGE RIGHT (default) to start merge from right with the first added matrix (mat[0]),
7388c8613bfSJakub Kruzik           MAT_COMPOSITE_MERGE_LEFT to start merge from left with the last added matrix (mat[nmat-1])
7398c8613bfSJakub Kruzik 
7408c8613bfSJakub Kruzik    Level: advanced
7418c8613bfSJakub Kruzik 
7428c8613bfSJakub Kruzik    Notes:
7438c8613bfSJakub Kruzik     The resulting matrix is the same regardles of the MergeType. Only the order of operation is changed.
7448c8613bfSJakub Kruzik     If set to MAT_COMPOSITE_MERGE_RIGHT the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0])))
7458c8613bfSJakub Kruzik     otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0].
7468c8613bfSJakub Kruzik 
747db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE`
7488c8613bfSJakub Kruzik 
7498c8613bfSJakub Kruzik @*/
7508c8613bfSJakub Kruzik PetscErrorCode MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type)
7518c8613bfSJakub Kruzik {
7528c8613bfSJakub Kruzik   PetscFunctionBegin;
7538c8613bfSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7548c8613bfSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat,type,2);
755cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeSetMergeType_C",(Mat,MatCompositeMergeType),(mat,type));
7568c8613bfSJakub Kruzik   PetscFunctionReturn(0);
7578c8613bfSJakub Kruzik }
7588c8613bfSJakub Kruzik 
759d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat)
760b52f573bSBarry Smith {
761b52f573bSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)mat->data;
7626d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head, prev = shell->tail;
7636d7c1e57SBarry Smith   Mat               tmat,newmat;
7641ba60692SJed Brown   Vec               left,right;
7651ba60692SJed Brown   PetscScalar       scale;
76603049c21SJunchao Zhang   PetscInt          i;
767b52f573bSBarry Smith 
768b52f573bSBarry Smith   PetscFunctionBegin;
76928b400f6SJacob Faibussowitsch   PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
77003049c21SJunchao Zhang   scale = shell->scale;
7716d7c1e57SBarry Smith   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
7728c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
77303049c21SJunchao Zhang       i = 0;
7749566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat));
7759566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i++]));
776b52f573bSBarry Smith       while ((next = next->next)) {
7779566063dSJacob Faibussowitsch         PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i++] : 1.0),next->mat,shell->structure));
778b52f573bSBarry Smith       }
7796d7c1e57SBarry Smith     } else {
78003049c21SJunchao Zhang       i = shell->nmat-1;
7819566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat));
7829566063dSJacob Faibussowitsch       if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i--]));
7838c8613bfSJakub Kruzik       while ((prev = prev->prev)) {
7849566063dSJacob Faibussowitsch         PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i--] : 1.0),prev->mat,shell->structure));
7858c8613bfSJakub Kruzik       }
7868c8613bfSJakub Kruzik     }
7878c8613bfSJakub Kruzik   } else {
7888c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
7899566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat));
790b6cfcaf5SJakub Kruzik       while ((next = next->next)) {
7919566063dSJacob Faibussowitsch         PetscCall(MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat));
7929566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
7936d7c1e57SBarry Smith         tmat = newmat;
7946d7c1e57SBarry Smith       }
79504d534ceSJakub Kruzik     } else {
7969566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat));
79704d534ceSJakub Kruzik       while ((prev = prev->prev)) {
7989566063dSJacob Faibussowitsch         PetscCall(MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat));
7999566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&tmat));
80004d534ceSJakub Kruzik         tmat = newmat;
80104d534ceSJakub Kruzik       }
80204d534ceSJakub Kruzik     }
80303049c21SJunchao Zhang     if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
8046d7c1e57SBarry Smith   }
8051ba60692SJed Brown 
8069566063dSJacob Faibussowitsch   if ((left = shell->left)) PetscCall(PetscObjectReference((PetscObject)left));
8079566063dSJacob Faibussowitsch   if ((right = shell->right)) PetscCall(PetscObjectReference((PetscObject)right));
8081ba60692SJed Brown 
8099566063dSJacob Faibussowitsch   PetscCall(MatHeaderReplace(mat,&tmat));
8101ba60692SJed Brown 
8119566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(mat,left,right));
8129566063dSJacob Faibussowitsch   PetscCall(MatScale(mat,scale));
8139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
8149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
815b52f573bSBarry Smith   PetscFunctionReturn(0);
816b52f573bSBarry Smith }
8176a7ede75SJakub Kruzik 
8186a7ede75SJakub Kruzik /*@
819d7f81bf2SJakub Kruzik    MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
8208bd739bdSJakub Kruzik      by summing or computing the product of all the matrices inside the composite matrix.
821d7f81bf2SJakub Kruzik 
822b2b245abSJakub Kruzik   Collective
823d7f81bf2SJakub Kruzik 
824f899ff85SJose E. Roman    Input Parameter:
825d7f81bf2SJakub Kruzik .  mat - the composite matrix
826d7f81bf2SJakub Kruzik 
8274b2558d6SJakub Kruzik    Options Database Keys:
828b28d0daaSJakub Kruzik +  -mat_composite_merge - merge in MatAssemblyEnd()
829b28d0daaSJakub Kruzik -  -mat_composite_merge_type - set merge direction
830d7f81bf2SJakub Kruzik 
831d7f81bf2SJakub Kruzik    Level: advanced
832d7f81bf2SJakub Kruzik 
833d7f81bf2SJakub Kruzik    Notes:
834d7f81bf2SJakub Kruzik       The MatType of the resulting matrix will be the same as the MatType of the FIRST
835d7f81bf2SJakub Kruzik     matrix in the composite matrix.
836d7f81bf2SJakub Kruzik 
837db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE`
838d7f81bf2SJakub Kruzik 
839d7f81bf2SJakub Kruzik @*/
840d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat)
841d7f81bf2SJakub Kruzik {
842d7f81bf2SJakub Kruzik   PetscFunctionBegin;
843d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
844cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));
845d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
846d7f81bf2SJakub Kruzik }
847d7f81bf2SJakub Kruzik 
8486d0add67SJakub Kruzik static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat,PetscInt *nmat)
849d7f81bf2SJakub Kruzik {
850d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
851d7f81bf2SJakub Kruzik 
852d7f81bf2SJakub Kruzik   PetscFunctionBegin;
853d7f81bf2SJakub Kruzik   *nmat = shell->nmat;
854d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
855d7f81bf2SJakub Kruzik }
856d7f81bf2SJakub Kruzik 
857d7f81bf2SJakub Kruzik /*@
8586d0add67SJakub Kruzik    MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix.
8596a7ede75SJakub Kruzik 
8606a7ede75SJakub Kruzik    Not Collective
8616a7ede75SJakub Kruzik 
8626a7ede75SJakub Kruzik    Input Parameter:
863d7f81bf2SJakub Kruzik .  mat - the composite matrix
8646a7ede75SJakub Kruzik 
8656a7ede75SJakub Kruzik    Output Parameter:
8666d0add67SJakub Kruzik .  nmat - number of matrices in the composite matrix
8676a7ede75SJakub Kruzik 
8688b5c3584SJakub Kruzik    Level: advanced
8696a7ede75SJakub Kruzik 
870db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE`
8716a7ede75SJakub Kruzik 
8726a7ede75SJakub Kruzik @*/
8736d0add67SJakub Kruzik PetscErrorCode MatCompositeGetNumberMat(Mat mat,PetscInt *nmat)
8746a7ede75SJakub Kruzik {
8756a7ede75SJakub Kruzik   PetscFunctionBegin;
876d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
877dadcf809SJacob Faibussowitsch   PetscValidIntPointer(nmat,2);
878cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeGetNumberMat_C",(Mat,PetscInt*),(mat,nmat));
879d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
880d7f81bf2SJakub Kruzik }
881d7f81bf2SJakub Kruzik 
882d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai)
883d7f81bf2SJakub Kruzik {
884d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
885d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink;
886d7f81bf2SJakub Kruzik   PetscInt          k;
887d7f81bf2SJakub Kruzik 
888d7f81bf2SJakub Kruzik   PetscFunctionBegin;
88908401ef6SPierre Jolivet   PetscCheck(i < shell->nmat,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT,i,shell->nmat);
890d7f81bf2SJakub Kruzik   ilink = shell->head;
891d7f81bf2SJakub Kruzik   for (k=0; k<i; k++) {
892d7f81bf2SJakub Kruzik     ilink = ilink->next;
893d7f81bf2SJakub Kruzik   }
894d7f81bf2SJakub Kruzik   *Ai = ilink->mat;
8956a7ede75SJakub Kruzik   PetscFunctionReturn(0);
8966a7ede75SJakub Kruzik }
8976a7ede75SJakub Kruzik 
8988b5c3584SJakub Kruzik /*@
8998bd739bdSJakub Kruzik    MatCompositeGetMat - Returns the ith matrix from the composite matrix.
9008b5c3584SJakub Kruzik 
9018b5c3584SJakub Kruzik    Logically Collective on Mat
9028b5c3584SJakub Kruzik 
903d8d19677SJose E. Roman    Input Parameters:
904d7f81bf2SJakub Kruzik +  mat - the composite matrix
9058b5c3584SJakub Kruzik -  i - the number of requested matrix
9068b5c3584SJakub Kruzik 
9078b5c3584SJakub Kruzik    Output Parameter:
9088b5c3584SJakub Kruzik .  Ai - ith matrix in composite
9098b5c3584SJakub Kruzik 
9108b5c3584SJakub Kruzik    Level: advanced
9118b5c3584SJakub Kruzik 
912db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE`
9138b5c3584SJakub Kruzik 
9148b5c3584SJakub Kruzik @*/
915d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai)
9168b5c3584SJakub Kruzik {
9178b5c3584SJakub Kruzik   PetscFunctionBegin;
918d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
919d7f81bf2SJakub Kruzik   PetscValidLogicalCollectiveInt(mat,i,2);
9208b5c3584SJakub Kruzik   PetscValidPointer(Ai,3);
921cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));
9228b5c3584SJakub Kruzik   PetscFunctionReturn(0);
9238b5c3584SJakub Kruzik }
9248b5c3584SJakub Kruzik 
92503049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings_Composite(Mat mat,const PetscScalar *scalings)
92603049c21SJunchao Zhang {
92703049c21SJunchao Zhang   Mat_Composite  *shell = (Mat_Composite*)mat->data;
92803049c21SJunchao Zhang   PetscInt       nmat;
92903049c21SJunchao Zhang 
93003049c21SJunchao Zhang   PetscFunctionBegin;
9319566063dSJacob Faibussowitsch   PetscCall(MatCompositeGetNumberMat(mat,&nmat));
9329566063dSJacob Faibussowitsch   if (!shell->scalings) PetscCall(PetscMalloc1(nmat,&shell->scalings));
9339566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(shell->scalings,scalings,nmat));
93403049c21SJunchao Zhang   PetscFunctionReturn(0);
93503049c21SJunchao Zhang }
93603049c21SJunchao Zhang 
93703049c21SJunchao Zhang /*@
93803049c21SJunchao Zhang    MatCompositeSetScalings - Sets separate scaling factors for component matrices.
93903049c21SJunchao Zhang 
94003049c21SJunchao Zhang    Logically Collective on Mat
94103049c21SJunchao Zhang 
942d8d19677SJose E. Roman    Input Parameters:
94303049c21SJunchao Zhang +  mat      - the composite matrix
94403049c21SJunchao Zhang -  scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat)
94503049c21SJunchao Zhang 
94603049c21SJunchao Zhang    Level: advanced
94703049c21SJunchao Zhang 
948db781477SPatrick Sanan .seealso: `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE`
94903049c21SJunchao Zhang 
95003049c21SJunchao Zhang @*/
95103049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings(Mat mat,const PetscScalar *scalings)
95203049c21SJunchao Zhang {
95303049c21SJunchao Zhang   PetscFunctionBegin;
95403049c21SJunchao Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
955dadcf809SJacob Faibussowitsch   PetscValidScalarPointer(scalings,2);
95603049c21SJunchao Zhang   PetscValidLogicalCollectiveScalar(mat,*scalings,2);
957cac4c232SBarry Smith   PetscUseMethod(mat,"MatCompositeSetScalings_C",(Mat,const PetscScalar*),(mat,scalings));
95803049c21SJunchao Zhang   PetscFunctionReturn(0);
95903049c21SJunchao Zhang }
96003049c21SJunchao Zhang 
961f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL,
962f4259b30SLisandro Dalcin                                        NULL,
963f4259b30SLisandro Dalcin                                        NULL,
96441cd0125SJakub Kruzik                                        MatMult_Composite,
9657bf3a71bSJakub Kruzik                                        MatMultAdd_Composite,
96641cd0125SJakub Kruzik                                /*  5*/ MatMultTranspose_Composite,
9677bf3a71bSJakub Kruzik                                        MatMultTransposeAdd_Composite,
968f4259b30SLisandro Dalcin                                        NULL,
969f4259b30SLisandro Dalcin                                        NULL,
970f4259b30SLisandro Dalcin                                        NULL,
971f4259b30SLisandro Dalcin                                /* 10*/ NULL,
972f4259b30SLisandro Dalcin                                        NULL,
973f4259b30SLisandro Dalcin                                        NULL,
974f4259b30SLisandro Dalcin                                        NULL,
975f4259b30SLisandro Dalcin                                        NULL,
976f4259b30SLisandro Dalcin                                /* 15*/ NULL,
977f4259b30SLisandro Dalcin                                        NULL,
97841cd0125SJakub Kruzik                                        MatGetDiagonal_Composite,
97941cd0125SJakub Kruzik                                        MatDiagonalScale_Composite,
980f4259b30SLisandro Dalcin                                        NULL,
981f4259b30SLisandro Dalcin                                /* 20*/ NULL,
98241cd0125SJakub Kruzik                                        MatAssemblyEnd_Composite,
983f4259b30SLisandro Dalcin                                        NULL,
984f4259b30SLisandro Dalcin                                        NULL,
985f4259b30SLisandro Dalcin                                /* 24*/ NULL,
986f4259b30SLisandro Dalcin                                        NULL,
987f4259b30SLisandro Dalcin                                        NULL,
988f4259b30SLisandro Dalcin                                        NULL,
989f4259b30SLisandro Dalcin                                        NULL,
990f4259b30SLisandro Dalcin                                /* 29*/ NULL,
991f4259b30SLisandro Dalcin                                        NULL,
992f4259b30SLisandro Dalcin                                        NULL,
993f4259b30SLisandro Dalcin                                        NULL,
994f4259b30SLisandro Dalcin                                        NULL,
995f4259b30SLisandro Dalcin                                /* 34*/ NULL,
996f4259b30SLisandro Dalcin                                        NULL,
997f4259b30SLisandro Dalcin                                        NULL,
998f4259b30SLisandro Dalcin                                        NULL,
999f4259b30SLisandro Dalcin                                        NULL,
1000f4259b30SLisandro Dalcin                                /* 39*/ NULL,
1001f4259b30SLisandro Dalcin                                        NULL,
1002f4259b30SLisandro Dalcin                                        NULL,
1003f4259b30SLisandro Dalcin                                        NULL,
1004f4259b30SLisandro Dalcin                                        NULL,
1005f4259b30SLisandro Dalcin                                /* 44*/ NULL,
100641cd0125SJakub Kruzik                                        MatScale_Composite,
100741cd0125SJakub Kruzik                                        MatShift_Basic,
1008f4259b30SLisandro Dalcin                                        NULL,
1009f4259b30SLisandro Dalcin                                        NULL,
1010f4259b30SLisandro Dalcin                                /* 49*/ NULL,
1011f4259b30SLisandro Dalcin                                        NULL,
1012f4259b30SLisandro Dalcin                                        NULL,
1013f4259b30SLisandro Dalcin                                        NULL,
1014f4259b30SLisandro Dalcin                                        NULL,
1015f4259b30SLisandro Dalcin                                /* 54*/ NULL,
1016f4259b30SLisandro Dalcin                                        NULL,
1017f4259b30SLisandro Dalcin                                        NULL,
1018f4259b30SLisandro Dalcin                                        NULL,
1019f4259b30SLisandro Dalcin                                        NULL,
1020f4259b30SLisandro Dalcin                                /* 59*/ NULL,
102141cd0125SJakub Kruzik                                        MatDestroy_Composite,
1022f4259b30SLisandro Dalcin                                        NULL,
1023f4259b30SLisandro Dalcin                                        NULL,
1024f4259b30SLisandro Dalcin                                        NULL,
1025f4259b30SLisandro Dalcin                                /* 64*/ NULL,
1026f4259b30SLisandro Dalcin                                        NULL,
1027f4259b30SLisandro Dalcin                                        NULL,
1028f4259b30SLisandro Dalcin                                        NULL,
1029f4259b30SLisandro Dalcin                                        NULL,
1030f4259b30SLisandro Dalcin                                /* 69*/ NULL,
1031f4259b30SLisandro Dalcin                                        NULL,
1032f4259b30SLisandro Dalcin                                        NULL,
1033f4259b30SLisandro Dalcin                                        NULL,
1034f4259b30SLisandro Dalcin                                        NULL,
1035f4259b30SLisandro Dalcin                                /* 74*/ NULL,
1036f4259b30SLisandro Dalcin                                        NULL,
10374b2558d6SJakub Kruzik                                        MatSetFromOptions_Composite,
1038f4259b30SLisandro Dalcin                                        NULL,
1039f4259b30SLisandro Dalcin                                        NULL,
1040f4259b30SLisandro Dalcin                                /* 79*/ NULL,
1041f4259b30SLisandro Dalcin                                        NULL,
1042f4259b30SLisandro Dalcin                                        NULL,
1043f4259b30SLisandro Dalcin                                        NULL,
1044f4259b30SLisandro Dalcin                                        NULL,
1045f4259b30SLisandro Dalcin                                /* 84*/ NULL,
1046f4259b30SLisandro Dalcin                                        NULL,
1047f4259b30SLisandro Dalcin                                        NULL,
1048f4259b30SLisandro Dalcin                                        NULL,
1049f4259b30SLisandro Dalcin                                        NULL,
1050f4259b30SLisandro Dalcin                                /* 89*/ NULL,
1051f4259b30SLisandro Dalcin                                        NULL,
1052f4259b30SLisandro Dalcin                                        NULL,
1053f4259b30SLisandro Dalcin                                        NULL,
1054f4259b30SLisandro Dalcin                                        NULL,
1055f4259b30SLisandro Dalcin                                /* 94*/ NULL,
1056f4259b30SLisandro Dalcin                                        NULL,
1057f4259b30SLisandro Dalcin                                        NULL,
1058f4259b30SLisandro Dalcin                                        NULL,
1059f4259b30SLisandro Dalcin                                        NULL,
1060f4259b30SLisandro Dalcin                                 /*99*/ NULL,
1061f4259b30SLisandro Dalcin                                        NULL,
1062f4259b30SLisandro Dalcin                                        NULL,
1063f4259b30SLisandro Dalcin                                        NULL,
1064f4259b30SLisandro Dalcin                                        NULL,
1065f4259b30SLisandro Dalcin                                /*104*/ NULL,
1066f4259b30SLisandro Dalcin                                        NULL,
1067f4259b30SLisandro Dalcin                                        NULL,
1068f4259b30SLisandro Dalcin                                        NULL,
1069f4259b30SLisandro Dalcin                                        NULL,
1070f4259b30SLisandro Dalcin                                /*109*/ NULL,
1071f4259b30SLisandro Dalcin                                        NULL,
1072f4259b30SLisandro Dalcin                                        NULL,
1073f4259b30SLisandro Dalcin                                        NULL,
1074f4259b30SLisandro Dalcin                                        NULL,
1075f4259b30SLisandro Dalcin                                /*114*/ NULL,
1076f4259b30SLisandro Dalcin                                        NULL,
1077f4259b30SLisandro Dalcin                                        NULL,
1078f4259b30SLisandro Dalcin                                        NULL,
1079f4259b30SLisandro Dalcin                                        NULL,
1080f4259b30SLisandro Dalcin                                /*119*/ NULL,
1081f4259b30SLisandro Dalcin                                        NULL,
1082f4259b30SLisandro Dalcin                                        NULL,
1083f4259b30SLisandro Dalcin                                        NULL,
1084f4259b30SLisandro Dalcin                                        NULL,
1085f4259b30SLisandro Dalcin                                /*124*/ NULL,
1086f4259b30SLisandro Dalcin                                        NULL,
1087f4259b30SLisandro Dalcin                                        NULL,
1088f4259b30SLisandro Dalcin                                        NULL,
1089f4259b30SLisandro Dalcin                                        NULL,
1090f4259b30SLisandro Dalcin                                /*129*/ NULL,
1091f4259b30SLisandro Dalcin                                        NULL,
1092f4259b30SLisandro Dalcin                                        NULL,
1093f4259b30SLisandro Dalcin                                        NULL,
1094f4259b30SLisandro Dalcin                                        NULL,
1095f4259b30SLisandro Dalcin                                /*134*/ NULL,
1096f4259b30SLisandro Dalcin                                        NULL,
1097f4259b30SLisandro Dalcin                                        NULL,
1098f4259b30SLisandro Dalcin                                        NULL,
1099f4259b30SLisandro Dalcin                                        NULL,
1100f4259b30SLisandro Dalcin                                /*139*/ NULL,
1101f4259b30SLisandro Dalcin                                        NULL,
1102d70f29a3SPierre Jolivet                                        NULL,
1103d70f29a3SPierre Jolivet                                        NULL,
1104d70f29a3SPierre Jolivet                                        NULL,
1105d70f29a3SPierre Jolivet                                 /*144*/NULL,
1106d70f29a3SPierre Jolivet                                        NULL,
1107d70f29a3SPierre Jolivet                                        NULL,
1108f4259b30SLisandro Dalcin                                        NULL
110941cd0125SJakub Kruzik };
111041cd0125SJakub Kruzik 
111141cd0125SJakub Kruzik /*MC
111241cd0125SJakub Kruzik    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices.
111341cd0125SJakub Kruzik     The matrices need to have a correct size and parallel layout for the sum or product to be valid.
111441cd0125SJakub Kruzik 
111541cd0125SJakub Kruzik    Notes:
111641cd0125SJakub Kruzik     to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE);
111741cd0125SJakub Kruzik 
111841cd0125SJakub Kruzik   Level: advanced
111941cd0125SJakub Kruzik 
1120db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`, `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()`
112141cd0125SJakub Kruzik M*/
112241cd0125SJakub Kruzik 
112341cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A)
112441cd0125SJakub Kruzik {
112541cd0125SJakub Kruzik   Mat_Composite  *b;
112641cd0125SJakub Kruzik 
112741cd0125SJakub Kruzik   PetscFunctionBegin;
11289566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A,&b));
112941cd0125SJakub Kruzik   A->data = (void*)b;
11309566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps)));
113141cd0125SJakub Kruzik 
11329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
11339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
113441cd0125SJakub Kruzik 
113541cd0125SJakub Kruzik   A->assembled    = PETSC_TRUE;
113641cd0125SJakub Kruzik   A->preallocated = PETSC_TRUE;
113741cd0125SJakub Kruzik   b->type         = MAT_COMPOSITE_ADDITIVE;
113841cd0125SJakub Kruzik   b->scale        = 1.0;
113941cd0125SJakub Kruzik   b->nmat         = 0;
11404b2558d6SJakub Kruzik   b->merge        = PETSC_FALSE;
11418c8613bfSJakub Kruzik   b->mergetype    = MAT_COMPOSITE_MERGE_RIGHT;
11423b35acafSJakub Kruzik   b->structure    = DIFFERENT_NONZERO_PATTERN;
114303049c21SJunchao Zhang   b->merge_mvctx  = PETSC_TRUE;
114403049c21SJunchao Zhang 
11459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite));
11469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite));
11479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite));
11489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeType_C",MatCompositeSetMergeType_Composite));
11499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMatStructure_C",MatCompositeSetMatStructure_Composite));
11509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMatStructure_C",MatCompositeGetMatStructure_Composite));
11519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite));
11529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNumberMat_C",MatCompositeGetNumberMat_Composite));
11539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite));
11549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetScalings_C",MatCompositeSetScalings_Composite));
115541cd0125SJakub Kruzik 
11569566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE));
115741cd0125SJakub Kruzik   PetscFunctionReturn(0);
115841cd0125SJakub Kruzik }
1159