xref: /petsc/src/mat/impls/composite/mcomposite.c (revision f4259b3095b7e7c300c2b3c3729548b071684bd7)
1793850ffSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
3793850ffSBarry Smith 
4*f4259b30SLisandro 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 {
37793850ffSBarry Smith   PetscErrorCode    ierr;
382c33bbaeSSatish Balay   Mat_Composite     *shell = (Mat_Composite*)mat->data;
396d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head,oldnext;
4003049c21SJunchao Zhang   PetscInt          i;
41793850ffSBarry Smith 
42793850ffSBarry Smith   PetscFunctionBegin;
43793850ffSBarry Smith   while (next) {
446bf464f9SBarry Smith     ierr = MatDestroy(&next->mat);CHKERRQ(ierr);
456d7c1e57SBarry Smith     if (next->work && (!next->next || next->work != next->next->work)) {
466bf464f9SBarry Smith       ierr = VecDestroy(&next->work);CHKERRQ(ierr);
476d7c1e57SBarry Smith     }
486d7c1e57SBarry Smith     oldnext = next;
49793850ffSBarry Smith     next    = next->next;
506d7c1e57SBarry Smith     ierr    = PetscFree(oldnext);CHKERRQ(ierr);
51793850ffSBarry Smith   }
526bf464f9SBarry Smith   ierr = VecDestroy(&shell->work);CHKERRQ(ierr);
536bf464f9SBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
546bf464f9SBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
556bf464f9SBarry Smith   ierr = VecDestroy(&shell->leftwork);CHKERRQ(ierr);
566bf464f9SBarry Smith   ierr = VecDestroy(&shell->rightwork);CHKERRQ(ierr);
5703049c21SJunchao Zhang   ierr = VecDestroy(&shell->leftwork2);CHKERRQ(ierr);
5803049c21SJunchao Zhang   ierr = VecDestroy(&shell->rightwork2);CHKERRQ(ierr);
5903049c21SJunchao Zhang 
6003049c21SJunchao Zhang   if (shell->Mvctx) {
6103049c21SJunchao Zhang     for (i=0; i<shell->nmat; i++) {ierr = VecDestroy(&shell->lvecs[i]);CHKERRQ(ierr);}
6203049c21SJunchao Zhang     ierr = PetscFree3(shell->location,shell->larray,shell->lvecs);CHKERRQ(ierr);
6303049c21SJunchao Zhang     ierr = PetscFree(shell->larray);CHKERRQ(ierr);
6403049c21SJunchao Zhang     ierr = VecDestroy(&shell->gvec);CHKERRQ(ierr);
6503049c21SJunchao Zhang     ierr = VecScatterDestroy(&shell->Mvctx);CHKERRQ(ierr);
6603049c21SJunchao Zhang   }
6703049c21SJunchao Zhang 
6803049c21SJunchao Zhang   ierr = PetscFree(shell->scalings);CHKERRQ(ierr);
696bf464f9SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
70793850ffSBarry Smith   PetscFunctionReturn(0);
71793850ffSBarry Smith }
72793850ffSBarry Smith 
736d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y)
746d7c1e57SBarry Smith {
756d7c1e57SBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
766d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head;
776d7c1e57SBarry Smith   PetscErrorCode    ierr;
786d7c1e57SBarry Smith   Vec               in,out;
7903049c21SJunchao Zhang   PetscScalar       scale;
8003049c21SJunchao Zhang   PetscInt          i;
816d7c1e57SBarry Smith 
826d7c1e57SBarry Smith   PetscFunctionBegin;
83e32f2f54SBarry Smith   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
846d7c1e57SBarry Smith   in = x;
85e4fc5df0SBarry Smith   if (shell->right) {
86e4fc5df0SBarry Smith     if (!shell->rightwork) {
87e4fc5df0SBarry Smith       ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr);
88e4fc5df0SBarry Smith     }
89e4fc5df0SBarry Smith     ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr);
90e4fc5df0SBarry Smith     in   = shell->rightwork;
91e4fc5df0SBarry Smith   }
926d7c1e57SBarry Smith   while (next->next) {
936d7c1e57SBarry Smith     if (!next->work) { /* should reuse previous work if the same size */
942a7a6963SBarry Smith       ierr = MatCreateVecs(next->mat,NULL,&next->work);CHKERRQ(ierr);
956d7c1e57SBarry Smith     }
966d7c1e57SBarry Smith     out  = next->work;
976d7c1e57SBarry Smith     ierr = MatMult(next->mat,in,out);CHKERRQ(ierr);
986d7c1e57SBarry Smith     in   = out;
996d7c1e57SBarry Smith     next = next->next;
1006d7c1e57SBarry Smith   }
1016d7c1e57SBarry Smith   ierr = MatMult(next->mat,in,y);CHKERRQ(ierr);
102e4fc5df0SBarry Smith   if (shell->left) {
103e4fc5df0SBarry Smith     ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr);
104e4fc5df0SBarry Smith   }
10503049c21SJunchao Zhang   scale = shell->scale;
10603049c21SJunchao Zhang   if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
10703049c21SJunchao Zhang   ierr = VecScale(y,scale);CHKERRQ(ierr);
1086d7c1e57SBarry Smith   PetscFunctionReturn(0);
1096d7c1e57SBarry Smith }
1106d7c1e57SBarry Smith 
1116d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y)
1126d7c1e57SBarry Smith {
1136d7c1e57SBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
1146d7c1e57SBarry Smith   Mat_CompositeLink tail   = shell->tail;
1156d7c1e57SBarry Smith   PetscErrorCode    ierr;
1166d7c1e57SBarry Smith   Vec               in,out;
11703049c21SJunchao Zhang   PetscScalar       scale;
11803049c21SJunchao Zhang   PetscInt          i;
1196d7c1e57SBarry Smith 
1206d7c1e57SBarry Smith   PetscFunctionBegin;
121e32f2f54SBarry Smith   if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
1226d7c1e57SBarry Smith   in = x;
123e4fc5df0SBarry Smith   if (shell->left) {
124e4fc5df0SBarry Smith     if (!shell->leftwork) {
125e4fc5df0SBarry Smith       ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr);
126e4fc5df0SBarry Smith     }
127e4fc5df0SBarry Smith     ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr);
128e4fc5df0SBarry Smith     in   = shell->leftwork;
129e4fc5df0SBarry Smith   }
1306d7c1e57SBarry Smith   while (tail->prev) {
1316d7c1e57SBarry Smith     if (!tail->prev->work) { /* should reuse previous work if the same size */
1322a7a6963SBarry Smith       ierr = MatCreateVecs(tail->mat,NULL,&tail->prev->work);CHKERRQ(ierr);
1336d7c1e57SBarry Smith     }
1346d7c1e57SBarry Smith     out  = tail->prev->work;
1356d7c1e57SBarry Smith     ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr);
1366d7c1e57SBarry Smith     in   = out;
1376d7c1e57SBarry Smith     tail = tail->prev;
1386d7c1e57SBarry Smith   }
1396d7c1e57SBarry Smith   ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr);
140e4fc5df0SBarry Smith   if (shell->right) {
141e4fc5df0SBarry Smith     ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr);
142e4fc5df0SBarry Smith   }
14303049c21SJunchao Zhang 
14403049c21SJunchao Zhang   scale = shell->scale;
14503049c21SJunchao Zhang   if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
14603049c21SJunchao Zhang   ierr = VecScale(y,scale);CHKERRQ(ierr);
1476d7c1e57SBarry Smith   PetscFunctionReturn(0);
1486d7c1e57SBarry Smith }
1496d7c1e57SBarry Smith 
15003049c21SJunchao Zhang PetscErrorCode MatMult_Composite(Mat mat,Vec x,Vec y)
151793850ffSBarry Smith {
152793850ffSBarry Smith   PetscErrorCode    ierr;
15303049c21SJunchao Zhang   Mat_Composite     *shell = (Mat_Composite*)mat->data;
15403049c21SJunchao Zhang   Mat_CompositeLink cur = shell->head;
15503049c21SJunchao Zhang   Vec               in,y2,xin;
15603049c21SJunchao Zhang   Mat               A,B;
15703049c21SJunchao Zhang   PetscInt          i,j,k,n,nuniq,lo,hi,mid,*gindices,*buf,*tmp,tot;
15803049c21SJunchao Zhang   const PetscScalar *vals;
15903049c21SJunchao Zhang   const PetscInt    *garray;
16003049c21SJunchao Zhang   IS                ix,iy;
16103049c21SJunchao Zhang   PetscBool         match;
162793850ffSBarry Smith 
163793850ffSBarry Smith   PetscFunctionBegin;
16403049c21SJunchao Zhang   if (!cur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
165e4fc5df0SBarry Smith   in = x;
166e4fc5df0SBarry Smith   if (shell->right) {
167e4fc5df0SBarry Smith     if (!shell->rightwork) {
168e4fc5df0SBarry Smith       ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr);
169793850ffSBarry Smith     }
170e4fc5df0SBarry Smith     ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr);
171e4fc5df0SBarry Smith     in   = shell->rightwork;
172e4fc5df0SBarry Smith   }
17303049c21SJunchao Zhang 
17403049c21SJunchao Zhang   /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time
17503049c21SJunchao Zhang      we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and
17603049c21SJunchao Zhang      it is legal to merge Mvctx, because all component matrices have the same size.
17703049c21SJunchao Zhang    */
17803049c21SJunchao Zhang   if (shell->merge_mvctx && !shell->Mvctx) {
17903049c21SJunchao Zhang     /* Currently only implemented for MATMPIAIJ */
18003049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
18103049c21SJunchao Zhang       ierr = PetscObjectTypeCompare((PetscObject)cur->mat,MATMPIAIJ,&match);CHKERRQ(ierr);
18203049c21SJunchao Zhang       if (!match) {
18303049c21SJunchao Zhang         shell->merge_mvctx = PETSC_FALSE;
18403049c21SJunchao Zhang         goto skip_merge_mvctx;
185e4fc5df0SBarry Smith       }
186e4fc5df0SBarry Smith     }
18703049c21SJunchao Zhang 
18803049c21SJunchao Zhang     /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */
18903049c21SJunchao Zhang     tot = 0;
19003049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
19103049c21SJunchao Zhang       ierr = MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,NULL);CHKERRQ(ierr);
19203049c21SJunchao Zhang       ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr);
19303049c21SJunchao Zhang       tot += n;
19403049c21SJunchao Zhang     }
19503049c21SJunchao Zhang     ierr = PetscMalloc3(tot,&shell->location,tot,&shell->larray,shell->nmat,&shell->lvecs);CHKERRQ(ierr);
19603049c21SJunchao Zhang     shell->len = tot;
19703049c21SJunchao Zhang 
19803049c21SJunchao Zhang     /* Go through matrices second time to sort off-diag columns and remove dups */
19903049c21SJunchao Zhang     ierr  = PetscMalloc1(tot,&gindices);CHKERRQ(ierr); /* No Malloc2() since we will give one to petsc and free the other */
20003049c21SJunchao Zhang     ierr  = PetscMalloc1(tot,&buf);CHKERRQ(ierr);
20103049c21SJunchao Zhang     nuniq = 0; /* Number of unique nonzero columns */
20203049c21SJunchao Zhang     for (cur=shell->head; cur; cur=cur->next) {
20303049c21SJunchao Zhang       ierr = MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray);CHKERRQ(ierr);
20403049c21SJunchao Zhang       ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr);
20503049c21SJunchao Zhang       /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */
20603049c21SJunchao Zhang       i = j = k = 0;
20703049c21SJunchao Zhang       while (i < n && j < nuniq) {
20803049c21SJunchao Zhang         if (garray[i] < gindices[j]) buf[k++] = garray[i++];
20903049c21SJunchao Zhang         else if (garray[i] > gindices[j]) buf[k++] = gindices[j++];
21003049c21SJunchao Zhang         else {buf[k++] = garray[i++]; j++;}
21103049c21SJunchao Zhang       }
21203049c21SJunchao Zhang       /* Copy leftover in garray[] or gindices[] */
21303049c21SJunchao Zhang       if (i < n) {
21403049c21SJunchao Zhang         ierr  = PetscArraycpy(buf+k,garray+i,n-i);CHKERRQ(ierr);
21503049c21SJunchao Zhang         nuniq = k + n-i;
21603049c21SJunchao Zhang       } else if (j < nuniq) {
21703049c21SJunchao Zhang         ierr  = PetscArraycpy(buf+k,gindices+j,nuniq-j);CHKERRQ(ierr);
21803049c21SJunchao Zhang         nuniq = k + nuniq-j;
21903049c21SJunchao Zhang       } else nuniq = k;
22003049c21SJunchao Zhang       /* Swap gindices and buf to merge garray of the next matrix */
22103049c21SJunchao Zhang       tmp      = gindices;
22203049c21SJunchao Zhang       gindices = buf;
22303049c21SJunchao Zhang       buf      = tmp;
22403049c21SJunchao Zhang     }
22503049c21SJunchao Zhang     ierr = PetscFree(buf);CHKERRQ(ierr);
22603049c21SJunchao Zhang 
22703049c21SJunchao Zhang     /* Go through matrices third time to build a map from gindices[] to garray[] */
22803049c21SJunchao Zhang     tot = 0;
22903049c21SJunchao Zhang     for (cur=shell->head,j=0; cur; cur=cur->next,j++) { /* j-th matrix */
23003049c21SJunchao Zhang       ierr = MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray);CHKERRQ(ierr);
23103049c21SJunchao Zhang       ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr);
23203049c21SJunchao Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&shell->lvecs[j]);CHKERRQ(ierr);
23303049c21SJunchao Zhang       /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */
23403049c21SJunchao Zhang       lo   = 0;
23503049c21SJunchao Zhang       for (i=0; i<n; i++) {
23603049c21SJunchao Zhang         hi = nuniq;
23703049c21SJunchao Zhang         while (hi - lo > 1) {
23803049c21SJunchao Zhang           mid = lo + (hi - lo)/2;
23903049c21SJunchao Zhang           if (garray[i] < gindices[mid]) hi = mid;
24003049c21SJunchao Zhang           else lo = mid;
24103049c21SJunchao Zhang         }
24203049c21SJunchao Zhang         shell->location[tot+i] = lo; /* gindices[lo] = garray[i] */
24303049c21SJunchao Zhang         lo++; /* Since garray[i+1] > garray[i], we can safely advance lo */
24403049c21SJunchao Zhang       }
24503049c21SJunchao Zhang       tot += n;
24603049c21SJunchao Zhang     }
24703049c21SJunchao Zhang 
24803049c21SJunchao Zhang     /* Build merged Mvctx */
24903049c21SJunchao Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,nuniq,gindices,PETSC_OWN_POINTER,&ix);CHKERRQ(ierr);
25003049c21SJunchao Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,nuniq,0,1,&iy);CHKERRQ(ierr);
25103049c21SJunchao Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&xin);CHKERRQ(ierr);
25203049c21SJunchao Zhang     ierr = VecCreateSeq(PETSC_COMM_SELF,nuniq,&shell->gvec);CHKERRQ(ierr);
25303049c21SJunchao Zhang     ierr = VecScatterCreate(xin,ix,shell->gvec,iy,&shell->Mvctx);CHKERRQ(ierr);
25403049c21SJunchao Zhang     ierr = VecDestroy(&xin);CHKERRQ(ierr);
25503049c21SJunchao Zhang     ierr = ISDestroy(&ix);CHKERRQ(ierr);
25603049c21SJunchao Zhang     ierr = ISDestroy(&iy);CHKERRQ(ierr);
25703049c21SJunchao Zhang   }
25803049c21SJunchao Zhang 
25903049c21SJunchao Zhang skip_merge_mvctx:
26003049c21SJunchao Zhang   ierr = VecSet(y,0);CHKERRQ(ierr);
26103049c21SJunchao Zhang   if (!shell->leftwork2) {ierr = VecDuplicate(y,&shell->leftwork2);CHKERRQ(ierr);}
26203049c21SJunchao Zhang   y2 = shell->leftwork2;
26303049c21SJunchao Zhang 
26403049c21SJunchao Zhang   if (shell->Mvctx) { /* Have a merged Mvctx */
26503049c21SJunchao 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
26603049c21SJunchao 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
26703049c21SJunchao Zhang        overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach.
26803049c21SJunchao Zhang      */
26903049c21SJunchao Zhang     ierr = VecScatterBegin(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
27003049c21SJunchao Zhang     ierr = VecScatterEnd(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
27103049c21SJunchao Zhang 
27203049c21SJunchao Zhang     ierr = VecGetArrayRead(shell->gvec,&vals);CHKERRQ(ierr);
27303049c21SJunchao Zhang     for (i=0; i<shell->len; i++) shell->larray[i] = vals[shell->location[i]];
27403049c21SJunchao Zhang     ierr = VecRestoreArrayRead(shell->gvec,&vals);CHKERRQ(ierr);
27503049c21SJunchao Zhang 
27603049c21SJunchao Zhang     for (cur=shell->head,tot=i=0; cur; cur=cur->next,i++) { /* i-th matrix */
27703049c21SJunchao Zhang       ierr = MatMPIAIJGetSeqAIJ(cur->mat,&A,&B,NULL);CHKERRQ(ierr);
27803049c21SJunchao Zhang       ierr = (*A->ops->mult)(A,in,y2);CHKERRQ(ierr);
27903049c21SJunchao Zhang       ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr);
28003049c21SJunchao Zhang       ierr = VecPlaceArray(shell->lvecs[i],&shell->larray[tot]);CHKERRQ(ierr);
28103049c21SJunchao Zhang       ierr = (*B->ops->multadd)(B,shell->lvecs[i],y2,y2);CHKERRQ(ierr);
28203049c21SJunchao Zhang       ierr = VecResetArray(shell->lvecs[i]);CHKERRQ(ierr);
28303049c21SJunchao Zhang       ierr = VecAXPY(y,(shell->scalings ? shell->scalings[i] : 1.0),y2);CHKERRQ(ierr);
28403049c21SJunchao Zhang       tot += n;
28503049c21SJunchao Zhang     }
28603049c21SJunchao Zhang   } else {
28703049c21SJunchao Zhang     if (shell->scalings) {
28803049c21SJunchao Zhang       for (cur=shell->head,i=0; cur; cur=cur->next,i++) {
28903049c21SJunchao Zhang         ierr = MatMult(cur->mat,in,y2);CHKERRQ(ierr);
29003049c21SJunchao Zhang         ierr = VecAXPY(y,shell->scalings[i],y2);CHKERRQ(ierr);
29103049c21SJunchao Zhang       }
29203049c21SJunchao Zhang     } else {
29303049c21SJunchao Zhang       for (cur=shell->head; cur; cur=cur->next) {ierr = MatMultAdd(cur->mat,in,y,y);CHKERRQ(ierr);}
29403049c21SJunchao Zhang     }
29503049c21SJunchao Zhang   }
29603049c21SJunchao Zhang 
29703049c21SJunchao Zhang   if (shell->left) {ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr);}
298e4fc5df0SBarry Smith   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
299793850ffSBarry Smith   PetscFunctionReturn(0);
300793850ffSBarry Smith }
301793850ffSBarry Smith 
302793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y)
303793850ffSBarry Smith {
304793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
305793850ffSBarry Smith   Mat_CompositeLink next   = shell->head;
306793850ffSBarry Smith   PetscErrorCode    ierr;
30703049c21SJunchao Zhang   Vec               in,y2 = NULL;
30803049c21SJunchao Zhang   PetscInt          i;
309793850ffSBarry Smith 
310793850ffSBarry Smith   PetscFunctionBegin;
311e32f2f54SBarry Smith   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
312e4fc5df0SBarry Smith   in = x;
313e4fc5df0SBarry Smith   if (shell->left) {
314e4fc5df0SBarry Smith     if (!shell->leftwork) {
315e4fc5df0SBarry Smith       ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr);
316793850ffSBarry Smith     }
317e4fc5df0SBarry Smith     ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr);
318e4fc5df0SBarry Smith     in   = shell->leftwork;
319e4fc5df0SBarry Smith   }
32003049c21SJunchao Zhang 
321e4fc5df0SBarry Smith   ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr);
32203049c21SJunchao Zhang   if (shell->scalings) {
32303049c21SJunchao Zhang     ierr = VecScale(y,shell->scalings[0]);CHKERRQ(ierr);
32403049c21SJunchao Zhang     if (!shell->rightwork2) {ierr = VecDuplicate(y,&shell->rightwork2);CHKERRQ(ierr);}
32503049c21SJunchao Zhang     y2 = shell->rightwork2;
32603049c21SJunchao Zhang   }
32703049c21SJunchao Zhang   i = 1;
328e4fc5df0SBarry Smith   while ((next = next->next)) {
32903049c21SJunchao Zhang     if (!shell->scalings) {ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr);}
33003049c21SJunchao Zhang     else {
33103049c21SJunchao Zhang       ierr = MatMultTranspose(next->mat,in,y2);CHKERRQ(ierr);
33203049c21SJunchao Zhang       ierr = VecAXPY(y,shell->scalings[i++],y2);CHKERRQ(ierr);
33303049c21SJunchao Zhang     }
334e4fc5df0SBarry Smith   }
335e4fc5df0SBarry Smith   if (shell->right) {
336e4fc5df0SBarry Smith     ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr);
337e4fc5df0SBarry Smith   }
338e4fc5df0SBarry Smith   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
339793850ffSBarry Smith   PetscFunctionReturn(0);
340793850ffSBarry Smith }
341793850ffSBarry Smith 
3427bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z)
3437bf3a71bSJakub Kruzik {
3447bf3a71bSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)A->data;
3457bf3a71bSJakub Kruzik   PetscErrorCode    ierr;
3467bf3a71bSJakub Kruzik 
3477bf3a71bSJakub Kruzik   PetscFunctionBegin;
3487bf3a71bSJakub Kruzik   if (y != z) {
3497bf3a71bSJakub Kruzik     ierr = MatMult(A,x,z);CHKERRQ(ierr);
3507bf3a71bSJakub Kruzik     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3517bf3a71bSJakub Kruzik   } else {
3527bf3a71bSJakub Kruzik     if (!shell->leftwork) {
3537bf3a71bSJakub Kruzik       ierr = VecDuplicate(z,&shell->leftwork);CHKERRQ(ierr);
3547bf3a71bSJakub Kruzik     }
3557bf3a71bSJakub Kruzik     ierr = MatMult(A,x,shell->leftwork);CHKERRQ(ierr);
3567bf3a71bSJakub Kruzik     ierr = VecCopy(y,z);CHKERRQ(ierr);
3577bf3a71bSJakub Kruzik     ierr = VecAXPY(z,1.0,shell->leftwork);CHKERRQ(ierr);
3587bf3a71bSJakub Kruzik   }
3597bf3a71bSJakub Kruzik   PetscFunctionReturn(0);
3607bf3a71bSJakub Kruzik }
3617bf3a71bSJakub Kruzik 
3627bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z)
3637bf3a71bSJakub Kruzik {
3647bf3a71bSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)A->data;
3657bf3a71bSJakub Kruzik   PetscErrorCode    ierr;
3667bf3a71bSJakub Kruzik 
3677bf3a71bSJakub Kruzik   PetscFunctionBegin;
3687bf3a71bSJakub Kruzik   if (y != z) {
3697bf3a71bSJakub Kruzik     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3707bf3a71bSJakub Kruzik     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3717bf3a71bSJakub Kruzik   } else {
3727bf3a71bSJakub Kruzik     if (!shell->rightwork) {
3737bf3a71bSJakub Kruzik       ierr = VecDuplicate(z,&shell->rightwork);CHKERRQ(ierr);
3747bf3a71bSJakub Kruzik     }
3757bf3a71bSJakub Kruzik     ierr = MatMultTranspose(A,x,shell->rightwork);CHKERRQ(ierr);
3767bf3a71bSJakub Kruzik     ierr = VecCopy(y,z);CHKERRQ(ierr);
3777bf3a71bSJakub Kruzik     ierr = VecAXPY(z,1.0,shell->rightwork);CHKERRQ(ierr);
3787bf3a71bSJakub Kruzik   }
3797bf3a71bSJakub Kruzik   PetscFunctionReturn(0);
3807bf3a71bSJakub Kruzik }
3817bf3a71bSJakub Kruzik 
382793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v)
383793850ffSBarry Smith {
384793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
385793850ffSBarry Smith   Mat_CompositeLink next   = shell->head;
386793850ffSBarry Smith   PetscErrorCode    ierr;
38703049c21SJunchao Zhang   PetscInt          i;
388793850ffSBarry Smith 
389793850ffSBarry Smith   PetscFunctionBegin;
390e32f2f54SBarry Smith   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
391e32f2f54SBarry Smith   if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling");
392e4fc5df0SBarry Smith 
393793850ffSBarry Smith   ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr);
39403049c21SJunchao Zhang   if (shell->scalings) {ierr = VecScale(v,shell->scalings[0]);CHKERRQ(ierr);}
39503049c21SJunchao Zhang 
396793850ffSBarry Smith   if (next->next && !shell->work) {
397793850ffSBarry Smith     ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr);
398793850ffSBarry Smith   }
39903049c21SJunchao Zhang   i = 1;
400793850ffSBarry Smith   while ((next = next->next)) {
401793850ffSBarry Smith     ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr);
40203049c21SJunchao Zhang     ierr = VecAXPY(v,(shell->scalings ? shell->scalings[i++] : 1.0),shell->work);CHKERRQ(ierr);
403793850ffSBarry Smith   }
404e4fc5df0SBarry Smith   ierr = VecScale(v,shell->scale);CHKERRQ(ierr);
405793850ffSBarry Smith   PetscFunctionReturn(0);
406793850ffSBarry Smith }
407793850ffSBarry Smith 
408793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)
409793850ffSBarry Smith {
4104b2558d6SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)Y->data;
411b52f573bSBarry Smith   PetscErrorCode    ierr;
412b52f573bSBarry Smith 
413793850ffSBarry Smith   PetscFunctionBegin;
4144b2558d6SJakub Kruzik   if (shell->merge) {
415b52f573bSBarry Smith     ierr = MatCompositeMerge(Y);CHKERRQ(ierr);
416b52f573bSBarry Smith   }
41703049c21SJunchao Zhang 
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   PetscErrorCode ierr;
434e4fc5df0SBarry Smith 
435e4fc5df0SBarry Smith   PetscFunctionBegin;
436e4fc5df0SBarry Smith   if (left) {
437321429dbSBarry Smith     if (!a->left) {
438e4fc5df0SBarry Smith       ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr);
439e4fc5df0SBarry Smith       ierr = VecCopy(left,a->left);CHKERRQ(ierr);
440321429dbSBarry Smith     } else {
441321429dbSBarry Smith       ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr);
442321429dbSBarry Smith     }
443e4fc5df0SBarry Smith   }
444e4fc5df0SBarry Smith   if (right) {
445321429dbSBarry Smith     if (!a->right) {
446e4fc5df0SBarry Smith       ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr);
447e4fc5df0SBarry Smith       ierr = VecCopy(right,a->right);CHKERRQ(ierr);
448321429dbSBarry Smith     } else {
449321429dbSBarry Smith       ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr);
450321429dbSBarry Smith     }
451e4fc5df0SBarry Smith   }
452e4fc5df0SBarry Smith   PetscFunctionReturn(0);
453e4fc5df0SBarry Smith }
454793850ffSBarry Smith 
4554b2558d6SJakub Kruzik PetscErrorCode MatSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,Mat A)
4564b2558d6SJakub Kruzik {
4574b2558d6SJakub Kruzik   Mat_Composite  *a = (Mat_Composite*)A->data;
4584b2558d6SJakub Kruzik   PetscErrorCode ierr;
4594b2558d6SJakub Kruzik 
4604b2558d6SJakub Kruzik   PetscFunctionBegin;
4614b2558d6SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"MATCOMPOSITE options");CHKERRQ(ierr);
4624b2558d6SJakub Kruzik   ierr = PetscOptionsBool("-mat_composite_merge","Merge at MatAssemblyEnd","MatCompositeMerge",a->merge,&a->merge,NULL);CHKERRQ(ierr);
4638c8613bfSJakub Kruzik   ierr = PetscOptionsEnum("-mat_composite_merge_type","Set composite merge direction","MatCompositeSetMergeType",MatCompositeMergeTypes,(PetscEnum)a->mergetype,(PetscEnum*)&a->mergetype,NULL);CHKERRQ(ierr);
46403049c21SJunchao Zhang   ierr = PetscOptionsBool("-mat_composite_merge_mvctx","Merge MatMult() vecscat contexts","MatCreateComposite",a->merge_mvctx,&a->merge_mvctx,NULL);CHKERRQ(ierr);
4654b2558d6SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
4664b2558d6SJakub Kruzik   PetscFunctionReturn(0);
4674b2558d6SJakub Kruzik }
4684b2558d6SJakub Kruzik 
4692c0821f3SBarry Smith /*@
4708bd739bdSJakub Kruzik    MatCreateComposite - Creates a matrix as the sum or product of one or more matrices
471793850ffSBarry Smith 
472d083f849SBarry Smith   Collective
473793850ffSBarry Smith 
474793850ffSBarry Smith    Input Parameters:
475793850ffSBarry Smith +  comm - MPI communicator
476793850ffSBarry Smith .  nmat - number of matrices to put in
477793850ffSBarry Smith -  mats - the matrices
478793850ffSBarry Smith 
479793850ffSBarry Smith    Output Parameter:
480db36af27SMatthew G. Knepley .  mat - the matrix
481793850ffSBarry Smith 
4824b2558d6SJakub Kruzik    Options Database Keys:
4834b2558d6SJakub Kruzik +  -mat_composite_merge         - merge in MatAssemblyEnd()
48403049c21SJunchao Zhang .  -mat_composite_merge_mvctx   - merge Mvctx of component matrices to optimize communication in MatMult() for ADDITIVE matrices
485b28d0daaSJakub Kruzik -  -mat_composite_merge_type    - set merge direction
4864b2558d6SJakub Kruzik 
487793850ffSBarry Smith    Level: advanced
488793850ffSBarry Smith 
489793850ffSBarry Smith    Notes:
490793850ffSBarry Smith      Alternative construction
491793850ffSBarry Smith $       MatCreate(comm,&mat);
492793850ffSBarry Smith $       MatSetSizes(mat,m,n,M,N);
493793850ffSBarry Smith $       MatSetType(mat,MATCOMPOSITE);
494793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[0]);
495793850ffSBarry Smith $       ....
496793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[nmat-1]);
497b52f573bSBarry Smith $       MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
498b52f573bSBarry Smith $       MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
499793850ffSBarry Smith 
5006d7c1e57SBarry Smith      For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0]
5016d7c1e57SBarry Smith 
5026dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeGetMat(), MatCompositeMerge(), MatCompositeSetType(), MATCOMPOSITE
503793850ffSBarry Smith 
504793850ffSBarry Smith @*/
5057087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat)
506793850ffSBarry Smith {
507793850ffSBarry Smith   PetscErrorCode ierr;
508793850ffSBarry Smith   PetscInt       m,n,M,N,i;
509793850ffSBarry Smith 
510793850ffSBarry Smith   PetscFunctionBegin;
511e32f2f54SBarry Smith   if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix");
512f3f84630SBarry Smith   PetscValidPointer(mat,3);
513793850ffSBarry Smith 
514c4f13f43SJakub Kruzik   ierr = MatGetLocalSize(mats[0],PETSC_IGNORE,&n);CHKERRQ(ierr);
515c4f13f43SJakub Kruzik   ierr = MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE);CHKERRQ(ierr);
516c4f13f43SJakub Kruzik   ierr = MatGetSize(mats[0],PETSC_IGNORE,&N);CHKERRQ(ierr);
517c4f13f43SJakub Kruzik   ierr = MatGetSize(mats[nmat-1],&M,PETSC_IGNORE);CHKERRQ(ierr);
518793850ffSBarry Smith   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
519793850ffSBarry Smith   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
520793850ffSBarry Smith   ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr);
521793850ffSBarry Smith   for (i=0; i<nmat; i++) {
522793850ffSBarry Smith     ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr);
523793850ffSBarry Smith   }
524b52f573bSBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
525b52f573bSBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
526793850ffSBarry Smith   PetscFunctionReturn(0);
527793850ffSBarry Smith }
528793850ffSBarry Smith 
529d7f81bf2SJakub Kruzik 
530d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat)
531d7f81bf2SJakub Kruzik {
532d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
533d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink,next = shell->head;
534d7f81bf2SJakub Kruzik   PetscErrorCode    ierr;
535d7f81bf2SJakub Kruzik 
536d7f81bf2SJakub Kruzik   PetscFunctionBegin;
537d7f81bf2SJakub Kruzik   ierr        = PetscNewLog(mat,&ilink);CHKERRQ(ierr);
538*f4259b30SLisandro Dalcin   ilink->next = NULL;
539d7f81bf2SJakub Kruzik   ierr        = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr);
540d7f81bf2SJakub Kruzik   ilink->mat  = smat;
541d7f81bf2SJakub Kruzik 
542d7f81bf2SJakub Kruzik   if (!next) shell->head = ilink;
543d7f81bf2SJakub Kruzik   else {
544d7f81bf2SJakub Kruzik     while (next->next) {
545d7f81bf2SJakub Kruzik       next = next->next;
546d7f81bf2SJakub Kruzik     }
547d7f81bf2SJakub Kruzik     next->next  = ilink;
548d7f81bf2SJakub Kruzik     ilink->prev = next;
549d7f81bf2SJakub Kruzik   }
550d7f81bf2SJakub Kruzik   shell->tail =  ilink;
551d7f81bf2SJakub Kruzik   shell->nmat += 1;
55203049c21SJunchao Zhang 
55303049c21SJunchao Zhang   /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */
55403049c21SJunchao Zhang   if (shell->scalings) {
55503049c21SJunchao Zhang     ierr = PetscRealloc(sizeof(PetscScalar)*shell->nmat,&shell->scalings);CHKERRQ(ierr);
55603049c21SJunchao Zhang     shell->scalings[shell->nmat-1] = 1.0;
55703049c21SJunchao Zhang   }
558d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
559d7f81bf2SJakub Kruzik }
560d7f81bf2SJakub Kruzik 
561793850ffSBarry Smith /*@
5628bd739bdSJakub Kruzik     MatCompositeAddMat - Add another matrix to a composite matrix.
563793850ffSBarry Smith 
564793850ffSBarry Smith    Collective on Mat
565793850ffSBarry Smith 
566793850ffSBarry Smith     Input Parameters:
567793850ffSBarry Smith +   mat - the composite matrix
568793850ffSBarry Smith -   smat - the partial matrix
569793850ffSBarry Smith 
570793850ffSBarry Smith    Level: advanced
571793850ffSBarry Smith 
5728bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE
573793850ffSBarry Smith @*/
5747087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat)
575793850ffSBarry Smith {
576793850ffSBarry Smith   PetscErrorCode    ierr;
577793850ffSBarry Smith 
578793850ffSBarry Smith   PetscFunctionBegin;
5790700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
5800700a824SBarry Smith   PetscValidHeaderSpecific(smat,MAT_CLASSID,2);
581d7f81bf2SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));CHKERRQ(ierr);
582d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
583d7f81bf2SJakub Kruzik }
584793850ffSBarry Smith 
585d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type)
586d7f81bf2SJakub Kruzik {
587d7f81bf2SJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
588d7f81bf2SJakub Kruzik 
589d7f81bf2SJakub Kruzik   PetscFunctionBegin;
590ced1ac25SJakub Kruzik   b->type = type;
591d7f81bf2SJakub Kruzik   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
592*f4259b30SLisandro Dalcin     mat->ops->getdiagonal   = NULL;
593d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite_Multiplicative;
594d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative;
59503049c21SJunchao Zhang     b->merge_mvctx          = PETSC_FALSE;
596d7f81bf2SJakub Kruzik   } else {
597d7f81bf2SJakub Kruzik     mat->ops->getdiagonal   = MatGetDiagonal_Composite;
598d7f81bf2SJakub Kruzik     mat->ops->mult          = MatMult_Composite;
599d7f81bf2SJakub Kruzik     mat->ops->multtranspose = MatMultTranspose_Composite;
600793850ffSBarry Smith   }
6016d7c1e57SBarry Smith   PetscFunctionReturn(0);
6026d7c1e57SBarry Smith }
6036d7c1e57SBarry Smith 
6042c0821f3SBarry Smith /*@
6058bd739bdSJakub Kruzik    MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product.
6066d7c1e57SBarry Smith 
607b2b245abSJakub Kruzik    Logically Collective on Mat
6086d7c1e57SBarry Smith 
6096d7c1e57SBarry Smith    Input Parameters:
6106d7c1e57SBarry Smith .  mat - the composite matrix
6116d7c1e57SBarry Smith 
6126d7c1e57SBarry Smith    Level: advanced
6136d7c1e57SBarry Smith 
6146dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeGetType(), MATCOMPOSITE
6156d7c1e57SBarry Smith 
6166d7c1e57SBarry Smith @*/
6177087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type)
6186d7c1e57SBarry Smith {
6196d7c1e57SBarry Smith   PetscErrorCode ierr;
6206d7c1e57SBarry Smith 
6216d7c1e57SBarry Smith   PetscFunctionBegin;
622d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
623b2b245abSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat,type,2);
624d7f81bf2SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));CHKERRQ(ierr);
625793850ffSBarry Smith   PetscFunctionReturn(0);
626793850ffSBarry Smith }
627793850ffSBarry Smith 
6286dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type)
6296dbc55e5SJakub Kruzik {
6306dbc55e5SJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
6316dbc55e5SJakub Kruzik 
6326dbc55e5SJakub Kruzik   PetscFunctionBegin;
6336dbc55e5SJakub Kruzik   *type = b->type;
6346dbc55e5SJakub Kruzik   PetscFunctionReturn(0);
6356dbc55e5SJakub Kruzik }
6366dbc55e5SJakub Kruzik 
6376dbc55e5SJakub Kruzik /*@
6386dbc55e5SJakub Kruzik    MatCompositeGetType - Returns type of composite.
6396dbc55e5SJakub Kruzik 
6406dbc55e5SJakub Kruzik    Not Collective
6416dbc55e5SJakub Kruzik 
6426dbc55e5SJakub Kruzik    Input Parameter:
6436dbc55e5SJakub Kruzik .  mat - the composite matrix
6446dbc55e5SJakub Kruzik 
6456dbc55e5SJakub Kruzik    Output Parameter:
6466dbc55e5SJakub Kruzik .  type - type of composite
6476dbc55e5SJakub Kruzik 
6486dbc55e5SJakub Kruzik    Level: advanced
6496dbc55e5SJakub Kruzik 
6506dbc55e5SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetType(), MATCOMPOSITE
6516dbc55e5SJakub Kruzik 
6526dbc55e5SJakub Kruzik @*/
6536dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type)
6546dbc55e5SJakub Kruzik {
6556dbc55e5SJakub Kruzik   PetscErrorCode ierr;
6566dbc55e5SJakub Kruzik 
6576dbc55e5SJakub Kruzik   PetscFunctionBegin;
6586dbc55e5SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6596dbc55e5SJakub Kruzik   PetscValidPointer(type,2);
6606dbc55e5SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type));CHKERRQ(ierr);
6616dbc55e5SJakub Kruzik   PetscFunctionReturn(0);
6626dbc55e5SJakub Kruzik }
6636dbc55e5SJakub Kruzik 
6643b35acafSJakub Kruzik static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat,MatStructure str)
6653b35acafSJakub Kruzik {
6663b35acafSJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
6673b35acafSJakub Kruzik 
6683b35acafSJakub Kruzik   PetscFunctionBegin;
6693b35acafSJakub Kruzik   b->structure = str;
6703b35acafSJakub Kruzik   PetscFunctionReturn(0);
6713b35acafSJakub Kruzik }
6723b35acafSJakub Kruzik 
6733b35acafSJakub Kruzik /*@
6743b35acafSJakub Kruzik    MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix.
6753b35acafSJakub Kruzik 
6763b35acafSJakub Kruzik    Not Collective
6773b35acafSJakub Kruzik 
6783b35acafSJakub Kruzik    Input Parameters:
6793b35acafSJakub Kruzik +  mat - the composite matrix
6803b35acafSJakub Kruzik -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN (default) or SUBSET_NONZERO_PATTERN
6813b35acafSJakub Kruzik 
6823b35acafSJakub Kruzik    Level: advanced
6833b35acafSJakub Kruzik 
6843b35acafSJakub Kruzik    Notes:
6853b35acafSJakub Kruzik     Information about the matrices structure is used in MatCompositeMerge() for additive composite matrix.
6863b35acafSJakub Kruzik 
6873b35acafSJakub Kruzik .seealso: MatAXPY(), MatCreateComposite(), MatCompositeMerge() MatCompositeGetMatStructure(), MATCOMPOSITE
6883b35acafSJakub Kruzik 
6893b35acafSJakub Kruzik @*/
6903b35acafSJakub Kruzik PetscErrorCode MatCompositeSetMatStructure(Mat mat,MatStructure str)
6913b35acafSJakub Kruzik {
6923b35acafSJakub Kruzik   PetscErrorCode ierr;
6933b35acafSJakub Kruzik 
6943b35acafSJakub Kruzik   PetscFunctionBegin;
6953b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6963b35acafSJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeSetMatStructure_C",(Mat,MatStructure),(mat,str));CHKERRQ(ierr);
6973b35acafSJakub Kruzik   PetscFunctionReturn(0);
6983b35acafSJakub Kruzik }
6993b35acafSJakub Kruzik 
7003b35acafSJakub Kruzik static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat,MatStructure *str)
7013b35acafSJakub Kruzik {
7023b35acafSJakub Kruzik   Mat_Composite  *b = (Mat_Composite*)mat->data;
7033b35acafSJakub Kruzik 
7043b35acafSJakub Kruzik   PetscFunctionBegin;
7053b35acafSJakub Kruzik   *str = b->structure;
7063b35acafSJakub Kruzik   PetscFunctionReturn(0);
7073b35acafSJakub Kruzik }
7083b35acafSJakub Kruzik 
7093b35acafSJakub Kruzik /*@
7103b35acafSJakub Kruzik    MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix.
7113b35acafSJakub Kruzik 
7123b35acafSJakub Kruzik    Not Collective
7133b35acafSJakub Kruzik 
7143b35acafSJakub Kruzik    Input Parameter:
7153b35acafSJakub Kruzik .  mat - the composite matrix
7163b35acafSJakub Kruzik 
7173b35acafSJakub Kruzik    Output Parameter:
7183b35acafSJakub Kruzik .  str - structure of the matrices
7193b35acafSJakub Kruzik 
7203b35acafSJakub Kruzik    Level: advanced
7213b35acafSJakub Kruzik 
7223b35acafSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetMatStructure(), MATCOMPOSITE
7233b35acafSJakub Kruzik 
7243b35acafSJakub Kruzik @*/
7253b35acafSJakub Kruzik PetscErrorCode MatCompositeGetMatStructure(Mat mat,MatStructure *str)
7263b35acafSJakub Kruzik {
7273b35acafSJakub Kruzik   PetscErrorCode ierr;
7283b35acafSJakub Kruzik 
7293b35acafSJakub Kruzik   PetscFunctionBegin;
7303b35acafSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7313b35acafSJakub Kruzik   PetscValidPointer(str,2);
7323b35acafSJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeGetMatStructure_C",(Mat,MatStructure*),(mat,str));CHKERRQ(ierr);
7333b35acafSJakub Kruzik   PetscFunctionReturn(0);
7343b35acafSJakub Kruzik }
7353b35acafSJakub Kruzik 
7368c8613bfSJakub Kruzik static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type)
7378c8613bfSJakub Kruzik {
7388c8613bfSJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
7398c8613bfSJakub Kruzik 
7408c8613bfSJakub Kruzik   PetscFunctionBegin;
7418c8613bfSJakub Kruzik   shell->mergetype = type;
7428c8613bfSJakub Kruzik   PetscFunctionReturn(0);
7438c8613bfSJakub Kruzik }
7448c8613bfSJakub Kruzik 
7458c8613bfSJakub Kruzik /*@
7468c8613bfSJakub Kruzik    MatCompositeSetMergeType - Sets order of MatCompositeMerge().
7478c8613bfSJakub Kruzik 
7488c8613bfSJakub Kruzik    Logically Collective on Mat
7498c8613bfSJakub Kruzik 
7508c8613bfSJakub Kruzik    Input Parameter:
7518c8613bfSJakub Kruzik +  mat - the composite matrix
7528c8613bfSJakub Kruzik -  type - MAT_COMPOSITE_MERGE RIGHT (default) to start merge from right with the first added matrix (mat[0]),
7538c8613bfSJakub Kruzik           MAT_COMPOSITE_MERGE_LEFT to start merge from left with the last added matrix (mat[nmat-1])
7548c8613bfSJakub Kruzik 
7558c8613bfSJakub Kruzik    Level: advanced
7568c8613bfSJakub Kruzik 
7578c8613bfSJakub Kruzik    Notes:
7588c8613bfSJakub Kruzik     The resulting matrix is the same regardles of the MergeType. Only the order of operation is changed.
7598c8613bfSJakub Kruzik     If set to MAT_COMPOSITE_MERGE_RIGHT the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0])))
7608c8613bfSJakub Kruzik     otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0].
7618c8613bfSJakub Kruzik 
7628c8613bfSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeMerge(), MATCOMPOSITE
7638c8613bfSJakub Kruzik 
7648c8613bfSJakub Kruzik @*/
7658c8613bfSJakub Kruzik PetscErrorCode MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type)
7668c8613bfSJakub Kruzik {
7678c8613bfSJakub Kruzik   PetscErrorCode ierr;
7688c8613bfSJakub Kruzik 
7698c8613bfSJakub Kruzik   PetscFunctionBegin;
7708c8613bfSJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7718c8613bfSJakub Kruzik   PetscValidLogicalCollectiveEnum(mat,type,2);
7728c8613bfSJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeSetMergeType_C",(Mat,MatCompositeMergeType),(mat,type));CHKERRQ(ierr);
7738c8613bfSJakub Kruzik   PetscFunctionReturn(0);
7748c8613bfSJakub Kruzik }
7758c8613bfSJakub Kruzik 
776d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat)
777b52f573bSBarry Smith {
778b52f573bSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)mat->data;
7796d7c1e57SBarry Smith   Mat_CompositeLink next   = shell->head, prev = shell->tail;
780b52f573bSBarry Smith   PetscErrorCode    ierr;
7816d7c1e57SBarry Smith   Mat               tmat,newmat;
7821ba60692SJed Brown   Vec               left,right;
7831ba60692SJed Brown   PetscScalar       scale;
78403049c21SJunchao Zhang   PetscInt          i;
785b52f573bSBarry Smith 
786b52f573bSBarry Smith   PetscFunctionBegin;
787e32f2f54SBarry Smith   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
788b52f573bSBarry Smith 
789b52f573bSBarry Smith   PetscFunctionBegin;
79003049c21SJunchao Zhang   scale = shell->scale;
7916d7c1e57SBarry Smith   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
7928c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
79303049c21SJunchao Zhang       i = 0;
794b52f573bSBarry Smith       ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
79503049c21SJunchao Zhang       if (shell->scalings) {ierr = MatScale(tmat,shell->scalings[i++]);CHKERRQ(ierr);}
796b52f573bSBarry Smith       while ((next = next->next)) {
79703049c21SJunchao Zhang         ierr = MatAXPY(tmat,(shell->scalings ? shell->scalings[i++] : 1.0),next->mat,shell->structure);CHKERRQ(ierr);
798b52f573bSBarry Smith       }
7996d7c1e57SBarry Smith     } else {
80003049c21SJunchao Zhang       i = shell->nmat-1;
8018c8613bfSJakub Kruzik       ierr = MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
80203049c21SJunchao Zhang       if (shell->scalings) {ierr = MatScale(tmat,shell->scalings[i--]);CHKERRQ(ierr);}
8038c8613bfSJakub Kruzik       while ((prev = prev->prev)) {
80403049c21SJunchao Zhang         ierr = MatAXPY(tmat,(shell->scalings ? shell->scalings[i--] : 1.0),prev->mat,shell->structure);CHKERRQ(ierr);
8058c8613bfSJakub Kruzik       }
8068c8613bfSJakub Kruzik     }
8078c8613bfSJakub Kruzik   } else {
8088c8613bfSJakub Kruzik     if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) {
8096d7c1e57SBarry Smith       ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
810b6cfcaf5SJakub Kruzik       while ((next = next->next)) {
811b6cfcaf5SJakub Kruzik         ierr = MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr);
8126bf464f9SBarry Smith         ierr = MatDestroy(&tmat);CHKERRQ(ierr);
8136d7c1e57SBarry Smith         tmat = newmat;
8146d7c1e57SBarry Smith       }
81504d534ceSJakub Kruzik     } else {
81604d534ceSJakub Kruzik       ierr = MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
81704d534ceSJakub Kruzik       while ((prev = prev->prev)) {
81804d534ceSJakub Kruzik         ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr);
81904d534ceSJakub Kruzik         ierr = MatDestroy(&tmat);CHKERRQ(ierr);
82004d534ceSJakub Kruzik         tmat = newmat;
82104d534ceSJakub Kruzik       }
82204d534ceSJakub Kruzik     }
82303049c21SJunchao Zhang     if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];}
8246d7c1e57SBarry Smith   }
8251ba60692SJed Brown 
8261ba60692SJed Brown   if ((left = shell->left)) {ierr = PetscObjectReference((PetscObject)left);CHKERRQ(ierr);}
8271ba60692SJed Brown   if ((right = shell->right)) {ierr = PetscObjectReference((PetscObject)right);CHKERRQ(ierr);}
8281ba60692SJed Brown 
82928be2f97SBarry Smith   ierr = MatHeaderReplace(mat,&tmat);CHKERRQ(ierr);
8301ba60692SJed Brown 
8311ba60692SJed Brown   ierr = MatDiagonalScale(mat,left,right);CHKERRQ(ierr);
8321ba60692SJed Brown   ierr = MatScale(mat,scale);CHKERRQ(ierr);
8331ba60692SJed Brown   ierr = VecDestroy(&left);CHKERRQ(ierr);
8341ba60692SJed Brown   ierr = VecDestroy(&right);CHKERRQ(ierr);
835b52f573bSBarry Smith   PetscFunctionReturn(0);
836b52f573bSBarry Smith }
8376a7ede75SJakub Kruzik 
8386a7ede75SJakub Kruzik /*@
839d7f81bf2SJakub Kruzik    MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
8408bd739bdSJakub Kruzik      by summing or computing the product of all the matrices inside the composite matrix.
841d7f81bf2SJakub Kruzik 
842b2b245abSJakub Kruzik   Collective
843d7f81bf2SJakub Kruzik 
844d7f81bf2SJakub Kruzik    Input Parameters:
845d7f81bf2SJakub Kruzik .  mat - the composite matrix
846d7f81bf2SJakub Kruzik 
847d7f81bf2SJakub Kruzik 
8484b2558d6SJakub Kruzik    Options Database Keys:
849b28d0daaSJakub Kruzik +  -mat_composite_merge - merge in MatAssemblyEnd()
850b28d0daaSJakub Kruzik -  -mat_composite_merge_type - set merge direction
851d7f81bf2SJakub Kruzik 
852d7f81bf2SJakub Kruzik    Level: advanced
853d7f81bf2SJakub Kruzik 
854d7f81bf2SJakub Kruzik    Notes:
855d7f81bf2SJakub Kruzik       The MatType of the resulting matrix will be the same as the MatType of the FIRST
856d7f81bf2SJakub Kruzik     matrix in the composite matrix.
857d7f81bf2SJakub Kruzik 
8583b35acafSJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeSetMatStructure(), MatCompositeSetMergeType(), MATCOMPOSITE
859d7f81bf2SJakub Kruzik 
860d7f81bf2SJakub Kruzik @*/
861d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat)
862d7f81bf2SJakub Kruzik {
863d7f81bf2SJakub Kruzik   PetscErrorCode ierr;
864d7f81bf2SJakub Kruzik 
865d7f81bf2SJakub Kruzik   PetscFunctionBegin;
866d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
867d7f81bf2SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));CHKERRQ(ierr);
868d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
869d7f81bf2SJakub Kruzik }
870d7f81bf2SJakub Kruzik 
8716d0add67SJakub Kruzik static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat,PetscInt *nmat)
872d7f81bf2SJakub Kruzik {
873d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
874d7f81bf2SJakub Kruzik 
875d7f81bf2SJakub Kruzik   PetscFunctionBegin;
876d7f81bf2SJakub Kruzik   *nmat = shell->nmat;
877d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
878d7f81bf2SJakub Kruzik }
879d7f81bf2SJakub Kruzik 
880d7f81bf2SJakub Kruzik /*@
8816d0add67SJakub Kruzik    MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix.
8826a7ede75SJakub Kruzik 
8836a7ede75SJakub Kruzik    Not Collective
8846a7ede75SJakub Kruzik 
8856a7ede75SJakub Kruzik    Input Parameter:
886d7f81bf2SJakub Kruzik .  mat - the composite matrix
8876a7ede75SJakub Kruzik 
8886a7ede75SJakub Kruzik    Output Parameter:
8896d0add67SJakub Kruzik .  nmat - number of matrices in the composite matrix
8906a7ede75SJakub Kruzik 
8918b5c3584SJakub Kruzik    Level: advanced
8926a7ede75SJakub Kruzik 
8938bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE
8946a7ede75SJakub Kruzik 
8956a7ede75SJakub Kruzik @*/
8966d0add67SJakub Kruzik PetscErrorCode MatCompositeGetNumberMat(Mat mat,PetscInt *nmat)
8976a7ede75SJakub Kruzik {
898d7f81bf2SJakub Kruzik   PetscErrorCode ierr;
8996a7ede75SJakub Kruzik 
9006a7ede75SJakub Kruzik   PetscFunctionBegin;
901d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9026a7ede75SJakub Kruzik   PetscValidPointer(nmat,2);
9036d0add67SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeGetNumberMat_C",(Mat,PetscInt*),(mat,nmat));CHKERRQ(ierr);
904d7f81bf2SJakub Kruzik   PetscFunctionReturn(0);
905d7f81bf2SJakub Kruzik }
906d7f81bf2SJakub Kruzik 
907d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai)
908d7f81bf2SJakub Kruzik {
909d7f81bf2SJakub Kruzik   Mat_Composite     *shell = (Mat_Composite*)mat->data;
910d7f81bf2SJakub Kruzik   Mat_CompositeLink ilink;
911d7f81bf2SJakub Kruzik   PetscInt          k;
912d7f81bf2SJakub Kruzik 
913d7f81bf2SJakub Kruzik   PetscFunctionBegin;
914d7f81bf2SJakub Kruzik   if (i >= shell->nmat) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %d >= %d",i,shell->nmat);
915d7f81bf2SJakub Kruzik   ilink = shell->head;
916d7f81bf2SJakub Kruzik   for (k=0; k<i; k++) {
917d7f81bf2SJakub Kruzik     ilink = ilink->next;
918d7f81bf2SJakub Kruzik   }
919d7f81bf2SJakub Kruzik   *Ai = ilink->mat;
9206a7ede75SJakub Kruzik   PetscFunctionReturn(0);
9216a7ede75SJakub Kruzik }
9226a7ede75SJakub Kruzik 
9238b5c3584SJakub Kruzik /*@
9248bd739bdSJakub Kruzik    MatCompositeGetMat - Returns the ith matrix from the composite matrix.
9258b5c3584SJakub Kruzik 
9268b5c3584SJakub Kruzik    Logically Collective on Mat
9278b5c3584SJakub Kruzik 
9288b5c3584SJakub Kruzik    Input Parameter:
929d7f81bf2SJakub Kruzik +  mat - the composite matrix
9308b5c3584SJakub Kruzik -  i - the number of requested matrix
9318b5c3584SJakub Kruzik 
9328b5c3584SJakub Kruzik    Output Parameter:
9338b5c3584SJakub Kruzik .  Ai - ith matrix in composite
9348b5c3584SJakub Kruzik 
9358b5c3584SJakub Kruzik    Level: advanced
9368b5c3584SJakub Kruzik 
9376d0add67SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNumberMat(), MatCompositeAddMat(), MATCOMPOSITE
9388b5c3584SJakub Kruzik 
9398b5c3584SJakub Kruzik @*/
940d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai)
9418b5c3584SJakub Kruzik {
942d7f81bf2SJakub Kruzik   PetscErrorCode ierr;
9438b5c3584SJakub Kruzik 
9448b5c3584SJakub Kruzik   PetscFunctionBegin;
945d7f81bf2SJakub Kruzik   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
946d7f81bf2SJakub Kruzik   PetscValidLogicalCollectiveInt(mat,i,2);
9478b5c3584SJakub Kruzik   PetscValidPointer(Ai,3);
948d7f81bf2SJakub Kruzik   ierr = PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));CHKERRQ(ierr);
9498b5c3584SJakub Kruzik   PetscFunctionReturn(0);
9508b5c3584SJakub Kruzik }
9518b5c3584SJakub Kruzik 
95203049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings_Composite(Mat mat,const PetscScalar *scalings)
95303049c21SJunchao Zhang {
95403049c21SJunchao Zhang   PetscErrorCode ierr;
95503049c21SJunchao Zhang   Mat_Composite  *shell = (Mat_Composite*)mat->data;
95603049c21SJunchao Zhang   PetscInt       nmat;
95703049c21SJunchao Zhang 
95803049c21SJunchao Zhang   PetscFunctionBegin;
95903049c21SJunchao Zhang   ierr = MatCompositeGetNumberMat(mat,&nmat);CHKERRQ(ierr);
96003049c21SJunchao Zhang   if (!shell->scalings) {ierr = PetscMalloc1(nmat,&shell->scalings);CHKERRQ(ierr);}
96103049c21SJunchao Zhang   ierr = PetscArraycpy(shell->scalings,scalings,nmat);CHKERRQ(ierr);
96203049c21SJunchao Zhang   PetscFunctionReturn(0);
96303049c21SJunchao Zhang }
96403049c21SJunchao Zhang 
96503049c21SJunchao Zhang /*@
96603049c21SJunchao Zhang    MatCompositeSetScalings - Sets separate scaling factors for component matrices.
96703049c21SJunchao Zhang 
96803049c21SJunchao Zhang    Logically Collective on Mat
96903049c21SJunchao Zhang 
97003049c21SJunchao Zhang    Input Parameter:
97103049c21SJunchao Zhang +  mat      - the composite matrix
97203049c21SJunchao Zhang -  scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat)
97303049c21SJunchao Zhang 
97403049c21SJunchao Zhang    Level: advanced
97503049c21SJunchao Zhang 
97603049c21SJunchao Zhang .seealso: MatScale(), MatDiagonalScale(), MATCOMPOSITE
97703049c21SJunchao Zhang 
97803049c21SJunchao Zhang @*/
97903049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings(Mat mat,const PetscScalar *scalings)
98003049c21SJunchao Zhang {
98103049c21SJunchao Zhang   PetscErrorCode ierr;
98203049c21SJunchao Zhang 
98303049c21SJunchao Zhang   PetscFunctionBegin;
98403049c21SJunchao Zhang   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
98503049c21SJunchao Zhang   PetscValidPointer(scalings,2);
98603049c21SJunchao Zhang   PetscValidLogicalCollectiveScalar(mat,*scalings,2);
98703049c21SJunchao Zhang   ierr = PetscUseMethod(mat,"MatCompositeSetScalings_C",(Mat,const PetscScalar*),(mat,scalings));CHKERRQ(ierr);
98803049c21SJunchao Zhang   PetscFunctionReturn(0);
98903049c21SJunchao Zhang }
99003049c21SJunchao Zhang 
991*f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL,
992*f4259b30SLisandro Dalcin                                        NULL,
993*f4259b30SLisandro Dalcin                                        NULL,
99441cd0125SJakub Kruzik                                        MatMult_Composite,
9957bf3a71bSJakub Kruzik                                        MatMultAdd_Composite,
99641cd0125SJakub Kruzik                                /*  5*/ MatMultTranspose_Composite,
9977bf3a71bSJakub Kruzik                                        MatMultTransposeAdd_Composite,
998*f4259b30SLisandro Dalcin                                        NULL,
999*f4259b30SLisandro Dalcin                                        NULL,
1000*f4259b30SLisandro Dalcin                                        NULL,
1001*f4259b30SLisandro Dalcin                                /* 10*/ NULL,
1002*f4259b30SLisandro Dalcin                                        NULL,
1003*f4259b30SLisandro Dalcin                                        NULL,
1004*f4259b30SLisandro Dalcin                                        NULL,
1005*f4259b30SLisandro Dalcin                                        NULL,
1006*f4259b30SLisandro Dalcin                                /* 15*/ NULL,
1007*f4259b30SLisandro Dalcin                                        NULL,
100841cd0125SJakub Kruzik                                        MatGetDiagonal_Composite,
100941cd0125SJakub Kruzik                                        MatDiagonalScale_Composite,
1010*f4259b30SLisandro Dalcin                                        NULL,
1011*f4259b30SLisandro Dalcin                                /* 20*/ NULL,
101241cd0125SJakub Kruzik                                        MatAssemblyEnd_Composite,
1013*f4259b30SLisandro Dalcin                                        NULL,
1014*f4259b30SLisandro Dalcin                                        NULL,
1015*f4259b30SLisandro Dalcin                                /* 24*/ NULL,
1016*f4259b30SLisandro Dalcin                                        NULL,
1017*f4259b30SLisandro Dalcin                                        NULL,
1018*f4259b30SLisandro Dalcin                                        NULL,
1019*f4259b30SLisandro Dalcin                                        NULL,
1020*f4259b30SLisandro Dalcin                                /* 29*/ NULL,
1021*f4259b30SLisandro Dalcin                                        NULL,
1022*f4259b30SLisandro Dalcin                                        NULL,
1023*f4259b30SLisandro Dalcin                                        NULL,
1024*f4259b30SLisandro Dalcin                                        NULL,
1025*f4259b30SLisandro Dalcin                                /* 34*/ NULL,
1026*f4259b30SLisandro Dalcin                                        NULL,
1027*f4259b30SLisandro Dalcin                                        NULL,
1028*f4259b30SLisandro Dalcin                                        NULL,
1029*f4259b30SLisandro Dalcin                                        NULL,
1030*f4259b30SLisandro Dalcin                                /* 39*/ NULL,
1031*f4259b30SLisandro Dalcin                                        NULL,
1032*f4259b30SLisandro Dalcin                                        NULL,
1033*f4259b30SLisandro Dalcin                                        NULL,
1034*f4259b30SLisandro Dalcin                                        NULL,
1035*f4259b30SLisandro Dalcin                                /* 44*/ NULL,
103641cd0125SJakub Kruzik                                        MatScale_Composite,
103741cd0125SJakub Kruzik                                        MatShift_Basic,
1038*f4259b30SLisandro Dalcin                                        NULL,
1039*f4259b30SLisandro Dalcin                                        NULL,
1040*f4259b30SLisandro Dalcin                                /* 49*/ NULL,
1041*f4259b30SLisandro Dalcin                                        NULL,
1042*f4259b30SLisandro Dalcin                                        NULL,
1043*f4259b30SLisandro Dalcin                                        NULL,
1044*f4259b30SLisandro Dalcin                                        NULL,
1045*f4259b30SLisandro Dalcin                                /* 54*/ NULL,
1046*f4259b30SLisandro Dalcin                                        NULL,
1047*f4259b30SLisandro Dalcin                                        NULL,
1048*f4259b30SLisandro Dalcin                                        NULL,
1049*f4259b30SLisandro Dalcin                                        NULL,
1050*f4259b30SLisandro Dalcin                                /* 59*/ NULL,
105141cd0125SJakub Kruzik                                        MatDestroy_Composite,
1052*f4259b30SLisandro Dalcin                                        NULL,
1053*f4259b30SLisandro Dalcin                                        NULL,
1054*f4259b30SLisandro Dalcin                                        NULL,
1055*f4259b30SLisandro Dalcin                                /* 64*/ NULL,
1056*f4259b30SLisandro Dalcin                                        NULL,
1057*f4259b30SLisandro Dalcin                                        NULL,
1058*f4259b30SLisandro Dalcin                                        NULL,
1059*f4259b30SLisandro Dalcin                                        NULL,
1060*f4259b30SLisandro Dalcin                                /* 69*/ NULL,
1061*f4259b30SLisandro Dalcin                                        NULL,
1062*f4259b30SLisandro Dalcin                                        NULL,
1063*f4259b30SLisandro Dalcin                                        NULL,
1064*f4259b30SLisandro Dalcin                                        NULL,
1065*f4259b30SLisandro Dalcin                                /* 74*/ NULL,
1066*f4259b30SLisandro Dalcin                                        NULL,
10674b2558d6SJakub Kruzik                                        MatSetFromOptions_Composite,
1068*f4259b30SLisandro Dalcin                                        NULL,
1069*f4259b30SLisandro Dalcin                                        NULL,
1070*f4259b30SLisandro Dalcin                                /* 79*/ NULL,
1071*f4259b30SLisandro Dalcin                                        NULL,
1072*f4259b30SLisandro Dalcin                                        NULL,
1073*f4259b30SLisandro Dalcin                                        NULL,
1074*f4259b30SLisandro Dalcin                                        NULL,
1075*f4259b30SLisandro Dalcin                                /* 84*/ NULL,
1076*f4259b30SLisandro Dalcin                                        NULL,
1077*f4259b30SLisandro Dalcin                                        NULL,
1078*f4259b30SLisandro Dalcin                                        NULL,
1079*f4259b30SLisandro Dalcin                                        NULL,
1080*f4259b30SLisandro Dalcin                                /* 89*/ NULL,
1081*f4259b30SLisandro Dalcin                                        NULL,
1082*f4259b30SLisandro Dalcin                                        NULL,
1083*f4259b30SLisandro Dalcin                                        NULL,
1084*f4259b30SLisandro Dalcin                                        NULL,
1085*f4259b30SLisandro Dalcin                                /* 94*/ NULL,
1086*f4259b30SLisandro Dalcin                                        NULL,
1087*f4259b30SLisandro Dalcin                                        NULL,
1088*f4259b30SLisandro Dalcin                                        NULL,
1089*f4259b30SLisandro Dalcin                                        NULL,
1090*f4259b30SLisandro Dalcin                                 /*99*/ NULL,
1091*f4259b30SLisandro Dalcin                                        NULL,
1092*f4259b30SLisandro Dalcin                                        NULL,
1093*f4259b30SLisandro Dalcin                                        NULL,
1094*f4259b30SLisandro Dalcin                                        NULL,
1095*f4259b30SLisandro Dalcin                                /*104*/ NULL,
1096*f4259b30SLisandro Dalcin                                        NULL,
1097*f4259b30SLisandro Dalcin                                        NULL,
1098*f4259b30SLisandro Dalcin                                        NULL,
1099*f4259b30SLisandro Dalcin                                        NULL,
1100*f4259b30SLisandro Dalcin                                /*109*/ NULL,
1101*f4259b30SLisandro Dalcin                                        NULL,
1102*f4259b30SLisandro Dalcin                                        NULL,
1103*f4259b30SLisandro Dalcin                                        NULL,
1104*f4259b30SLisandro Dalcin                                        NULL,
1105*f4259b30SLisandro Dalcin                                /*114*/ NULL,
1106*f4259b30SLisandro Dalcin                                        NULL,
1107*f4259b30SLisandro Dalcin                                        NULL,
1108*f4259b30SLisandro Dalcin                                        NULL,
1109*f4259b30SLisandro Dalcin                                        NULL,
1110*f4259b30SLisandro Dalcin                                /*119*/ NULL,
1111*f4259b30SLisandro Dalcin                                        NULL,
1112*f4259b30SLisandro Dalcin                                        NULL,
1113*f4259b30SLisandro Dalcin                                        NULL,
1114*f4259b30SLisandro Dalcin                                        NULL,
1115*f4259b30SLisandro Dalcin                                /*124*/ NULL,
1116*f4259b30SLisandro Dalcin                                        NULL,
1117*f4259b30SLisandro Dalcin                                        NULL,
1118*f4259b30SLisandro Dalcin                                        NULL,
1119*f4259b30SLisandro Dalcin                                        NULL,
1120*f4259b30SLisandro Dalcin                                /*129*/ NULL,
1121*f4259b30SLisandro Dalcin                                        NULL,
1122*f4259b30SLisandro Dalcin                                        NULL,
1123*f4259b30SLisandro Dalcin                                        NULL,
1124*f4259b30SLisandro Dalcin                                        NULL,
1125*f4259b30SLisandro Dalcin                                /*134*/ NULL,
1126*f4259b30SLisandro Dalcin                                        NULL,
1127*f4259b30SLisandro Dalcin                                        NULL,
1128*f4259b30SLisandro Dalcin                                        NULL,
1129*f4259b30SLisandro Dalcin                                        NULL,
1130*f4259b30SLisandro Dalcin                                /*139*/ NULL,
1131*f4259b30SLisandro Dalcin                                        NULL,
1132*f4259b30SLisandro Dalcin                                        NULL
113341cd0125SJakub Kruzik };
113441cd0125SJakub Kruzik 
113541cd0125SJakub Kruzik /*MC
113641cd0125SJakub Kruzik    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices.
113741cd0125SJakub Kruzik     The matrices need to have a correct size and parallel layout for the sum or product to be valid.
113841cd0125SJakub Kruzik 
113941cd0125SJakub Kruzik    Notes:
114041cd0125SJakub Kruzik     to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE);
114141cd0125SJakub Kruzik 
114241cd0125SJakub Kruzik   Level: advanced
114341cd0125SJakub Kruzik 
114403049c21SJunchao Zhang .seealso: MatCreateComposite(), MatCompositeSetScalings(), MatCompositeAddMat(), MatSetType(), MatCompositeSetType(), MatCompositeGetType(), MatCompositeSetMatStructure(), MatCompositeGetMatStructure(), MatCompositeMerge(), MatCompositeSetMergeType(), MatCompositeGetNumberMat(), MatCompositeGetMat()
114541cd0125SJakub Kruzik M*/
114641cd0125SJakub Kruzik 
114741cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A)
114841cd0125SJakub Kruzik {
114941cd0125SJakub Kruzik   Mat_Composite  *b;
115041cd0125SJakub Kruzik   PetscErrorCode ierr;
115141cd0125SJakub Kruzik 
115241cd0125SJakub Kruzik   PetscFunctionBegin;
115341cd0125SJakub Kruzik   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
115441cd0125SJakub Kruzik   A->data = (void*)b;
115541cd0125SJakub Kruzik   ierr    = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
115641cd0125SJakub Kruzik 
115741cd0125SJakub Kruzik   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
115841cd0125SJakub Kruzik   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
115941cd0125SJakub Kruzik 
116041cd0125SJakub Kruzik   A->assembled    = PETSC_TRUE;
116141cd0125SJakub Kruzik   A->preallocated = PETSC_TRUE;
116241cd0125SJakub Kruzik   b->type         = MAT_COMPOSITE_ADDITIVE;
116341cd0125SJakub Kruzik   b->scale        = 1.0;
116441cd0125SJakub Kruzik   b->nmat         = 0;
11654b2558d6SJakub Kruzik   b->merge        = PETSC_FALSE;
11668c8613bfSJakub Kruzik   b->mergetype    = MAT_COMPOSITE_MERGE_RIGHT;
11673b35acafSJakub Kruzik   b->structure    = DIFFERENT_NONZERO_PATTERN;
116803049c21SJunchao Zhang   b->merge_mvctx  = PETSC_TRUE;
116903049c21SJunchao Zhang 
117003049c21SJunchao Zhang 
117141cd0125SJakub Kruzik 
117241cd0125SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite);CHKERRQ(ierr);
117341cd0125SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite);CHKERRQ(ierr);
11746dbc55e5SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite);CHKERRQ(ierr);
11758c8613bfSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeType_C",MatCompositeSetMergeType_Composite);CHKERRQ(ierr);
11763b35acafSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMatStructure_C",MatCompositeSetMatStructure_Composite);CHKERRQ(ierr);
11773b35acafSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMatStructure_C",MatCompositeGetMatStructure_Composite);CHKERRQ(ierr);
117841cd0125SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite);CHKERRQ(ierr);
11796d0add67SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNumberMat_C",MatCompositeGetNumberMat_Composite);CHKERRQ(ierr);
118041cd0125SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite);CHKERRQ(ierr);
118103049c21SJunchao Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetScalings_C",MatCompositeSetScalings_Composite);CHKERRQ(ierr);
118241cd0125SJakub Kruzik 
118341cd0125SJakub Kruzik   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr);
118441cd0125SJakub Kruzik   PetscFunctionReturn(0);
118541cd0125SJakub Kruzik }
118641cd0125SJakub Kruzik 
1187