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)); 682e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeAddMat_C",NULL)); 692e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetType_C",NULL)); 702e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetType_C",NULL)); 712e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetMergeType_C",NULL)); 722e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeSetMatStructure_C",NULL)); 732e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetMatStructure_C",NULL)); 742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeMerge_C",NULL)); 752e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetNumberMat_C",NULL)); 762e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatCompositeGetMat_C",NULL)); 772e956fe4SStefano 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)); 1101baa6e33SBarry Smith if (shell->left) PetscCall(VecPointwiseMult(y,shell->left,y)); 11103049c21SJunchao Zhang scale = shell->scale; 11203049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 1139566063dSJacob Faibussowitsch PetscCall(VecScale(y,scale)); 1146d7c1e57SBarry Smith PetscFunctionReturn(0); 1156d7c1e57SBarry Smith } 1166d7c1e57SBarry Smith 1176d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 1186d7c1e57SBarry Smith { 1196d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 1206d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 1216d7c1e57SBarry Smith Vec in,out; 12203049c21SJunchao Zhang PetscScalar scale; 12303049c21SJunchao Zhang PetscInt i; 1246d7c1e57SBarry Smith 1256d7c1e57SBarry Smith PetscFunctionBegin; 12628b400f6SJacob Faibussowitsch PetscCheck(tail,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 1276d7c1e57SBarry Smith in = x; 128e4fc5df0SBarry Smith if (shell->left) { 129e4fc5df0SBarry Smith if (!shell->leftwork) { 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left,&shell->leftwork)); 131e4fc5df0SBarry Smith } 1329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in)); 133e4fc5df0SBarry Smith in = shell->leftwork; 134e4fc5df0SBarry Smith } 1356d7c1e57SBarry Smith while (tail->prev) { 1366d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1379566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(tail->mat,NULL,&tail->prev->work)); 1386d7c1e57SBarry Smith } 1396d7c1e57SBarry Smith out = tail->prev->work; 1409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat,in,out)); 1416d7c1e57SBarry Smith in = out; 1426d7c1e57SBarry Smith tail = tail->prev; 1436d7c1e57SBarry Smith } 1449566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat,in,y)); 1451baa6e33SBarry Smith if (shell->right) PetscCall(VecPointwiseMult(y,shell->right,y)); 14603049c21SJunchao Zhang 14703049c21SJunchao Zhang scale = shell->scale; 14803049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 1499566063dSJacob Faibussowitsch PetscCall(VecScale(y,scale)); 1506d7c1e57SBarry Smith PetscFunctionReturn(0); 1516d7c1e57SBarry Smith } 1526d7c1e57SBarry Smith 15303049c21SJunchao Zhang PetscErrorCode MatMult_Composite(Mat mat,Vec x,Vec y) 154793850ffSBarry Smith { 15503049c21SJunchao Zhang Mat_Composite *shell = (Mat_Composite*)mat->data; 15603049c21SJunchao Zhang Mat_CompositeLink cur = shell->head; 15703049c21SJunchao Zhang Vec in,y2,xin; 15803049c21SJunchao Zhang Mat A,B; 15903049c21SJunchao Zhang PetscInt i,j,k,n,nuniq,lo,hi,mid,*gindices,*buf,*tmp,tot; 16003049c21SJunchao Zhang const PetscScalar *vals; 16103049c21SJunchao Zhang const PetscInt *garray; 16203049c21SJunchao Zhang IS ix,iy; 16303049c21SJunchao Zhang PetscBool match; 164793850ffSBarry Smith 165793850ffSBarry Smith PetscFunctionBegin; 16628b400f6SJacob Faibussowitsch PetscCheck(cur,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 167e4fc5df0SBarry Smith in = x; 168e4fc5df0SBarry Smith if (shell->right) { 169e4fc5df0SBarry Smith if (!shell->rightwork) { 1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->right,&shell->rightwork)); 171793850ffSBarry Smith } 1729566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->rightwork,shell->right,in)); 173e4fc5df0SBarry Smith in = shell->rightwork; 174e4fc5df0SBarry Smith } 17503049c21SJunchao Zhang 17603049c21SJunchao Zhang /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time 17703049c21SJunchao Zhang we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and 17803049c21SJunchao Zhang it is legal to merge Mvctx, because all component matrices have the same size. 17903049c21SJunchao Zhang */ 18003049c21SJunchao Zhang if (shell->merge_mvctx && !shell->Mvctx) { 18103049c21SJunchao Zhang /* Currently only implemented for MATMPIAIJ */ 18203049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat,MATMPIAIJ,&match)); 18403049c21SJunchao Zhang if (!match) { 18503049c21SJunchao Zhang shell->merge_mvctx = PETSC_FALSE; 18603049c21SJunchao Zhang goto skip_merge_mvctx; 187e4fc5df0SBarry Smith } 188e4fc5df0SBarry Smith } 18903049c21SJunchao Zhang 19003049c21SJunchao Zhang /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */ 19103049c21SJunchao Zhang tot = 0; 19203049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 1939566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,NULL)); 1949566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 19503049c21SJunchao Zhang tot += n; 19603049c21SJunchao Zhang } 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tot,&shell->location,tot,&shell->larray,shell->nmat,&shell->lvecs)); 19803049c21SJunchao Zhang shell->len = tot; 19903049c21SJunchao Zhang 20003049c21SJunchao Zhang /* Go through matrices second time to sort off-diag columns and remove dups */ 2019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot,&gindices)); /* No Malloc2() since we will give one to petsc and free the other */ 2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot,&buf)); 20303049c21SJunchao Zhang nuniq = 0; /* Number of unique nonzero columns */ 20403049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 2059566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray)); 2069566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 20703049c21SJunchao Zhang /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */ 20803049c21SJunchao Zhang i = j = k = 0; 20903049c21SJunchao Zhang while (i < n && j < nuniq) { 21003049c21SJunchao Zhang if (garray[i] < gindices[j]) buf[k++] = garray[i++]; 21103049c21SJunchao Zhang else if (garray[i] > gindices[j]) buf[k++] = gindices[j++]; 21203049c21SJunchao Zhang else {buf[k++] = garray[i++]; j++;} 21303049c21SJunchao Zhang } 21403049c21SJunchao Zhang /* Copy leftover in garray[] or gindices[] */ 21503049c21SJunchao Zhang if (i < n) { 2169566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf+k,garray+i,n-i)); 21703049c21SJunchao Zhang nuniq = k + n-i; 21803049c21SJunchao Zhang } else if (j < nuniq) { 2199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf+k,gindices+j,nuniq-j)); 22003049c21SJunchao Zhang nuniq = k + nuniq-j; 22103049c21SJunchao Zhang } else nuniq = k; 22203049c21SJunchao Zhang /* Swap gindices and buf to merge garray of the next matrix */ 22303049c21SJunchao Zhang tmp = gindices; 22403049c21SJunchao Zhang gindices = buf; 22503049c21SJunchao Zhang buf = tmp; 22603049c21SJunchao Zhang } 2279566063dSJacob Faibussowitsch PetscCall(PetscFree(buf)); 22803049c21SJunchao Zhang 22903049c21SJunchao Zhang /* Go through matrices third time to build a map from gindices[] to garray[] */ 23003049c21SJunchao Zhang tot = 0; 23103049c21SJunchao Zhang for (cur=shell->head,j=0; cur; cur=cur->next,j++) { /* j-th matrix */ 2329566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray)); 2339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 2349566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&shell->lvecs[j])); 23503049c21SJunchao Zhang /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */ 23603049c21SJunchao Zhang lo = 0; 23703049c21SJunchao Zhang for (i=0; i<n; i++) { 23803049c21SJunchao Zhang hi = nuniq; 23903049c21SJunchao Zhang while (hi - lo > 1) { 24003049c21SJunchao Zhang mid = lo + (hi - lo)/2; 24103049c21SJunchao Zhang if (garray[i] < gindices[mid]) hi = mid; 24203049c21SJunchao Zhang else lo = mid; 24303049c21SJunchao Zhang } 24403049c21SJunchao Zhang shell->location[tot+i] = lo; /* gindices[lo] = garray[i] */ 24503049c21SJunchao Zhang lo++; /* Since garray[i+1] > garray[i], we can safely advance lo */ 24603049c21SJunchao Zhang } 24703049c21SJunchao Zhang tot += n; 24803049c21SJunchao Zhang } 24903049c21SJunchao Zhang 25003049c21SJunchao Zhang /* Build merged Mvctx */ 2519566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nuniq,gindices,PETSC_OWN_POINTER,&ix)); 2529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,nuniq,0,1,&iy)); 2539566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&xin)); 2549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,nuniq,&shell->gvec)); 2559566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xin,ix,shell->gvec,iy,&shell->Mvctx)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xin)); 2579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ix)); 2589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iy)); 25903049c21SJunchao Zhang } 26003049c21SJunchao Zhang 26103049c21SJunchao Zhang skip_merge_mvctx: 2629566063dSJacob Faibussowitsch PetscCall(VecSet(y,0)); 2639566063dSJacob Faibussowitsch if (!shell->leftwork2) PetscCall(VecDuplicate(y,&shell->leftwork2)); 26403049c21SJunchao Zhang y2 = shell->leftwork2; 26503049c21SJunchao Zhang 26603049c21SJunchao Zhang if (shell->Mvctx) { /* Have a merged Mvctx */ 26703049c21SJunchao 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 26803049c21SJunchao 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 26903049c21SJunchao Zhang overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach. 27003049c21SJunchao Zhang */ 2719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD)); 2729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD)); 27303049c21SJunchao Zhang 2749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(shell->gvec,&vals)); 27503049c21SJunchao Zhang for (i=0; i<shell->len; i++) shell->larray[i] = vals[shell->location[i]]; 2769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(shell->gvec,&vals)); 27703049c21SJunchao Zhang 27803049c21SJunchao Zhang for (cur=shell->head,tot=i=0; cur; cur=cur->next,i++) { /* i-th matrix */ 2799566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,&A,&B,NULL)); 280*dbbe0bcdSBarry Smith PetscUseTypeMethod(A,mult ,in,y2); 2819566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 2829566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(shell->lvecs[i],&shell->larray[tot])); 2839566063dSJacob Faibussowitsch PetscCall((*B->ops->multadd)(B,shell->lvecs[i],y2,y2)); 2849566063dSJacob Faibussowitsch PetscCall(VecResetArray(shell->lvecs[i])); 2859566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,(shell->scalings ? shell->scalings[i] : 1.0),y2)); 28603049c21SJunchao Zhang tot += n; 28703049c21SJunchao Zhang } 28803049c21SJunchao Zhang } else { 28903049c21SJunchao Zhang if (shell->scalings) { 29003049c21SJunchao Zhang for (cur=shell->head,i=0; cur; cur=cur->next,i++) { 2919566063dSJacob Faibussowitsch PetscCall(MatMult(cur->mat,in,y2)); 2929566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->scalings[i],y2)); 29303049c21SJunchao Zhang } 29403049c21SJunchao Zhang } else { 2959566063dSJacob Faibussowitsch for (cur=shell->head; cur; cur=cur->next) PetscCall(MatMultAdd(cur->mat,in,y,y)); 29603049c21SJunchao Zhang } 29703049c21SJunchao Zhang } 29803049c21SJunchao Zhang 2999566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(y,shell->left,y)); 3009566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scale)); 301793850ffSBarry Smith PetscFunctionReturn(0); 302793850ffSBarry Smith } 303793850ffSBarry Smith 304793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 305793850ffSBarry Smith { 306793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 307793850ffSBarry Smith Mat_CompositeLink next = shell->head; 30803049c21SJunchao Zhang Vec in,y2 = NULL; 30903049c21SJunchao Zhang PetscInt i; 310793850ffSBarry Smith 311793850ffSBarry Smith PetscFunctionBegin; 31228b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 313e4fc5df0SBarry Smith in = x; 314e4fc5df0SBarry Smith if (shell->left) { 315e4fc5df0SBarry Smith if (!shell->leftwork) { 3169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left,&shell->leftwork)); 317793850ffSBarry Smith } 3189566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in)); 319e4fc5df0SBarry Smith in = shell->leftwork; 320e4fc5df0SBarry Smith } 32103049c21SJunchao Zhang 3229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat,in,y)); 32303049c21SJunchao Zhang if (shell->scalings) { 3249566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scalings[0])); 3259566063dSJacob Faibussowitsch if (!shell->rightwork2) PetscCall(VecDuplicate(y,&shell->rightwork2)); 32603049c21SJunchao Zhang y2 = shell->rightwork2; 32703049c21SJunchao Zhang } 32803049c21SJunchao Zhang i = 1; 329e4fc5df0SBarry Smith while ((next = next->next)) { 3309566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat,in,y,y)); 33103049c21SJunchao Zhang else { 3329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat,in,y2)); 3339566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->scalings[i++],y2)); 33403049c21SJunchao Zhang } 335e4fc5df0SBarry Smith } 3361baa6e33SBarry Smith if (shell->right) PetscCall(VecPointwiseMult(y,shell->right,y)); 3379566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scale)); 338793850ffSBarry Smith PetscFunctionReturn(0); 339793850ffSBarry Smith } 340793850ffSBarry Smith 3417bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z) 3427bf3a71bSJakub Kruzik { 3437bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 3447bf3a71bSJakub Kruzik 3457bf3a71bSJakub Kruzik PetscFunctionBegin; 3467bf3a71bSJakub Kruzik if (y != z) { 3479566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,z)); 3489566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 3497bf3a71bSJakub Kruzik } else { 3507bf3a71bSJakub Kruzik if (!shell->leftwork) { 3519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(z,&shell->leftwork)); 3527bf3a71bSJakub Kruzik } 3539566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,shell->leftwork)); 3549566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 3559566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->leftwork)); 3567bf3a71bSJakub Kruzik } 3577bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3587bf3a71bSJakub Kruzik } 3597bf3a71bSJakub Kruzik 3607bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z) 3617bf3a71bSJakub Kruzik { 3627bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 3637bf3a71bSJakub Kruzik 3647bf3a71bSJakub Kruzik PetscFunctionBegin; 3657bf3a71bSJakub Kruzik if (y != z) { 3669566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,z)); 3679566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 3687bf3a71bSJakub Kruzik } else { 3697bf3a71bSJakub Kruzik if (!shell->rightwork) { 3709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(z,&shell->rightwork)); 3717bf3a71bSJakub Kruzik } 3729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,shell->rightwork)); 3739566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 3749566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->rightwork)); 3757bf3a71bSJakub Kruzik } 3767bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3777bf3a71bSJakub Kruzik } 3787bf3a71bSJakub Kruzik 379793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 380793850ffSBarry Smith { 381793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 382793850ffSBarry Smith Mat_CompositeLink next = shell->head; 38303049c21SJunchao Zhang PetscInt i; 384793850ffSBarry Smith 385793850ffSBarry Smith PetscFunctionBegin; 38628b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 38708401ef6SPierre Jolivet PetscCheck(!shell->right && !shell->left,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 388e4fc5df0SBarry Smith 3899566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat,v)); 3909566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(VecScale(v,shell->scalings[0])); 39103049c21SJunchao Zhang 392793850ffSBarry Smith if (next->next && !shell->work) { 3939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v,&shell->work)); 394793850ffSBarry Smith } 39503049c21SJunchao Zhang i = 1; 396793850ffSBarry Smith while ((next = next->next)) { 3979566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat,shell->work)); 3989566063dSJacob Faibussowitsch PetscCall(VecAXPY(v,(shell->scalings ? shell->scalings[i++] : 1.0),shell->work)); 399793850ffSBarry Smith } 4009566063dSJacob Faibussowitsch PetscCall(VecScale(v,shell->scale)); 401793850ffSBarry Smith PetscFunctionReturn(0); 402793850ffSBarry Smith } 403793850ffSBarry Smith 404793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 405793850ffSBarry Smith { 4064b2558d6SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)Y->data; 407b52f573bSBarry Smith 408793850ffSBarry Smith PetscFunctionBegin; 4091baa6e33SBarry Smith if (shell->merge) PetscCall(MatCompositeMerge(Y)); 410793850ffSBarry Smith PetscFunctionReturn(0); 411793850ffSBarry Smith } 412793850ffSBarry Smith 413e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 414e4fc5df0SBarry Smith { 415e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 416e4fc5df0SBarry Smith 417e4fc5df0SBarry Smith PetscFunctionBegin; 418321429dbSBarry Smith a->scale *= alpha; 419e4fc5df0SBarry Smith PetscFunctionReturn(0); 420e4fc5df0SBarry Smith } 421e4fc5df0SBarry Smith 422e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 423e4fc5df0SBarry Smith { 424e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 425e4fc5df0SBarry Smith 426e4fc5df0SBarry Smith PetscFunctionBegin; 427e4fc5df0SBarry Smith if (left) { 428321429dbSBarry Smith if (!a->left) { 4299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&a->left)); 4309566063dSJacob Faibussowitsch PetscCall(VecCopy(left,a->left)); 431321429dbSBarry Smith } else { 4329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left,left,a->left)); 433321429dbSBarry Smith } 434e4fc5df0SBarry Smith } 435e4fc5df0SBarry Smith if (right) { 436321429dbSBarry Smith if (!a->right) { 4379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&a->right)); 4389566063dSJacob Faibussowitsch PetscCall(VecCopy(right,a->right)); 439321429dbSBarry Smith } else { 4409566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right,right,a->right)); 441321429dbSBarry Smith } 442e4fc5df0SBarry Smith } 443e4fc5df0SBarry Smith PetscFunctionReturn(0); 444e4fc5df0SBarry Smith } 445793850ffSBarry Smith 446*dbbe0bcdSBarry Smith PetscErrorCode MatSetFromOptions_Composite(Mat A,PetscOptionItems *PetscOptionsObject) 4474b2558d6SJakub Kruzik { 4484b2558d6SJakub Kruzik Mat_Composite *a = (Mat_Composite*)A->data; 4494b2558d6SJakub Kruzik 4504b2558d6SJakub Kruzik PetscFunctionBegin; 451d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"MATCOMPOSITE options"); 4529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge","Merge at MatAssemblyEnd","MatCompositeMerge",a->merge,&a->merge,NULL)); 4539566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_composite_merge_type","Set composite merge direction","MatCompositeSetMergeType",MatCompositeMergeTypes,(PetscEnum)a->mergetype,(PetscEnum*)&a->mergetype,NULL)); 4549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx","Merge MatMult() vecscat contexts","MatCreateComposite",a->merge_mvctx,&a->merge_mvctx,NULL)); 455d0609cedSBarry Smith PetscOptionsHeadEnd(); 4564b2558d6SJakub Kruzik PetscFunctionReturn(0); 4574b2558d6SJakub Kruzik } 4584b2558d6SJakub Kruzik 4592c0821f3SBarry Smith /*@ 4608bd739bdSJakub Kruzik MatCreateComposite - Creates a matrix as the sum or product of one or more matrices 461793850ffSBarry Smith 462d083f849SBarry Smith Collective 463793850ffSBarry Smith 464793850ffSBarry Smith Input Parameters: 465793850ffSBarry Smith + comm - MPI communicator 466793850ffSBarry Smith . nmat - number of matrices to put in 467793850ffSBarry Smith - mats - the matrices 468793850ffSBarry Smith 469793850ffSBarry Smith Output Parameter: 470db36af27SMatthew G. Knepley . mat - the matrix 471793850ffSBarry Smith 4724b2558d6SJakub Kruzik Options Database Keys: 4734b2558d6SJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 47403049c21SJunchao Zhang . -mat_composite_merge_mvctx - merge Mvctx of component matrices to optimize communication in MatMult() for ADDITIVE matrices 475b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 4764b2558d6SJakub Kruzik 477793850ffSBarry Smith Level: advanced 478793850ffSBarry Smith 479793850ffSBarry Smith Notes: 480793850ffSBarry Smith Alternative construction 481793850ffSBarry Smith $ MatCreate(comm,&mat); 482793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 483793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 484793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 485793850ffSBarry Smith $ .... 486793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 487b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 488b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 489793850ffSBarry Smith 4906d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4916d7c1e57SBarry Smith 492db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCompositeGetMat()`, `MatCompositeMerge()`, `MatCompositeSetType()`, `MATCOMPOSITE` 493793850ffSBarry Smith 494793850ffSBarry Smith @*/ 4957087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 496793850ffSBarry Smith { 497793850ffSBarry Smith PetscInt m,n,M,N,i; 498793850ffSBarry Smith 499793850ffSBarry Smith PetscFunctionBegin; 50008401ef6SPierre Jolivet PetscCheck(nmat >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 501064a246eSJacob Faibussowitsch PetscValidPointer(mat,4); 502793850ffSBarry Smith 5039566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[0],PETSC_IGNORE,&n)); 5049566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE)); 5059566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[0],PETSC_IGNORE,&N)); 5069566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[nmat-1],&M,PETSC_IGNORE)); 5079566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,mat)); 5089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*mat,m,n,M,N)); 5099566063dSJacob Faibussowitsch PetscCall(MatSetType(*mat,MATCOMPOSITE)); 510793850ffSBarry Smith for (i=0; i<nmat; i++) { 5119566063dSJacob Faibussowitsch PetscCall(MatCompositeAddMat(*mat,mats[i])); 512793850ffSBarry Smith } 5139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY)); 5149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY)); 515793850ffSBarry Smith PetscFunctionReturn(0); 516793850ffSBarry Smith } 517793850ffSBarry Smith 518d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 519d7f81bf2SJakub Kruzik { 520d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 521d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 522d7f81bf2SJakub Kruzik 523d7f81bf2SJakub Kruzik PetscFunctionBegin; 5249566063dSJacob Faibussowitsch PetscCall(PetscNewLog(mat,&ilink)); 525f4259b30SLisandro Dalcin ilink->next = NULL; 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)smat)); 527d7f81bf2SJakub Kruzik ilink->mat = smat; 528d7f81bf2SJakub Kruzik 529d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 530d7f81bf2SJakub Kruzik else { 531d7f81bf2SJakub Kruzik while (next->next) { 532d7f81bf2SJakub Kruzik next = next->next; 533d7f81bf2SJakub Kruzik } 534d7f81bf2SJakub Kruzik next->next = ilink; 535d7f81bf2SJakub Kruzik ilink->prev = next; 536d7f81bf2SJakub Kruzik } 537d7f81bf2SJakub Kruzik shell->tail = ilink; 538d7f81bf2SJakub Kruzik shell->nmat += 1; 53903049c21SJunchao Zhang 54003049c21SJunchao Zhang /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */ 54103049c21SJunchao Zhang if (shell->scalings) { 5429566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(PetscScalar)*shell->nmat,&shell->scalings)); 54303049c21SJunchao Zhang shell->scalings[shell->nmat-1] = 1.0; 54403049c21SJunchao Zhang } 545d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 546d7f81bf2SJakub Kruzik } 547d7f81bf2SJakub Kruzik 548793850ffSBarry Smith /*@ 5498bd739bdSJakub Kruzik MatCompositeAddMat - Add another matrix to a composite matrix. 550793850ffSBarry Smith 551793850ffSBarry Smith Collective on Mat 552793850ffSBarry Smith 553793850ffSBarry Smith Input Parameters: 554793850ffSBarry Smith + mat - the composite matrix 555793850ffSBarry Smith - smat - the partial matrix 556793850ffSBarry Smith 557793850ffSBarry Smith Level: advanced 558793850ffSBarry Smith 559db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE` 560793850ffSBarry Smith @*/ 5617087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 562793850ffSBarry Smith { 563793850ffSBarry Smith PetscFunctionBegin; 5640700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5650700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 566cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat)); 567d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 568d7f81bf2SJakub Kruzik } 569793850ffSBarry Smith 570d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 571d7f81bf2SJakub Kruzik { 572d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 573d7f81bf2SJakub Kruzik 574d7f81bf2SJakub Kruzik PetscFunctionBegin; 575ced1ac25SJakub Kruzik b->type = type; 576d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 577f4259b30SLisandro Dalcin mat->ops->getdiagonal = NULL; 578d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 579d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 58003049c21SJunchao Zhang b->merge_mvctx = PETSC_FALSE; 581d7f81bf2SJakub Kruzik } else { 582d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 583d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 584d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 585793850ffSBarry Smith } 5866d7c1e57SBarry Smith PetscFunctionReturn(0); 5876d7c1e57SBarry Smith } 5886d7c1e57SBarry Smith 5892c0821f3SBarry Smith /*@ 5908bd739bdSJakub Kruzik MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product. 5916d7c1e57SBarry Smith 592b2b245abSJakub Kruzik Logically Collective on Mat 5936d7c1e57SBarry Smith 5946d7c1e57SBarry Smith Input Parameters: 5956d7c1e57SBarry Smith . mat - the composite matrix 5966d7c1e57SBarry Smith 5976d7c1e57SBarry Smith Level: advanced 5986d7c1e57SBarry Smith 599db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE` 6006d7c1e57SBarry Smith 6016d7c1e57SBarry Smith @*/ 6027087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 6036d7c1e57SBarry Smith { 6046d7c1e57SBarry Smith PetscFunctionBegin; 605d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 606b2b245abSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 607cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type)); 608793850ffSBarry Smith PetscFunctionReturn(0); 609793850ffSBarry Smith } 610793850ffSBarry Smith 6116dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type) 6126dbc55e5SJakub Kruzik { 6136dbc55e5SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6146dbc55e5SJakub Kruzik 6156dbc55e5SJakub Kruzik PetscFunctionBegin; 6166dbc55e5SJakub Kruzik *type = b->type; 6176dbc55e5SJakub Kruzik PetscFunctionReturn(0); 6186dbc55e5SJakub Kruzik } 6196dbc55e5SJakub Kruzik 6206dbc55e5SJakub Kruzik /*@ 6216dbc55e5SJakub Kruzik MatCompositeGetType - Returns type of composite. 6226dbc55e5SJakub Kruzik 6236dbc55e5SJakub Kruzik Not Collective 6246dbc55e5SJakub Kruzik 6256dbc55e5SJakub Kruzik Input Parameter: 6266dbc55e5SJakub Kruzik . mat - the composite matrix 6276dbc55e5SJakub Kruzik 6286dbc55e5SJakub Kruzik Output Parameter: 6296dbc55e5SJakub Kruzik . type - type of composite 6306dbc55e5SJakub Kruzik 6316dbc55e5SJakub Kruzik Level: advanced 6326dbc55e5SJakub Kruzik 633db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE` 6346dbc55e5SJakub Kruzik 6356dbc55e5SJakub Kruzik @*/ 6366dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type) 6376dbc55e5SJakub Kruzik { 6386dbc55e5SJakub Kruzik PetscFunctionBegin; 6396dbc55e5SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6406dbc55e5SJakub Kruzik PetscValidPointer(type,2); 641cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type)); 6426dbc55e5SJakub Kruzik PetscFunctionReturn(0); 6436dbc55e5SJakub Kruzik } 6446dbc55e5SJakub Kruzik 6453b35acafSJakub Kruzik static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat,MatStructure str) 6463b35acafSJakub Kruzik { 6473b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6483b35acafSJakub Kruzik 6493b35acafSJakub Kruzik PetscFunctionBegin; 6503b35acafSJakub Kruzik b->structure = str; 6513b35acafSJakub Kruzik PetscFunctionReturn(0); 6523b35acafSJakub Kruzik } 6533b35acafSJakub Kruzik 6543b35acafSJakub Kruzik /*@ 6553b35acafSJakub Kruzik MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix. 6563b35acafSJakub Kruzik 6573b35acafSJakub Kruzik Not Collective 6583b35acafSJakub Kruzik 6593b35acafSJakub Kruzik Input Parameters: 6603b35acafSJakub Kruzik + mat - the composite matrix 6613b35acafSJakub Kruzik - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN (default) or SUBSET_NONZERO_PATTERN 6623b35acafSJakub Kruzik 6633b35acafSJakub Kruzik Level: advanced 6643b35acafSJakub Kruzik 6653b35acafSJakub Kruzik Notes: 6663b35acafSJakub Kruzik Information about the matrices structure is used in MatCompositeMerge() for additive composite matrix. 6673b35acafSJakub Kruzik 668db781477SPatrick Sanan .seealso: `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE` 6693b35acafSJakub Kruzik 6703b35acafSJakub Kruzik @*/ 6713b35acafSJakub Kruzik PetscErrorCode MatCompositeSetMatStructure(Mat mat,MatStructure str) 6723b35acafSJakub Kruzik { 6733b35acafSJakub Kruzik PetscFunctionBegin; 6743b35acafSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 675cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetMatStructure_C",(Mat,MatStructure),(mat,str)); 6763b35acafSJakub Kruzik PetscFunctionReturn(0); 6773b35acafSJakub Kruzik } 6783b35acafSJakub Kruzik 6793b35acafSJakub Kruzik static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat,MatStructure *str) 6803b35acafSJakub Kruzik { 6813b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6823b35acafSJakub Kruzik 6833b35acafSJakub Kruzik PetscFunctionBegin; 6843b35acafSJakub Kruzik *str = b->structure; 6853b35acafSJakub Kruzik PetscFunctionReturn(0); 6863b35acafSJakub Kruzik } 6873b35acafSJakub Kruzik 6883b35acafSJakub Kruzik /*@ 6893b35acafSJakub Kruzik MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix. 6903b35acafSJakub Kruzik 6913b35acafSJakub Kruzik Not Collective 6923b35acafSJakub Kruzik 6933b35acafSJakub Kruzik Input Parameter: 6943b35acafSJakub Kruzik . mat - the composite matrix 6953b35acafSJakub Kruzik 6963b35acafSJakub Kruzik Output Parameter: 6973b35acafSJakub Kruzik . str - structure of the matrices 6983b35acafSJakub Kruzik 6993b35acafSJakub Kruzik Level: advanced 7003b35acafSJakub Kruzik 701db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE` 7023b35acafSJakub Kruzik 7033b35acafSJakub Kruzik @*/ 7043b35acafSJakub Kruzik PetscErrorCode MatCompositeGetMatStructure(Mat mat,MatStructure *str) 7053b35acafSJakub Kruzik { 7063b35acafSJakub Kruzik PetscFunctionBegin; 7073b35acafSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7083b35acafSJakub Kruzik PetscValidPointer(str,2); 709cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetMatStructure_C",(Mat,MatStructure*),(mat,str)); 7103b35acafSJakub Kruzik PetscFunctionReturn(0); 7113b35acafSJakub Kruzik } 7123b35acafSJakub Kruzik 7138c8613bfSJakub Kruzik static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type) 7148c8613bfSJakub Kruzik { 7158c8613bfSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 7168c8613bfSJakub Kruzik 7178c8613bfSJakub Kruzik PetscFunctionBegin; 7188c8613bfSJakub Kruzik shell->mergetype = type; 7198c8613bfSJakub Kruzik PetscFunctionReturn(0); 7208c8613bfSJakub Kruzik } 7218c8613bfSJakub Kruzik 7228c8613bfSJakub Kruzik /*@ 7238c8613bfSJakub Kruzik MatCompositeSetMergeType - Sets order of MatCompositeMerge(). 7248c8613bfSJakub Kruzik 7258c8613bfSJakub Kruzik Logically Collective on Mat 7268c8613bfSJakub Kruzik 727d8d19677SJose E. Roman Input Parameters: 7288c8613bfSJakub Kruzik + mat - the composite matrix 7298c8613bfSJakub Kruzik - type - MAT_COMPOSITE_MERGE RIGHT (default) to start merge from right with the first added matrix (mat[0]), 7308c8613bfSJakub Kruzik MAT_COMPOSITE_MERGE_LEFT to start merge from left with the last added matrix (mat[nmat-1]) 7318c8613bfSJakub Kruzik 7328c8613bfSJakub Kruzik Level: advanced 7338c8613bfSJakub Kruzik 7348c8613bfSJakub Kruzik Notes: 7358c8613bfSJakub Kruzik The resulting matrix is the same regardles of the MergeType. Only the order of operation is changed. 7368c8613bfSJakub Kruzik If set to MAT_COMPOSITE_MERGE_RIGHT the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0]))) 7378c8613bfSJakub Kruzik otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0]. 7388c8613bfSJakub Kruzik 739db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE` 7408c8613bfSJakub Kruzik 7418c8613bfSJakub Kruzik @*/ 7428c8613bfSJakub Kruzik PetscErrorCode MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type) 7438c8613bfSJakub Kruzik { 7448c8613bfSJakub Kruzik PetscFunctionBegin; 7458c8613bfSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7468c8613bfSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 747cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetMergeType_C",(Mat,MatCompositeMergeType),(mat,type)); 7488c8613bfSJakub Kruzik PetscFunctionReturn(0); 7498c8613bfSJakub Kruzik } 7508c8613bfSJakub Kruzik 751d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 752b52f573bSBarry Smith { 753b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 7546d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 7556d7c1e57SBarry Smith Mat tmat,newmat; 7561ba60692SJed Brown Vec left,right; 7571ba60692SJed Brown PetscScalar scale; 75803049c21SJunchao Zhang PetscInt i; 759b52f573bSBarry Smith 760b52f573bSBarry Smith PetscFunctionBegin; 76128b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 76203049c21SJunchao Zhang scale = shell->scale; 7636d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 7648c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 76503049c21SJunchao Zhang i = 0; 7669566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat)); 7679566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i++])); 768b52f573bSBarry Smith while ((next = next->next)) { 7699566063dSJacob Faibussowitsch PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i++] : 1.0),next->mat,shell->structure)); 770b52f573bSBarry Smith } 7716d7c1e57SBarry Smith } else { 77203049c21SJunchao Zhang i = shell->nmat-1; 7739566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat)); 7749566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i--])); 7758c8613bfSJakub Kruzik while ((prev = prev->prev)) { 7769566063dSJacob Faibussowitsch PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i--] : 1.0),prev->mat,shell->structure)); 7778c8613bfSJakub Kruzik } 7788c8613bfSJakub Kruzik } 7798c8613bfSJakub Kruzik } else { 7808c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 7819566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat)); 782b6cfcaf5SJakub Kruzik while ((next = next->next)) { 7839566063dSJacob Faibussowitsch PetscCall(MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat)); 7849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 7856d7c1e57SBarry Smith tmat = newmat; 7866d7c1e57SBarry Smith } 78704d534ceSJakub Kruzik } else { 7889566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat)); 78904d534ceSJakub Kruzik while ((prev = prev->prev)) { 7909566063dSJacob Faibussowitsch PetscCall(MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat)); 7919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 79204d534ceSJakub Kruzik tmat = newmat; 79304d534ceSJakub Kruzik } 79404d534ceSJakub Kruzik } 79503049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 7966d7c1e57SBarry Smith } 7971ba60692SJed Brown 7989566063dSJacob Faibussowitsch if ((left = shell->left)) PetscCall(PetscObjectReference((PetscObject)left)); 7999566063dSJacob Faibussowitsch if ((right = shell->right)) PetscCall(PetscObjectReference((PetscObject)right)); 8001ba60692SJed Brown 8019566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(mat,&tmat)); 8021ba60692SJed Brown 8039566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mat,left,right)); 8049566063dSJacob Faibussowitsch PetscCall(MatScale(mat,scale)); 8059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 8069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 807b52f573bSBarry Smith PetscFunctionReturn(0); 808b52f573bSBarry Smith } 8096a7ede75SJakub Kruzik 8106a7ede75SJakub Kruzik /*@ 811d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 8128bd739bdSJakub Kruzik by summing or computing the product of all the matrices inside the composite matrix. 813d7f81bf2SJakub Kruzik 814b2b245abSJakub Kruzik Collective 815d7f81bf2SJakub Kruzik 816f899ff85SJose E. Roman Input Parameter: 817d7f81bf2SJakub Kruzik . mat - the composite matrix 818d7f81bf2SJakub Kruzik 8194b2558d6SJakub Kruzik Options Database Keys: 820b28d0daaSJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 821b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 822d7f81bf2SJakub Kruzik 823d7f81bf2SJakub Kruzik Level: advanced 824d7f81bf2SJakub Kruzik 825d7f81bf2SJakub Kruzik Notes: 826d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 827d7f81bf2SJakub Kruzik matrix in the composite matrix. 828d7f81bf2SJakub Kruzik 829db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE` 830d7f81bf2SJakub Kruzik 831d7f81bf2SJakub Kruzik @*/ 832d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 833d7f81bf2SJakub Kruzik { 834d7f81bf2SJakub Kruzik PetscFunctionBegin; 835d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 836cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat)); 837d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 838d7f81bf2SJakub Kruzik } 839d7f81bf2SJakub Kruzik 8406d0add67SJakub Kruzik static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat,PetscInt *nmat) 841d7f81bf2SJakub Kruzik { 842d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 843d7f81bf2SJakub Kruzik 844d7f81bf2SJakub Kruzik PetscFunctionBegin; 845d7f81bf2SJakub Kruzik *nmat = shell->nmat; 846d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 847d7f81bf2SJakub Kruzik } 848d7f81bf2SJakub Kruzik 849d7f81bf2SJakub Kruzik /*@ 8506d0add67SJakub Kruzik MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix. 8516a7ede75SJakub Kruzik 8526a7ede75SJakub Kruzik Not Collective 8536a7ede75SJakub Kruzik 8546a7ede75SJakub Kruzik Input Parameter: 855d7f81bf2SJakub Kruzik . mat - the composite matrix 8566a7ede75SJakub Kruzik 8576a7ede75SJakub Kruzik Output Parameter: 8586d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix 8596a7ede75SJakub Kruzik 8608b5c3584SJakub Kruzik Level: advanced 8616a7ede75SJakub Kruzik 862db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE` 8636a7ede75SJakub Kruzik 8646a7ede75SJakub Kruzik @*/ 8656d0add67SJakub Kruzik PetscErrorCode MatCompositeGetNumberMat(Mat mat,PetscInt *nmat) 8666a7ede75SJakub Kruzik { 8676a7ede75SJakub Kruzik PetscFunctionBegin; 868d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 869dadcf809SJacob Faibussowitsch PetscValidIntPointer(nmat,2); 870cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetNumberMat_C",(Mat,PetscInt*),(mat,nmat)); 871d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 872d7f81bf2SJakub Kruzik } 873d7f81bf2SJakub Kruzik 874d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 875d7f81bf2SJakub Kruzik { 876d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 877d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 878d7f81bf2SJakub Kruzik PetscInt k; 879d7f81bf2SJakub Kruzik 880d7f81bf2SJakub Kruzik PetscFunctionBegin; 88108401ef6SPierre Jolivet PetscCheck(i < shell->nmat,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT,i,shell->nmat); 882d7f81bf2SJakub Kruzik ilink = shell->head; 883d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 884d7f81bf2SJakub Kruzik ilink = ilink->next; 885d7f81bf2SJakub Kruzik } 886d7f81bf2SJakub Kruzik *Ai = ilink->mat; 8876a7ede75SJakub Kruzik PetscFunctionReturn(0); 8886a7ede75SJakub Kruzik } 8896a7ede75SJakub Kruzik 8908b5c3584SJakub Kruzik /*@ 8918bd739bdSJakub Kruzik MatCompositeGetMat - Returns the ith matrix from the composite matrix. 8928b5c3584SJakub Kruzik 8938b5c3584SJakub Kruzik Logically Collective on Mat 8948b5c3584SJakub Kruzik 895d8d19677SJose E. Roman Input Parameters: 896d7f81bf2SJakub Kruzik + mat - the composite matrix 8978b5c3584SJakub Kruzik - i - the number of requested matrix 8988b5c3584SJakub Kruzik 8998b5c3584SJakub Kruzik Output Parameter: 9008b5c3584SJakub Kruzik . Ai - ith matrix in composite 9018b5c3584SJakub Kruzik 9028b5c3584SJakub Kruzik Level: advanced 9038b5c3584SJakub Kruzik 904db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE` 9058b5c3584SJakub Kruzik 9068b5c3584SJakub Kruzik @*/ 907d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 9088b5c3584SJakub Kruzik { 9098b5c3584SJakub Kruzik PetscFunctionBegin; 910d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 911d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 9128b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 913cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai)); 9148b5c3584SJakub Kruzik PetscFunctionReturn(0); 9158b5c3584SJakub Kruzik } 9168b5c3584SJakub Kruzik 91703049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings_Composite(Mat mat,const PetscScalar *scalings) 91803049c21SJunchao Zhang { 91903049c21SJunchao Zhang Mat_Composite *shell = (Mat_Composite*)mat->data; 92003049c21SJunchao Zhang PetscInt nmat; 92103049c21SJunchao Zhang 92203049c21SJunchao Zhang PetscFunctionBegin; 9239566063dSJacob Faibussowitsch PetscCall(MatCompositeGetNumberMat(mat,&nmat)); 9249566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(PetscMalloc1(nmat,&shell->scalings)); 9259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(shell->scalings,scalings,nmat)); 92603049c21SJunchao Zhang PetscFunctionReturn(0); 92703049c21SJunchao Zhang } 92803049c21SJunchao Zhang 92903049c21SJunchao Zhang /*@ 93003049c21SJunchao Zhang MatCompositeSetScalings - Sets separate scaling factors for component matrices. 93103049c21SJunchao Zhang 93203049c21SJunchao Zhang Logically Collective on Mat 93303049c21SJunchao Zhang 934d8d19677SJose E. Roman Input Parameters: 93503049c21SJunchao Zhang + mat - the composite matrix 93603049c21SJunchao Zhang - scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat) 93703049c21SJunchao Zhang 93803049c21SJunchao Zhang Level: advanced 93903049c21SJunchao Zhang 940db781477SPatrick Sanan .seealso: `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE` 94103049c21SJunchao Zhang 94203049c21SJunchao Zhang @*/ 94303049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings(Mat mat,const PetscScalar *scalings) 94403049c21SJunchao Zhang { 94503049c21SJunchao Zhang PetscFunctionBegin; 94603049c21SJunchao Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 947dadcf809SJacob Faibussowitsch PetscValidScalarPointer(scalings,2); 94803049c21SJunchao Zhang PetscValidLogicalCollectiveScalar(mat,*scalings,2); 949cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetScalings_C",(Mat,const PetscScalar*),(mat,scalings)); 95003049c21SJunchao Zhang PetscFunctionReturn(0); 95103049c21SJunchao Zhang } 95203049c21SJunchao Zhang 953f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 954f4259b30SLisandro Dalcin NULL, 955f4259b30SLisandro Dalcin NULL, 95641cd0125SJakub Kruzik MatMult_Composite, 9577bf3a71bSJakub Kruzik MatMultAdd_Composite, 95841cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 9597bf3a71bSJakub Kruzik MatMultTransposeAdd_Composite, 960f4259b30SLisandro Dalcin NULL, 961f4259b30SLisandro Dalcin NULL, 962f4259b30SLisandro Dalcin NULL, 963f4259b30SLisandro Dalcin /* 10*/ NULL, 964f4259b30SLisandro Dalcin NULL, 965f4259b30SLisandro Dalcin NULL, 966f4259b30SLisandro Dalcin NULL, 967f4259b30SLisandro Dalcin NULL, 968f4259b30SLisandro Dalcin /* 15*/ NULL, 969f4259b30SLisandro Dalcin NULL, 97041cd0125SJakub Kruzik MatGetDiagonal_Composite, 97141cd0125SJakub Kruzik MatDiagonalScale_Composite, 972f4259b30SLisandro Dalcin NULL, 973f4259b30SLisandro Dalcin /* 20*/ NULL, 97441cd0125SJakub Kruzik MatAssemblyEnd_Composite, 975f4259b30SLisandro Dalcin NULL, 976f4259b30SLisandro Dalcin NULL, 977f4259b30SLisandro Dalcin /* 24*/ NULL, 978f4259b30SLisandro Dalcin NULL, 979f4259b30SLisandro Dalcin NULL, 980f4259b30SLisandro Dalcin NULL, 981f4259b30SLisandro Dalcin NULL, 982f4259b30SLisandro Dalcin /* 29*/ NULL, 983f4259b30SLisandro Dalcin NULL, 984f4259b30SLisandro Dalcin NULL, 985f4259b30SLisandro Dalcin NULL, 986f4259b30SLisandro Dalcin NULL, 987f4259b30SLisandro Dalcin /* 34*/ NULL, 988f4259b30SLisandro Dalcin NULL, 989f4259b30SLisandro Dalcin NULL, 990f4259b30SLisandro Dalcin NULL, 991f4259b30SLisandro Dalcin NULL, 992f4259b30SLisandro Dalcin /* 39*/ NULL, 993f4259b30SLisandro Dalcin NULL, 994f4259b30SLisandro Dalcin NULL, 995f4259b30SLisandro Dalcin NULL, 996f4259b30SLisandro Dalcin NULL, 997f4259b30SLisandro Dalcin /* 44*/ NULL, 99841cd0125SJakub Kruzik MatScale_Composite, 99941cd0125SJakub Kruzik MatShift_Basic, 1000f4259b30SLisandro Dalcin NULL, 1001f4259b30SLisandro Dalcin NULL, 1002f4259b30SLisandro Dalcin /* 49*/ NULL, 1003f4259b30SLisandro Dalcin NULL, 1004f4259b30SLisandro Dalcin NULL, 1005f4259b30SLisandro Dalcin NULL, 1006f4259b30SLisandro Dalcin NULL, 1007f4259b30SLisandro Dalcin /* 54*/ NULL, 1008f4259b30SLisandro Dalcin NULL, 1009f4259b30SLisandro Dalcin NULL, 1010f4259b30SLisandro Dalcin NULL, 1011f4259b30SLisandro Dalcin NULL, 1012f4259b30SLisandro Dalcin /* 59*/ NULL, 101341cd0125SJakub Kruzik MatDestroy_Composite, 1014f4259b30SLisandro Dalcin NULL, 1015f4259b30SLisandro Dalcin NULL, 1016f4259b30SLisandro Dalcin NULL, 1017f4259b30SLisandro Dalcin /* 64*/ NULL, 1018f4259b30SLisandro Dalcin NULL, 1019f4259b30SLisandro Dalcin NULL, 1020f4259b30SLisandro Dalcin NULL, 1021f4259b30SLisandro Dalcin NULL, 1022f4259b30SLisandro Dalcin /* 69*/ NULL, 1023f4259b30SLisandro Dalcin NULL, 1024f4259b30SLisandro Dalcin NULL, 1025f4259b30SLisandro Dalcin NULL, 1026f4259b30SLisandro Dalcin NULL, 1027f4259b30SLisandro Dalcin /* 74*/ NULL, 1028f4259b30SLisandro Dalcin NULL, 10294b2558d6SJakub Kruzik MatSetFromOptions_Composite, 1030f4259b30SLisandro Dalcin NULL, 1031f4259b30SLisandro Dalcin NULL, 1032f4259b30SLisandro Dalcin /* 79*/ NULL, 1033f4259b30SLisandro Dalcin NULL, 1034f4259b30SLisandro Dalcin NULL, 1035f4259b30SLisandro Dalcin NULL, 1036f4259b30SLisandro Dalcin NULL, 1037f4259b30SLisandro Dalcin /* 84*/ NULL, 1038f4259b30SLisandro Dalcin NULL, 1039f4259b30SLisandro Dalcin NULL, 1040f4259b30SLisandro Dalcin NULL, 1041f4259b30SLisandro Dalcin NULL, 1042f4259b30SLisandro Dalcin /* 89*/ NULL, 1043f4259b30SLisandro Dalcin NULL, 1044f4259b30SLisandro Dalcin NULL, 1045f4259b30SLisandro Dalcin NULL, 1046f4259b30SLisandro Dalcin NULL, 1047f4259b30SLisandro Dalcin /* 94*/ NULL, 1048f4259b30SLisandro Dalcin NULL, 1049f4259b30SLisandro Dalcin NULL, 1050f4259b30SLisandro Dalcin NULL, 1051f4259b30SLisandro Dalcin NULL, 1052f4259b30SLisandro Dalcin /*99*/ NULL, 1053f4259b30SLisandro Dalcin NULL, 1054f4259b30SLisandro Dalcin NULL, 1055f4259b30SLisandro Dalcin NULL, 1056f4259b30SLisandro Dalcin NULL, 1057f4259b30SLisandro Dalcin /*104*/ NULL, 1058f4259b30SLisandro Dalcin NULL, 1059f4259b30SLisandro Dalcin NULL, 1060f4259b30SLisandro Dalcin NULL, 1061f4259b30SLisandro Dalcin NULL, 1062f4259b30SLisandro Dalcin /*109*/ NULL, 1063f4259b30SLisandro Dalcin NULL, 1064f4259b30SLisandro Dalcin NULL, 1065f4259b30SLisandro Dalcin NULL, 1066f4259b30SLisandro Dalcin NULL, 1067f4259b30SLisandro Dalcin /*114*/ NULL, 1068f4259b30SLisandro Dalcin NULL, 1069f4259b30SLisandro Dalcin NULL, 1070f4259b30SLisandro Dalcin NULL, 1071f4259b30SLisandro Dalcin NULL, 1072f4259b30SLisandro Dalcin /*119*/ NULL, 1073f4259b30SLisandro Dalcin NULL, 1074f4259b30SLisandro Dalcin NULL, 1075f4259b30SLisandro Dalcin NULL, 1076f4259b30SLisandro Dalcin NULL, 1077f4259b30SLisandro Dalcin /*124*/ NULL, 1078f4259b30SLisandro Dalcin NULL, 1079f4259b30SLisandro Dalcin NULL, 1080f4259b30SLisandro Dalcin NULL, 1081f4259b30SLisandro Dalcin NULL, 1082f4259b30SLisandro Dalcin /*129*/ NULL, 1083f4259b30SLisandro Dalcin NULL, 1084f4259b30SLisandro Dalcin NULL, 1085f4259b30SLisandro Dalcin NULL, 1086f4259b30SLisandro Dalcin NULL, 1087f4259b30SLisandro Dalcin /*134*/ NULL, 1088f4259b30SLisandro Dalcin NULL, 1089f4259b30SLisandro Dalcin NULL, 1090f4259b30SLisandro Dalcin NULL, 1091f4259b30SLisandro Dalcin NULL, 1092f4259b30SLisandro Dalcin /*139*/ NULL, 1093f4259b30SLisandro Dalcin NULL, 1094d70f29a3SPierre Jolivet NULL, 1095d70f29a3SPierre Jolivet NULL, 1096d70f29a3SPierre Jolivet NULL, 1097d70f29a3SPierre Jolivet /*144*/ NULL, 1098d70f29a3SPierre Jolivet NULL, 1099d70f29a3SPierre Jolivet NULL, 110099a7f59eSMark Adams NULL, 110199a7f59eSMark Adams NULL, 11027fb60732SBarry Smith NULL, 11037fb60732SBarry Smith /*150*/ NULL 110441cd0125SJakub Kruzik }; 110541cd0125SJakub Kruzik 110641cd0125SJakub Kruzik /*MC 110741cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 110841cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 110941cd0125SJakub Kruzik 111041cd0125SJakub Kruzik Notes: 111141cd0125SJakub Kruzik to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 111241cd0125SJakub Kruzik 111341cd0125SJakub Kruzik Level: advanced 111441cd0125SJakub Kruzik 1115db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`, `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()` 111641cd0125SJakub Kruzik M*/ 111741cd0125SJakub Kruzik 111841cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 111941cd0125SJakub Kruzik { 112041cd0125SJakub Kruzik Mat_Composite *b; 112141cd0125SJakub Kruzik 112241cd0125SJakub Kruzik PetscFunctionBegin; 11239566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&b)); 112441cd0125SJakub Kruzik A->data = (void*)b; 11259566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps))); 112641cd0125SJakub Kruzik 11279566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 11289566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 112941cd0125SJakub Kruzik 113041cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 113141cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 113241cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 113341cd0125SJakub Kruzik b->scale = 1.0; 113441cd0125SJakub Kruzik b->nmat = 0; 11354b2558d6SJakub Kruzik b->merge = PETSC_FALSE; 11368c8613bfSJakub Kruzik b->mergetype = MAT_COMPOSITE_MERGE_RIGHT; 11373b35acafSJakub Kruzik b->structure = DIFFERENT_NONZERO_PATTERN; 113803049c21SJunchao Zhang b->merge_mvctx = PETSC_TRUE; 113903049c21SJunchao Zhang 11409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite)); 11419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite)); 11429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite)); 11439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeType_C",MatCompositeSetMergeType_Composite)); 11449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMatStructure_C",MatCompositeSetMatStructure_Composite)); 11459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMatStructure_C",MatCompositeGetMatStructure_Composite)); 11469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite)); 11479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNumberMat_C",MatCompositeGetNumberMat_Composite)); 11489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite)); 11499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetScalings_C",MatCompositeSetScalings_Composite)); 115041cd0125SJakub Kruzik 11519566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE)); 115241cd0125SJakub Kruzik PetscFunctionReturn(0); 115341cd0125SJakub Kruzik } 1154