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)); 689566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 69793850ffSBarry Smith PetscFunctionReturn(0); 70793850ffSBarry Smith } 71793850ffSBarry Smith 726d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 736d7c1e57SBarry Smith { 746d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 756d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 766d7c1e57SBarry Smith Vec in,out; 7703049c21SJunchao Zhang PetscScalar scale; 7803049c21SJunchao Zhang PetscInt i; 796d7c1e57SBarry Smith 806d7c1e57SBarry Smith PetscFunctionBegin; 8128b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 826d7c1e57SBarry Smith in = x; 83e4fc5df0SBarry Smith if (shell->right) { 84e4fc5df0SBarry Smith if (!shell->rightwork) { 859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->right,&shell->rightwork)); 86e4fc5df0SBarry Smith } 879566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->rightwork,shell->right,in)); 88e4fc5df0SBarry Smith in = shell->rightwork; 89e4fc5df0SBarry Smith } 906d7c1e57SBarry Smith while (next->next) { 916d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 929566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(next->mat,NULL,&next->work)); 936d7c1e57SBarry Smith } 946d7c1e57SBarry Smith out = next->work; 959566063dSJacob Faibussowitsch PetscCall(MatMult(next->mat,in,out)); 966d7c1e57SBarry Smith in = out; 976d7c1e57SBarry Smith next = next->next; 986d7c1e57SBarry Smith } 999566063dSJacob Faibussowitsch PetscCall(MatMult(next->mat,in,y)); 100e4fc5df0SBarry Smith if (shell->left) { 1019566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,shell->left,y)); 102e4fc5df0SBarry Smith } 10303049c21SJunchao Zhang scale = shell->scale; 10403049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 1059566063dSJacob Faibussowitsch PetscCall(VecScale(y,scale)); 1066d7c1e57SBarry Smith PetscFunctionReturn(0); 1076d7c1e57SBarry Smith } 1086d7c1e57SBarry Smith 1096d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 1106d7c1e57SBarry Smith { 1116d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 1126d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 1136d7c1e57SBarry Smith Vec in,out; 11403049c21SJunchao Zhang PetscScalar scale; 11503049c21SJunchao Zhang PetscInt i; 1166d7c1e57SBarry Smith 1176d7c1e57SBarry Smith PetscFunctionBegin; 11828b400f6SJacob Faibussowitsch PetscCheck(tail,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 1196d7c1e57SBarry Smith in = x; 120e4fc5df0SBarry Smith if (shell->left) { 121e4fc5df0SBarry Smith if (!shell->leftwork) { 1229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left,&shell->leftwork)); 123e4fc5df0SBarry Smith } 1249566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in)); 125e4fc5df0SBarry Smith in = shell->leftwork; 126e4fc5df0SBarry Smith } 1276d7c1e57SBarry Smith while (tail->prev) { 1286d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(tail->mat,NULL,&tail->prev->work)); 1306d7c1e57SBarry Smith } 1316d7c1e57SBarry Smith out = tail->prev->work; 1329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat,in,out)); 1336d7c1e57SBarry Smith in = out; 1346d7c1e57SBarry Smith tail = tail->prev; 1356d7c1e57SBarry Smith } 1369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat,in,y)); 137e4fc5df0SBarry Smith if (shell->right) { 1389566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,shell->right,y)); 139e4fc5df0SBarry Smith } 14003049c21SJunchao Zhang 14103049c21SJunchao Zhang scale = shell->scale; 14203049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 1439566063dSJacob Faibussowitsch PetscCall(VecScale(y,scale)); 1446d7c1e57SBarry Smith PetscFunctionReturn(0); 1456d7c1e57SBarry Smith } 1466d7c1e57SBarry Smith 14703049c21SJunchao Zhang PetscErrorCode MatMult_Composite(Mat mat,Vec x,Vec y) 148793850ffSBarry Smith { 14903049c21SJunchao Zhang Mat_Composite *shell = (Mat_Composite*)mat->data; 15003049c21SJunchao Zhang Mat_CompositeLink cur = shell->head; 15103049c21SJunchao Zhang Vec in,y2,xin; 15203049c21SJunchao Zhang Mat A,B; 15303049c21SJunchao Zhang PetscInt i,j,k,n,nuniq,lo,hi,mid,*gindices,*buf,*tmp,tot; 15403049c21SJunchao Zhang const PetscScalar *vals; 15503049c21SJunchao Zhang const PetscInt *garray; 15603049c21SJunchao Zhang IS ix,iy; 15703049c21SJunchao Zhang PetscBool match; 158793850ffSBarry Smith 159793850ffSBarry Smith PetscFunctionBegin; 16028b400f6SJacob Faibussowitsch PetscCheck(cur,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 161e4fc5df0SBarry Smith in = x; 162e4fc5df0SBarry Smith if (shell->right) { 163e4fc5df0SBarry Smith if (!shell->rightwork) { 1649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->right,&shell->rightwork)); 165793850ffSBarry Smith } 1669566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->rightwork,shell->right,in)); 167e4fc5df0SBarry Smith in = shell->rightwork; 168e4fc5df0SBarry Smith } 16903049c21SJunchao Zhang 17003049c21SJunchao Zhang /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time 17103049c21SJunchao Zhang we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and 17203049c21SJunchao Zhang it is legal to merge Mvctx, because all component matrices have the same size. 17303049c21SJunchao Zhang */ 17403049c21SJunchao Zhang if (shell->merge_mvctx && !shell->Mvctx) { 17503049c21SJunchao Zhang /* Currently only implemented for MATMPIAIJ */ 17603049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat,MATMPIAIJ,&match)); 17803049c21SJunchao Zhang if (!match) { 17903049c21SJunchao Zhang shell->merge_mvctx = PETSC_FALSE; 18003049c21SJunchao Zhang goto skip_merge_mvctx; 181e4fc5df0SBarry Smith } 182e4fc5df0SBarry Smith } 18303049c21SJunchao Zhang 18403049c21SJunchao Zhang /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */ 18503049c21SJunchao Zhang tot = 0; 18603049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 1879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,NULL)); 1889566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 18903049c21SJunchao Zhang tot += n; 19003049c21SJunchao Zhang } 1919566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tot,&shell->location,tot,&shell->larray,shell->nmat,&shell->lvecs)); 19203049c21SJunchao Zhang shell->len = tot; 19303049c21SJunchao Zhang 19403049c21SJunchao Zhang /* Go through matrices second time to sort off-diag columns and remove dups */ 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot,&gindices)); /* No Malloc2() since we will give one to petsc and free the other */ 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot,&buf)); 19703049c21SJunchao Zhang nuniq = 0; /* Number of unique nonzero columns */ 19803049c21SJunchao Zhang for (cur=shell->head; cur; cur=cur->next) { 1999566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray)); 2009566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 20103049c21SJunchao Zhang /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */ 20203049c21SJunchao Zhang i = j = k = 0; 20303049c21SJunchao Zhang while (i < n && j < nuniq) { 20403049c21SJunchao Zhang if (garray[i] < gindices[j]) buf[k++] = garray[i++]; 20503049c21SJunchao Zhang else if (garray[i] > gindices[j]) buf[k++] = gindices[j++]; 20603049c21SJunchao Zhang else {buf[k++] = garray[i++]; j++;} 20703049c21SJunchao Zhang } 20803049c21SJunchao Zhang /* Copy leftover in garray[] or gindices[] */ 20903049c21SJunchao Zhang if (i < n) { 2109566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf+k,garray+i,n-i)); 21103049c21SJunchao Zhang nuniq = k + n-i; 21203049c21SJunchao Zhang } else if (j < nuniq) { 2139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf+k,gindices+j,nuniq-j)); 21403049c21SJunchao Zhang nuniq = k + nuniq-j; 21503049c21SJunchao Zhang } else nuniq = k; 21603049c21SJunchao Zhang /* Swap gindices and buf to merge garray of the next matrix */ 21703049c21SJunchao Zhang tmp = gindices; 21803049c21SJunchao Zhang gindices = buf; 21903049c21SJunchao Zhang buf = tmp; 22003049c21SJunchao Zhang } 2219566063dSJacob Faibussowitsch PetscCall(PetscFree(buf)); 22203049c21SJunchao Zhang 22303049c21SJunchao Zhang /* Go through matrices third time to build a map from gindices[] to garray[] */ 22403049c21SJunchao Zhang tot = 0; 22503049c21SJunchao Zhang for (cur=shell->head,j=0; cur; cur=cur->next,j++) { /* j-th matrix */ 2269566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,NULL,&B,&garray)); 2279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 2289566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&shell->lvecs[j])); 22903049c21SJunchao Zhang /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */ 23003049c21SJunchao Zhang lo = 0; 23103049c21SJunchao Zhang for (i=0; i<n; i++) { 23203049c21SJunchao Zhang hi = nuniq; 23303049c21SJunchao Zhang while (hi - lo > 1) { 23403049c21SJunchao Zhang mid = lo + (hi - lo)/2; 23503049c21SJunchao Zhang if (garray[i] < gindices[mid]) hi = mid; 23603049c21SJunchao Zhang else lo = mid; 23703049c21SJunchao Zhang } 23803049c21SJunchao Zhang shell->location[tot+i] = lo; /* gindices[lo] = garray[i] */ 23903049c21SJunchao Zhang lo++; /* Since garray[i+1] > garray[i], we can safely advance lo */ 24003049c21SJunchao Zhang } 24103049c21SJunchao Zhang tot += n; 24203049c21SJunchao Zhang } 24303049c21SJunchao Zhang 24403049c21SJunchao Zhang /* Build merged Mvctx */ 2459566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nuniq,gindices,PETSC_OWN_POINTER,&ix)); 2469566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,nuniq,0,1,&iy)); 2479566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&xin)); 2489566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,nuniq,&shell->gvec)); 2499566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xin,ix,shell->gvec,iy,&shell->Mvctx)); 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xin)); 2519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ix)); 2529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iy)); 25303049c21SJunchao Zhang } 25403049c21SJunchao Zhang 25503049c21SJunchao Zhang skip_merge_mvctx: 2569566063dSJacob Faibussowitsch PetscCall(VecSet(y,0)); 2579566063dSJacob Faibussowitsch if (!shell->leftwork2) PetscCall(VecDuplicate(y,&shell->leftwork2)); 25803049c21SJunchao Zhang y2 = shell->leftwork2; 25903049c21SJunchao Zhang 26003049c21SJunchao Zhang if (shell->Mvctx) { /* Have a merged Mvctx */ 26103049c21SJunchao 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 26203049c21SJunchao 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 26303049c21SJunchao Zhang overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach. 26403049c21SJunchao Zhang */ 2659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD)); 2669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->Mvctx,in,shell->gvec,INSERT_VALUES,SCATTER_FORWARD)); 26703049c21SJunchao Zhang 2689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(shell->gvec,&vals)); 26903049c21SJunchao Zhang for (i=0; i<shell->len; i++) shell->larray[i] = vals[shell->location[i]]; 2709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(shell->gvec,&vals)); 27103049c21SJunchao Zhang 27203049c21SJunchao Zhang for (cur=shell->head,tot=i=0; cur; cur=cur->next,i++) { /* i-th matrix */ 2739566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat,&A,&B,NULL)); 2749566063dSJacob Faibussowitsch PetscCall((*A->ops->mult)(A,in,y2)); 2759566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 2769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(shell->lvecs[i],&shell->larray[tot])); 2779566063dSJacob Faibussowitsch PetscCall((*B->ops->multadd)(B,shell->lvecs[i],y2,y2)); 2789566063dSJacob Faibussowitsch PetscCall(VecResetArray(shell->lvecs[i])); 2799566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,(shell->scalings ? shell->scalings[i] : 1.0),y2)); 28003049c21SJunchao Zhang tot += n; 28103049c21SJunchao Zhang } 28203049c21SJunchao Zhang } else { 28303049c21SJunchao Zhang if (shell->scalings) { 28403049c21SJunchao Zhang for (cur=shell->head,i=0; cur; cur=cur->next,i++) { 2859566063dSJacob Faibussowitsch PetscCall(MatMult(cur->mat,in,y2)); 2869566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->scalings[i],y2)); 28703049c21SJunchao Zhang } 28803049c21SJunchao Zhang } else { 2899566063dSJacob Faibussowitsch for (cur=shell->head; cur; cur=cur->next) PetscCall(MatMultAdd(cur->mat,in,y,y)); 29003049c21SJunchao Zhang } 29103049c21SJunchao Zhang } 29203049c21SJunchao Zhang 2939566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(y,shell->left,y)); 2949566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scale)); 295793850ffSBarry Smith PetscFunctionReturn(0); 296793850ffSBarry Smith } 297793850ffSBarry Smith 298793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 299793850ffSBarry Smith { 300793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 301793850ffSBarry Smith Mat_CompositeLink next = shell->head; 30203049c21SJunchao Zhang Vec in,y2 = NULL; 30303049c21SJunchao Zhang PetscInt i; 304793850ffSBarry Smith 305793850ffSBarry Smith PetscFunctionBegin; 30628b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 307e4fc5df0SBarry Smith in = x; 308e4fc5df0SBarry Smith if (shell->left) { 309e4fc5df0SBarry Smith if (!shell->leftwork) { 3109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left,&shell->leftwork)); 311793850ffSBarry Smith } 3129566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork,shell->left,in)); 313e4fc5df0SBarry Smith in = shell->leftwork; 314e4fc5df0SBarry Smith } 31503049c21SJunchao Zhang 3169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat,in,y)); 31703049c21SJunchao Zhang if (shell->scalings) { 3189566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scalings[0])); 3199566063dSJacob Faibussowitsch if (!shell->rightwork2) PetscCall(VecDuplicate(y,&shell->rightwork2)); 32003049c21SJunchao Zhang y2 = shell->rightwork2; 32103049c21SJunchao Zhang } 32203049c21SJunchao Zhang i = 1; 323e4fc5df0SBarry Smith while ((next = next->next)) { 3249566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat,in,y,y)); 32503049c21SJunchao Zhang else { 3269566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat,in,y2)); 3279566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->scalings[i++],y2)); 32803049c21SJunchao Zhang } 329e4fc5df0SBarry Smith } 330e4fc5df0SBarry Smith if (shell->right) { 3319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,shell->right,y)); 332e4fc5df0SBarry Smith } 3339566063dSJacob Faibussowitsch PetscCall(VecScale(y,shell->scale)); 334793850ffSBarry Smith PetscFunctionReturn(0); 335793850ffSBarry Smith } 336793850ffSBarry Smith 3377bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z) 3387bf3a71bSJakub Kruzik { 3397bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 3407bf3a71bSJakub Kruzik 3417bf3a71bSJakub Kruzik PetscFunctionBegin; 3427bf3a71bSJakub Kruzik if (y != z) { 3439566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,z)); 3449566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 3457bf3a71bSJakub Kruzik } else { 3467bf3a71bSJakub Kruzik if (!shell->leftwork) { 3479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(z,&shell->leftwork)); 3487bf3a71bSJakub Kruzik } 3499566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,shell->leftwork)); 3509566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 3519566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->leftwork)); 3527bf3a71bSJakub Kruzik } 3537bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3547bf3a71bSJakub Kruzik } 3557bf3a71bSJakub Kruzik 3567bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z) 3577bf3a71bSJakub Kruzik { 3587bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 3597bf3a71bSJakub Kruzik 3607bf3a71bSJakub Kruzik PetscFunctionBegin; 3617bf3a71bSJakub Kruzik if (y != z) { 3629566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,z)); 3639566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 3647bf3a71bSJakub Kruzik } else { 3657bf3a71bSJakub Kruzik if (!shell->rightwork) { 3669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(z,&shell->rightwork)); 3677bf3a71bSJakub Kruzik } 3689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,shell->rightwork)); 3699566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 3709566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->rightwork)); 3717bf3a71bSJakub Kruzik } 3727bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3737bf3a71bSJakub Kruzik } 3747bf3a71bSJakub Kruzik 375793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 376793850ffSBarry Smith { 377793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 378793850ffSBarry Smith Mat_CompositeLink next = shell->head; 37903049c21SJunchao Zhang PetscInt i; 380793850ffSBarry Smith 381793850ffSBarry Smith PetscFunctionBegin; 38228b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 3832c71b3e2SJacob Faibussowitsch PetscCheckFalse(shell->right || shell->left,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 384e4fc5df0SBarry Smith 3859566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat,v)); 3869566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(VecScale(v,shell->scalings[0])); 38703049c21SJunchao Zhang 388793850ffSBarry Smith if (next->next && !shell->work) { 3899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v,&shell->work)); 390793850ffSBarry Smith } 39103049c21SJunchao Zhang i = 1; 392793850ffSBarry Smith while ((next = next->next)) { 3939566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat,shell->work)); 3949566063dSJacob Faibussowitsch PetscCall(VecAXPY(v,(shell->scalings ? shell->scalings[i++] : 1.0),shell->work)); 395793850ffSBarry Smith } 3969566063dSJacob Faibussowitsch PetscCall(VecScale(v,shell->scale)); 397793850ffSBarry Smith PetscFunctionReturn(0); 398793850ffSBarry Smith } 399793850ffSBarry Smith 400793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 401793850ffSBarry Smith { 4024b2558d6SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)Y->data; 403b52f573bSBarry Smith 404793850ffSBarry Smith PetscFunctionBegin; 4054b2558d6SJakub Kruzik if (shell->merge) { 4069566063dSJacob Faibussowitsch PetscCall(MatCompositeMerge(Y)); 407b52f573bSBarry Smith } 408793850ffSBarry Smith PetscFunctionReturn(0); 409793850ffSBarry Smith } 410793850ffSBarry Smith 411e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 412e4fc5df0SBarry Smith { 413e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 414e4fc5df0SBarry Smith 415e4fc5df0SBarry Smith PetscFunctionBegin; 416321429dbSBarry Smith a->scale *= alpha; 417e4fc5df0SBarry Smith PetscFunctionReturn(0); 418e4fc5df0SBarry Smith } 419e4fc5df0SBarry Smith 420e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 421e4fc5df0SBarry Smith { 422e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 423e4fc5df0SBarry Smith 424e4fc5df0SBarry Smith PetscFunctionBegin; 425e4fc5df0SBarry Smith if (left) { 426321429dbSBarry Smith if (!a->left) { 4279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&a->left)); 4289566063dSJacob Faibussowitsch PetscCall(VecCopy(left,a->left)); 429321429dbSBarry Smith } else { 4309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left,left,a->left)); 431321429dbSBarry Smith } 432e4fc5df0SBarry Smith } 433e4fc5df0SBarry Smith if (right) { 434321429dbSBarry Smith if (!a->right) { 4359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&a->right)); 4369566063dSJacob Faibussowitsch PetscCall(VecCopy(right,a->right)); 437321429dbSBarry Smith } else { 4389566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right,right,a->right)); 439321429dbSBarry Smith } 440e4fc5df0SBarry Smith } 441e4fc5df0SBarry Smith PetscFunctionReturn(0); 442e4fc5df0SBarry Smith } 443793850ffSBarry Smith 4444b2558d6SJakub Kruzik PetscErrorCode MatSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,Mat A) 4454b2558d6SJakub Kruzik { 4464b2558d6SJakub Kruzik Mat_Composite *a = (Mat_Composite*)A->data; 4474b2558d6SJakub Kruzik 4484b2558d6SJakub Kruzik PetscFunctionBegin; 4499566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"MATCOMPOSITE options")); 4509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge","Merge at MatAssemblyEnd","MatCompositeMerge",a->merge,&a->merge,NULL)); 4519566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_composite_merge_type","Set composite merge direction","MatCompositeSetMergeType",MatCompositeMergeTypes,(PetscEnum)a->mergetype,(PetscEnum*)&a->mergetype,NULL)); 4529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx","Merge MatMult() vecscat contexts","MatCreateComposite",a->merge_mvctx,&a->merge_mvctx,NULL)); 4539566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 4544b2558d6SJakub Kruzik PetscFunctionReturn(0); 4554b2558d6SJakub Kruzik } 4564b2558d6SJakub Kruzik 4572c0821f3SBarry Smith /*@ 4588bd739bdSJakub Kruzik MatCreateComposite - Creates a matrix as the sum or product of one or more matrices 459793850ffSBarry Smith 460d083f849SBarry Smith Collective 461793850ffSBarry Smith 462793850ffSBarry Smith Input Parameters: 463793850ffSBarry Smith + comm - MPI communicator 464793850ffSBarry Smith . nmat - number of matrices to put in 465793850ffSBarry Smith - mats - the matrices 466793850ffSBarry Smith 467793850ffSBarry Smith Output Parameter: 468db36af27SMatthew G. Knepley . mat - the matrix 469793850ffSBarry Smith 4704b2558d6SJakub Kruzik Options Database Keys: 4714b2558d6SJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 47203049c21SJunchao Zhang . -mat_composite_merge_mvctx - merge Mvctx of component matrices to optimize communication in MatMult() for ADDITIVE matrices 473b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 4744b2558d6SJakub Kruzik 475793850ffSBarry Smith Level: advanced 476793850ffSBarry Smith 477793850ffSBarry Smith Notes: 478793850ffSBarry Smith Alternative construction 479793850ffSBarry Smith $ MatCreate(comm,&mat); 480793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 481793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 482793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 483793850ffSBarry Smith $ .... 484793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 485b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 486b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 487793850ffSBarry Smith 4886d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4896d7c1e57SBarry Smith 4906dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeGetMat(), MatCompositeMerge(), MatCompositeSetType(), MATCOMPOSITE 491793850ffSBarry Smith 492793850ffSBarry Smith @*/ 4937087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 494793850ffSBarry Smith { 495793850ffSBarry Smith PetscInt m,n,M,N,i; 496793850ffSBarry Smith 497793850ffSBarry Smith PetscFunctionBegin; 4982c71b3e2SJacob Faibussowitsch PetscCheckFalse(nmat < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 499064a246eSJacob Faibussowitsch PetscValidPointer(mat,4); 500793850ffSBarry Smith 5019566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[0],PETSC_IGNORE,&n)); 5029566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE)); 5039566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[0],PETSC_IGNORE,&N)); 5049566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[nmat-1],&M,PETSC_IGNORE)); 5059566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,mat)); 5069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*mat,m,n,M,N)); 5079566063dSJacob Faibussowitsch PetscCall(MatSetType(*mat,MATCOMPOSITE)); 508793850ffSBarry Smith for (i=0; i<nmat; i++) { 5099566063dSJacob Faibussowitsch PetscCall(MatCompositeAddMat(*mat,mats[i])); 510793850ffSBarry Smith } 5119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY)); 5129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY)); 513793850ffSBarry Smith PetscFunctionReturn(0); 514793850ffSBarry Smith } 515793850ffSBarry Smith 516d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 517d7f81bf2SJakub Kruzik { 518d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 519d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 520d7f81bf2SJakub Kruzik 521d7f81bf2SJakub Kruzik PetscFunctionBegin; 5229566063dSJacob Faibussowitsch PetscCall(PetscNewLog(mat,&ilink)); 523f4259b30SLisandro Dalcin ilink->next = NULL; 5249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)smat)); 525d7f81bf2SJakub Kruzik ilink->mat = smat; 526d7f81bf2SJakub Kruzik 527d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 528d7f81bf2SJakub Kruzik else { 529d7f81bf2SJakub Kruzik while (next->next) { 530d7f81bf2SJakub Kruzik next = next->next; 531d7f81bf2SJakub Kruzik } 532d7f81bf2SJakub Kruzik next->next = ilink; 533d7f81bf2SJakub Kruzik ilink->prev = next; 534d7f81bf2SJakub Kruzik } 535d7f81bf2SJakub Kruzik shell->tail = ilink; 536d7f81bf2SJakub Kruzik shell->nmat += 1; 53703049c21SJunchao Zhang 53803049c21SJunchao Zhang /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */ 53903049c21SJunchao Zhang if (shell->scalings) { 5409566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(PetscScalar)*shell->nmat,&shell->scalings)); 54103049c21SJunchao Zhang shell->scalings[shell->nmat-1] = 1.0; 54203049c21SJunchao Zhang } 543d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 544d7f81bf2SJakub Kruzik } 545d7f81bf2SJakub Kruzik 546793850ffSBarry Smith /*@ 5478bd739bdSJakub Kruzik MatCompositeAddMat - Add another matrix to a composite matrix. 548793850ffSBarry Smith 549793850ffSBarry Smith Collective on Mat 550793850ffSBarry Smith 551793850ffSBarry Smith Input Parameters: 552793850ffSBarry Smith + mat - the composite matrix 553793850ffSBarry Smith - smat - the partial matrix 554793850ffSBarry Smith 555793850ffSBarry Smith Level: advanced 556793850ffSBarry Smith 5578bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE 558793850ffSBarry Smith @*/ 5597087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 560793850ffSBarry Smith { 561793850ffSBarry Smith PetscFunctionBegin; 5620700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5630700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 564*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat)); 565d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 566d7f81bf2SJakub Kruzik } 567793850ffSBarry Smith 568d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 569d7f81bf2SJakub Kruzik { 570d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 571d7f81bf2SJakub Kruzik 572d7f81bf2SJakub Kruzik PetscFunctionBegin; 573ced1ac25SJakub Kruzik b->type = type; 574d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 575f4259b30SLisandro Dalcin mat->ops->getdiagonal = NULL; 576d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 577d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 57803049c21SJunchao Zhang b->merge_mvctx = PETSC_FALSE; 579d7f81bf2SJakub Kruzik } else { 580d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 581d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 582d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 583793850ffSBarry Smith } 5846d7c1e57SBarry Smith PetscFunctionReturn(0); 5856d7c1e57SBarry Smith } 5866d7c1e57SBarry Smith 5872c0821f3SBarry Smith /*@ 5888bd739bdSJakub Kruzik MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product. 5896d7c1e57SBarry Smith 590b2b245abSJakub Kruzik Logically Collective on Mat 5916d7c1e57SBarry Smith 5926d7c1e57SBarry Smith Input Parameters: 5936d7c1e57SBarry Smith . mat - the composite matrix 5946d7c1e57SBarry Smith 5956d7c1e57SBarry Smith Level: advanced 5966d7c1e57SBarry Smith 5976dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeGetType(), MATCOMPOSITE 5986d7c1e57SBarry Smith 5996d7c1e57SBarry Smith @*/ 6007087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 6016d7c1e57SBarry Smith { 6026d7c1e57SBarry Smith PetscFunctionBegin; 603d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 604b2b245abSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 605*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type)); 606793850ffSBarry Smith PetscFunctionReturn(0); 607793850ffSBarry Smith } 608793850ffSBarry Smith 6096dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type) 6106dbc55e5SJakub Kruzik { 6116dbc55e5SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6126dbc55e5SJakub Kruzik 6136dbc55e5SJakub Kruzik PetscFunctionBegin; 6146dbc55e5SJakub Kruzik *type = b->type; 6156dbc55e5SJakub Kruzik PetscFunctionReturn(0); 6166dbc55e5SJakub Kruzik } 6176dbc55e5SJakub Kruzik 6186dbc55e5SJakub Kruzik /*@ 6196dbc55e5SJakub Kruzik MatCompositeGetType - Returns type of composite. 6206dbc55e5SJakub Kruzik 6216dbc55e5SJakub Kruzik Not Collective 6226dbc55e5SJakub Kruzik 6236dbc55e5SJakub Kruzik Input Parameter: 6246dbc55e5SJakub Kruzik . mat - the composite matrix 6256dbc55e5SJakub Kruzik 6266dbc55e5SJakub Kruzik Output Parameter: 6276dbc55e5SJakub Kruzik . type - type of composite 6286dbc55e5SJakub Kruzik 6296dbc55e5SJakub Kruzik Level: advanced 6306dbc55e5SJakub Kruzik 6316dbc55e5SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetType(), MATCOMPOSITE 6326dbc55e5SJakub Kruzik 6336dbc55e5SJakub Kruzik @*/ 6346dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type) 6356dbc55e5SJakub Kruzik { 6366dbc55e5SJakub Kruzik PetscFunctionBegin; 6376dbc55e5SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6386dbc55e5SJakub Kruzik PetscValidPointer(type,2); 639*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type)); 6406dbc55e5SJakub Kruzik PetscFunctionReturn(0); 6416dbc55e5SJakub Kruzik } 6426dbc55e5SJakub Kruzik 6433b35acafSJakub Kruzik static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat,MatStructure str) 6443b35acafSJakub Kruzik { 6453b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6463b35acafSJakub Kruzik 6473b35acafSJakub Kruzik PetscFunctionBegin; 6483b35acafSJakub Kruzik b->structure = str; 6493b35acafSJakub Kruzik PetscFunctionReturn(0); 6503b35acafSJakub Kruzik } 6513b35acafSJakub Kruzik 6523b35acafSJakub Kruzik /*@ 6533b35acafSJakub Kruzik MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix. 6543b35acafSJakub Kruzik 6553b35acafSJakub Kruzik Not Collective 6563b35acafSJakub Kruzik 6573b35acafSJakub Kruzik Input Parameters: 6583b35acafSJakub Kruzik + mat - the composite matrix 6593b35acafSJakub Kruzik - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN (default) or SUBSET_NONZERO_PATTERN 6603b35acafSJakub Kruzik 6613b35acafSJakub Kruzik Level: advanced 6623b35acafSJakub Kruzik 6633b35acafSJakub Kruzik Notes: 6643b35acafSJakub Kruzik Information about the matrices structure is used in MatCompositeMerge() for additive composite matrix. 6653b35acafSJakub Kruzik 6663b35acafSJakub Kruzik .seealso: MatAXPY(), MatCreateComposite(), MatCompositeMerge() MatCompositeGetMatStructure(), MATCOMPOSITE 6673b35acafSJakub Kruzik 6683b35acafSJakub Kruzik @*/ 6693b35acafSJakub Kruzik PetscErrorCode MatCompositeSetMatStructure(Mat mat,MatStructure str) 6703b35acafSJakub Kruzik { 6713b35acafSJakub Kruzik PetscFunctionBegin; 6723b35acafSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 673*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetMatStructure_C",(Mat,MatStructure),(mat,str)); 6743b35acafSJakub Kruzik PetscFunctionReturn(0); 6753b35acafSJakub Kruzik } 6763b35acafSJakub Kruzik 6773b35acafSJakub Kruzik static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat,MatStructure *str) 6783b35acafSJakub Kruzik { 6793b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 6803b35acafSJakub Kruzik 6813b35acafSJakub Kruzik PetscFunctionBegin; 6823b35acafSJakub Kruzik *str = b->structure; 6833b35acafSJakub Kruzik PetscFunctionReturn(0); 6843b35acafSJakub Kruzik } 6853b35acafSJakub Kruzik 6863b35acafSJakub Kruzik /*@ 6873b35acafSJakub Kruzik MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix. 6883b35acafSJakub Kruzik 6893b35acafSJakub Kruzik Not Collective 6903b35acafSJakub Kruzik 6913b35acafSJakub Kruzik Input Parameter: 6923b35acafSJakub Kruzik . mat - the composite matrix 6933b35acafSJakub Kruzik 6943b35acafSJakub Kruzik Output Parameter: 6953b35acafSJakub Kruzik . str - structure of the matrices 6963b35acafSJakub Kruzik 6973b35acafSJakub Kruzik Level: advanced 6983b35acafSJakub Kruzik 6993b35acafSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetMatStructure(), MATCOMPOSITE 7003b35acafSJakub Kruzik 7013b35acafSJakub Kruzik @*/ 7023b35acafSJakub Kruzik PetscErrorCode MatCompositeGetMatStructure(Mat mat,MatStructure *str) 7033b35acafSJakub Kruzik { 7043b35acafSJakub Kruzik PetscFunctionBegin; 7053b35acafSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7063b35acafSJakub Kruzik PetscValidPointer(str,2); 707*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetMatStructure_C",(Mat,MatStructure*),(mat,str)); 7083b35acafSJakub Kruzik PetscFunctionReturn(0); 7093b35acafSJakub Kruzik } 7103b35acafSJakub Kruzik 7118c8613bfSJakub Kruzik static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type) 7128c8613bfSJakub Kruzik { 7138c8613bfSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 7148c8613bfSJakub Kruzik 7158c8613bfSJakub Kruzik PetscFunctionBegin; 7168c8613bfSJakub Kruzik shell->mergetype = type; 7178c8613bfSJakub Kruzik PetscFunctionReturn(0); 7188c8613bfSJakub Kruzik } 7198c8613bfSJakub Kruzik 7208c8613bfSJakub Kruzik /*@ 7218c8613bfSJakub Kruzik MatCompositeSetMergeType - Sets order of MatCompositeMerge(). 7228c8613bfSJakub Kruzik 7238c8613bfSJakub Kruzik Logically Collective on Mat 7248c8613bfSJakub Kruzik 725d8d19677SJose E. Roman Input Parameters: 7268c8613bfSJakub Kruzik + mat - the composite matrix 7278c8613bfSJakub Kruzik - type - MAT_COMPOSITE_MERGE RIGHT (default) to start merge from right with the first added matrix (mat[0]), 7288c8613bfSJakub Kruzik MAT_COMPOSITE_MERGE_LEFT to start merge from left with the last added matrix (mat[nmat-1]) 7298c8613bfSJakub Kruzik 7308c8613bfSJakub Kruzik Level: advanced 7318c8613bfSJakub Kruzik 7328c8613bfSJakub Kruzik Notes: 7338c8613bfSJakub Kruzik The resulting matrix is the same regardles of the MergeType. Only the order of operation is changed. 7348c8613bfSJakub Kruzik If set to MAT_COMPOSITE_MERGE_RIGHT the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0]))) 7358c8613bfSJakub Kruzik otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0]. 7368c8613bfSJakub Kruzik 7378c8613bfSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeMerge(), MATCOMPOSITE 7388c8613bfSJakub Kruzik 7398c8613bfSJakub Kruzik @*/ 7408c8613bfSJakub Kruzik PetscErrorCode MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type) 7418c8613bfSJakub Kruzik { 7428c8613bfSJakub Kruzik PetscFunctionBegin; 7438c8613bfSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7448c8613bfSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 745*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetMergeType_C",(Mat,MatCompositeMergeType),(mat,type)); 7468c8613bfSJakub Kruzik PetscFunctionReturn(0); 7478c8613bfSJakub Kruzik } 7488c8613bfSJakub Kruzik 749d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 750b52f573bSBarry Smith { 751b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 7526d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 7536d7c1e57SBarry Smith Mat tmat,newmat; 7541ba60692SJed Brown Vec left,right; 7551ba60692SJed Brown PetscScalar scale; 75603049c21SJunchao Zhang PetscInt i; 757b52f573bSBarry Smith 758b52f573bSBarry Smith PetscFunctionBegin; 75928b400f6SJacob Faibussowitsch PetscCheck(next,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 76003049c21SJunchao Zhang scale = shell->scale; 7616d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 7628c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 76303049c21SJunchao Zhang i = 0; 7649566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat)); 7659566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i++])); 766b52f573bSBarry Smith while ((next = next->next)) { 7679566063dSJacob Faibussowitsch PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i++] : 1.0),next->mat,shell->structure)); 768b52f573bSBarry Smith } 7696d7c1e57SBarry Smith } else { 77003049c21SJunchao Zhang i = shell->nmat-1; 7719566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat)); 7729566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat,shell->scalings[i--])); 7738c8613bfSJakub Kruzik while ((prev = prev->prev)) { 7749566063dSJacob Faibussowitsch PetscCall(MatAXPY(tmat,(shell->scalings ? shell->scalings[i--] : 1.0),prev->mat,shell->structure)); 7758c8613bfSJakub Kruzik } 7768c8613bfSJakub Kruzik } 7778c8613bfSJakub Kruzik } else { 7788c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 7799566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat)); 780b6cfcaf5SJakub Kruzik while ((next = next->next)) { 7819566063dSJacob Faibussowitsch PetscCall(MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat)); 7829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 7836d7c1e57SBarry Smith tmat = newmat; 7846d7c1e57SBarry Smith } 78504d534ceSJakub Kruzik } else { 7869566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat)); 78704d534ceSJakub Kruzik while ((prev = prev->prev)) { 7889566063dSJacob Faibussowitsch PetscCall(MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat)); 7899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 79004d534ceSJakub Kruzik tmat = newmat; 79104d534ceSJakub Kruzik } 79204d534ceSJakub Kruzik } 79303049c21SJunchao Zhang if (shell->scalings) {for (i=0; i<shell->nmat; i++) scale *= shell->scalings[i];} 7946d7c1e57SBarry Smith } 7951ba60692SJed Brown 7969566063dSJacob Faibussowitsch if ((left = shell->left)) PetscCall(PetscObjectReference((PetscObject)left)); 7979566063dSJacob Faibussowitsch if ((right = shell->right)) PetscCall(PetscObjectReference((PetscObject)right)); 7981ba60692SJed Brown 7999566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(mat,&tmat)); 8001ba60692SJed Brown 8019566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mat,left,right)); 8029566063dSJacob Faibussowitsch PetscCall(MatScale(mat,scale)); 8039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 8049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 805b52f573bSBarry Smith PetscFunctionReturn(0); 806b52f573bSBarry Smith } 8076a7ede75SJakub Kruzik 8086a7ede75SJakub Kruzik /*@ 809d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 8108bd739bdSJakub Kruzik by summing or computing the product of all the matrices inside the composite matrix. 811d7f81bf2SJakub Kruzik 812b2b245abSJakub Kruzik Collective 813d7f81bf2SJakub Kruzik 814f899ff85SJose E. Roman Input Parameter: 815d7f81bf2SJakub Kruzik . mat - the composite matrix 816d7f81bf2SJakub Kruzik 8174b2558d6SJakub Kruzik Options Database Keys: 818b28d0daaSJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 819b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 820d7f81bf2SJakub Kruzik 821d7f81bf2SJakub Kruzik Level: advanced 822d7f81bf2SJakub Kruzik 823d7f81bf2SJakub Kruzik Notes: 824d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 825d7f81bf2SJakub Kruzik matrix in the composite matrix. 826d7f81bf2SJakub Kruzik 8273b35acafSJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeSetMatStructure(), MatCompositeSetMergeType(), MATCOMPOSITE 828d7f81bf2SJakub Kruzik 829d7f81bf2SJakub Kruzik @*/ 830d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 831d7f81bf2SJakub Kruzik { 832d7f81bf2SJakub Kruzik PetscFunctionBegin; 833d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 834*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat)); 835d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 836d7f81bf2SJakub Kruzik } 837d7f81bf2SJakub Kruzik 8386d0add67SJakub Kruzik static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat,PetscInt *nmat) 839d7f81bf2SJakub Kruzik { 840d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 841d7f81bf2SJakub Kruzik 842d7f81bf2SJakub Kruzik PetscFunctionBegin; 843d7f81bf2SJakub Kruzik *nmat = shell->nmat; 844d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 845d7f81bf2SJakub Kruzik } 846d7f81bf2SJakub Kruzik 847d7f81bf2SJakub Kruzik /*@ 8486d0add67SJakub Kruzik MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix. 8496a7ede75SJakub Kruzik 8506a7ede75SJakub Kruzik Not Collective 8516a7ede75SJakub Kruzik 8526a7ede75SJakub Kruzik Input Parameter: 853d7f81bf2SJakub Kruzik . mat - the composite matrix 8546a7ede75SJakub Kruzik 8556a7ede75SJakub Kruzik Output Parameter: 8566d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix 8576a7ede75SJakub Kruzik 8588b5c3584SJakub Kruzik Level: advanced 8596a7ede75SJakub Kruzik 8608bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE 8616a7ede75SJakub Kruzik 8626a7ede75SJakub Kruzik @*/ 8636d0add67SJakub Kruzik PetscErrorCode MatCompositeGetNumberMat(Mat mat,PetscInt *nmat) 8646a7ede75SJakub Kruzik { 8656a7ede75SJakub Kruzik PetscFunctionBegin; 866d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 867dadcf809SJacob Faibussowitsch PetscValidIntPointer(nmat,2); 868*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetNumberMat_C",(Mat,PetscInt*),(mat,nmat)); 869d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 870d7f81bf2SJakub Kruzik } 871d7f81bf2SJakub Kruzik 872d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 873d7f81bf2SJakub Kruzik { 874d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 875d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 876d7f81bf2SJakub Kruzik PetscInt k; 877d7f81bf2SJakub Kruzik 878d7f81bf2SJakub Kruzik PetscFunctionBegin; 8792c71b3e2SJacob Faibussowitsch PetscCheckFalse(i >= shell->nmat,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT,i,shell->nmat); 880d7f81bf2SJakub Kruzik ilink = shell->head; 881d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 882d7f81bf2SJakub Kruzik ilink = ilink->next; 883d7f81bf2SJakub Kruzik } 884d7f81bf2SJakub Kruzik *Ai = ilink->mat; 8856a7ede75SJakub Kruzik PetscFunctionReturn(0); 8866a7ede75SJakub Kruzik } 8876a7ede75SJakub Kruzik 8888b5c3584SJakub Kruzik /*@ 8898bd739bdSJakub Kruzik MatCompositeGetMat - Returns the ith matrix from the composite matrix. 8908b5c3584SJakub Kruzik 8918b5c3584SJakub Kruzik Logically Collective on Mat 8928b5c3584SJakub Kruzik 893d8d19677SJose E. Roman Input Parameters: 894d7f81bf2SJakub Kruzik + mat - the composite matrix 8958b5c3584SJakub Kruzik - i - the number of requested matrix 8968b5c3584SJakub Kruzik 8978b5c3584SJakub Kruzik Output Parameter: 8988b5c3584SJakub Kruzik . Ai - ith matrix in composite 8998b5c3584SJakub Kruzik 9008b5c3584SJakub Kruzik Level: advanced 9018b5c3584SJakub Kruzik 9026d0add67SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNumberMat(), MatCompositeAddMat(), MATCOMPOSITE 9038b5c3584SJakub Kruzik 9048b5c3584SJakub Kruzik @*/ 905d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 9068b5c3584SJakub Kruzik { 9078b5c3584SJakub Kruzik PetscFunctionBegin; 908d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 909d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 9108b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 911*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai)); 9128b5c3584SJakub Kruzik PetscFunctionReturn(0); 9138b5c3584SJakub Kruzik } 9148b5c3584SJakub Kruzik 91503049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings_Composite(Mat mat,const PetscScalar *scalings) 91603049c21SJunchao Zhang { 91703049c21SJunchao Zhang Mat_Composite *shell = (Mat_Composite*)mat->data; 91803049c21SJunchao Zhang PetscInt nmat; 91903049c21SJunchao Zhang 92003049c21SJunchao Zhang PetscFunctionBegin; 9219566063dSJacob Faibussowitsch PetscCall(MatCompositeGetNumberMat(mat,&nmat)); 9229566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(PetscMalloc1(nmat,&shell->scalings)); 9239566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(shell->scalings,scalings,nmat)); 92403049c21SJunchao Zhang PetscFunctionReturn(0); 92503049c21SJunchao Zhang } 92603049c21SJunchao Zhang 92703049c21SJunchao Zhang /*@ 92803049c21SJunchao Zhang MatCompositeSetScalings - Sets separate scaling factors for component matrices. 92903049c21SJunchao Zhang 93003049c21SJunchao Zhang Logically Collective on Mat 93103049c21SJunchao Zhang 932d8d19677SJose E. Roman Input Parameters: 93303049c21SJunchao Zhang + mat - the composite matrix 93403049c21SJunchao Zhang - scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat) 93503049c21SJunchao Zhang 93603049c21SJunchao Zhang Level: advanced 93703049c21SJunchao Zhang 93803049c21SJunchao Zhang .seealso: MatScale(), MatDiagonalScale(), MATCOMPOSITE 93903049c21SJunchao Zhang 94003049c21SJunchao Zhang @*/ 94103049c21SJunchao Zhang PetscErrorCode MatCompositeSetScalings(Mat mat,const PetscScalar *scalings) 94203049c21SJunchao Zhang { 94303049c21SJunchao Zhang PetscFunctionBegin; 94403049c21SJunchao Zhang PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 945dadcf809SJacob Faibussowitsch PetscValidScalarPointer(scalings,2); 94603049c21SJunchao Zhang PetscValidLogicalCollectiveScalar(mat,*scalings,2); 947*cac4c232SBarry Smith PetscUseMethod(mat,"MatCompositeSetScalings_C",(Mat,const PetscScalar*),(mat,scalings)); 94803049c21SJunchao Zhang PetscFunctionReturn(0); 94903049c21SJunchao Zhang } 95003049c21SJunchao Zhang 951f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 952f4259b30SLisandro Dalcin NULL, 953f4259b30SLisandro Dalcin NULL, 95441cd0125SJakub Kruzik MatMult_Composite, 9557bf3a71bSJakub Kruzik MatMultAdd_Composite, 95641cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 9577bf3a71bSJakub Kruzik MatMultTransposeAdd_Composite, 958f4259b30SLisandro Dalcin NULL, 959f4259b30SLisandro Dalcin NULL, 960f4259b30SLisandro Dalcin NULL, 961f4259b30SLisandro Dalcin /* 10*/ NULL, 962f4259b30SLisandro Dalcin NULL, 963f4259b30SLisandro Dalcin NULL, 964f4259b30SLisandro Dalcin NULL, 965f4259b30SLisandro Dalcin NULL, 966f4259b30SLisandro Dalcin /* 15*/ NULL, 967f4259b30SLisandro Dalcin NULL, 96841cd0125SJakub Kruzik MatGetDiagonal_Composite, 96941cd0125SJakub Kruzik MatDiagonalScale_Composite, 970f4259b30SLisandro Dalcin NULL, 971f4259b30SLisandro Dalcin /* 20*/ NULL, 97241cd0125SJakub Kruzik MatAssemblyEnd_Composite, 973f4259b30SLisandro Dalcin NULL, 974f4259b30SLisandro Dalcin NULL, 975f4259b30SLisandro Dalcin /* 24*/ NULL, 976f4259b30SLisandro Dalcin NULL, 977f4259b30SLisandro Dalcin NULL, 978f4259b30SLisandro Dalcin NULL, 979f4259b30SLisandro Dalcin NULL, 980f4259b30SLisandro Dalcin /* 29*/ NULL, 981f4259b30SLisandro Dalcin NULL, 982f4259b30SLisandro Dalcin NULL, 983f4259b30SLisandro Dalcin NULL, 984f4259b30SLisandro Dalcin NULL, 985f4259b30SLisandro Dalcin /* 34*/ NULL, 986f4259b30SLisandro Dalcin NULL, 987f4259b30SLisandro Dalcin NULL, 988f4259b30SLisandro Dalcin NULL, 989f4259b30SLisandro Dalcin NULL, 990f4259b30SLisandro Dalcin /* 39*/ NULL, 991f4259b30SLisandro Dalcin NULL, 992f4259b30SLisandro Dalcin NULL, 993f4259b30SLisandro Dalcin NULL, 994f4259b30SLisandro Dalcin NULL, 995f4259b30SLisandro Dalcin /* 44*/ NULL, 99641cd0125SJakub Kruzik MatScale_Composite, 99741cd0125SJakub Kruzik MatShift_Basic, 998f4259b30SLisandro Dalcin NULL, 999f4259b30SLisandro Dalcin NULL, 1000f4259b30SLisandro Dalcin /* 49*/ NULL, 1001f4259b30SLisandro Dalcin NULL, 1002f4259b30SLisandro Dalcin NULL, 1003f4259b30SLisandro Dalcin NULL, 1004f4259b30SLisandro Dalcin NULL, 1005f4259b30SLisandro Dalcin /* 54*/ NULL, 1006f4259b30SLisandro Dalcin NULL, 1007f4259b30SLisandro Dalcin NULL, 1008f4259b30SLisandro Dalcin NULL, 1009f4259b30SLisandro Dalcin NULL, 1010f4259b30SLisandro Dalcin /* 59*/ NULL, 101141cd0125SJakub Kruzik MatDestroy_Composite, 1012f4259b30SLisandro Dalcin NULL, 1013f4259b30SLisandro Dalcin NULL, 1014f4259b30SLisandro Dalcin NULL, 1015f4259b30SLisandro Dalcin /* 64*/ NULL, 1016f4259b30SLisandro Dalcin NULL, 1017f4259b30SLisandro Dalcin NULL, 1018f4259b30SLisandro Dalcin NULL, 1019f4259b30SLisandro Dalcin NULL, 1020f4259b30SLisandro Dalcin /* 69*/ NULL, 1021f4259b30SLisandro Dalcin NULL, 1022f4259b30SLisandro Dalcin NULL, 1023f4259b30SLisandro Dalcin NULL, 1024f4259b30SLisandro Dalcin NULL, 1025f4259b30SLisandro Dalcin /* 74*/ NULL, 1026f4259b30SLisandro Dalcin NULL, 10274b2558d6SJakub Kruzik MatSetFromOptions_Composite, 1028f4259b30SLisandro Dalcin NULL, 1029f4259b30SLisandro Dalcin NULL, 1030f4259b30SLisandro Dalcin /* 79*/ NULL, 1031f4259b30SLisandro Dalcin NULL, 1032f4259b30SLisandro Dalcin NULL, 1033f4259b30SLisandro Dalcin NULL, 1034f4259b30SLisandro Dalcin NULL, 1035f4259b30SLisandro Dalcin /* 84*/ NULL, 1036f4259b30SLisandro Dalcin NULL, 1037f4259b30SLisandro Dalcin NULL, 1038f4259b30SLisandro Dalcin NULL, 1039f4259b30SLisandro Dalcin NULL, 1040f4259b30SLisandro Dalcin /* 89*/ NULL, 1041f4259b30SLisandro Dalcin NULL, 1042f4259b30SLisandro Dalcin NULL, 1043f4259b30SLisandro Dalcin NULL, 1044f4259b30SLisandro Dalcin NULL, 1045f4259b30SLisandro Dalcin /* 94*/ NULL, 1046f4259b30SLisandro Dalcin NULL, 1047f4259b30SLisandro Dalcin NULL, 1048f4259b30SLisandro Dalcin NULL, 1049f4259b30SLisandro Dalcin NULL, 1050f4259b30SLisandro Dalcin /*99*/ NULL, 1051f4259b30SLisandro Dalcin NULL, 1052f4259b30SLisandro Dalcin NULL, 1053f4259b30SLisandro Dalcin NULL, 1054f4259b30SLisandro Dalcin NULL, 1055f4259b30SLisandro Dalcin /*104*/ NULL, 1056f4259b30SLisandro Dalcin NULL, 1057f4259b30SLisandro Dalcin NULL, 1058f4259b30SLisandro Dalcin NULL, 1059f4259b30SLisandro Dalcin NULL, 1060f4259b30SLisandro Dalcin /*109*/ NULL, 1061f4259b30SLisandro Dalcin NULL, 1062f4259b30SLisandro Dalcin NULL, 1063f4259b30SLisandro Dalcin NULL, 1064f4259b30SLisandro Dalcin NULL, 1065f4259b30SLisandro Dalcin /*114*/ NULL, 1066f4259b30SLisandro Dalcin NULL, 1067f4259b30SLisandro Dalcin NULL, 1068f4259b30SLisandro Dalcin NULL, 1069f4259b30SLisandro Dalcin NULL, 1070f4259b30SLisandro Dalcin /*119*/ NULL, 1071f4259b30SLisandro Dalcin NULL, 1072f4259b30SLisandro Dalcin NULL, 1073f4259b30SLisandro Dalcin NULL, 1074f4259b30SLisandro Dalcin NULL, 1075f4259b30SLisandro Dalcin /*124*/ NULL, 1076f4259b30SLisandro Dalcin NULL, 1077f4259b30SLisandro Dalcin NULL, 1078f4259b30SLisandro Dalcin NULL, 1079f4259b30SLisandro Dalcin NULL, 1080f4259b30SLisandro Dalcin /*129*/ NULL, 1081f4259b30SLisandro Dalcin NULL, 1082f4259b30SLisandro Dalcin NULL, 1083f4259b30SLisandro Dalcin NULL, 1084f4259b30SLisandro Dalcin NULL, 1085f4259b30SLisandro Dalcin /*134*/ NULL, 1086f4259b30SLisandro Dalcin NULL, 1087f4259b30SLisandro Dalcin NULL, 1088f4259b30SLisandro Dalcin NULL, 1089f4259b30SLisandro Dalcin NULL, 1090f4259b30SLisandro Dalcin /*139*/ NULL, 1091f4259b30SLisandro Dalcin NULL, 1092d70f29a3SPierre Jolivet NULL, 1093d70f29a3SPierre Jolivet NULL, 1094d70f29a3SPierre Jolivet NULL, 1095d70f29a3SPierre Jolivet /*144*/NULL, 1096d70f29a3SPierre Jolivet NULL, 1097d70f29a3SPierre Jolivet NULL, 1098f4259b30SLisandro Dalcin NULL 109941cd0125SJakub Kruzik }; 110041cd0125SJakub Kruzik 110141cd0125SJakub Kruzik /*MC 110241cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 110341cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 110441cd0125SJakub Kruzik 110541cd0125SJakub Kruzik Notes: 110641cd0125SJakub Kruzik to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 110741cd0125SJakub Kruzik 110841cd0125SJakub Kruzik Level: advanced 110941cd0125SJakub Kruzik 111003049c21SJunchao Zhang .seealso: MatCreateComposite(), MatCompositeSetScalings(), MatCompositeAddMat(), MatSetType(), MatCompositeSetType(), MatCompositeGetType(), MatCompositeSetMatStructure(), MatCompositeGetMatStructure(), MatCompositeMerge(), MatCompositeSetMergeType(), MatCompositeGetNumberMat(), MatCompositeGetMat() 111141cd0125SJakub Kruzik M*/ 111241cd0125SJakub Kruzik 111341cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 111441cd0125SJakub Kruzik { 111541cd0125SJakub Kruzik Mat_Composite *b; 111641cd0125SJakub Kruzik 111741cd0125SJakub Kruzik PetscFunctionBegin; 11189566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&b)); 111941cd0125SJakub Kruzik A->data = (void*)b; 11209566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps))); 112141cd0125SJakub Kruzik 11229566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 11239566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 112441cd0125SJakub Kruzik 112541cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 112641cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 112741cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 112841cd0125SJakub Kruzik b->scale = 1.0; 112941cd0125SJakub Kruzik b->nmat = 0; 11304b2558d6SJakub Kruzik b->merge = PETSC_FALSE; 11318c8613bfSJakub Kruzik b->mergetype = MAT_COMPOSITE_MERGE_RIGHT; 11323b35acafSJakub Kruzik b->structure = DIFFERENT_NONZERO_PATTERN; 113303049c21SJunchao Zhang b->merge_mvctx = PETSC_TRUE; 113403049c21SJunchao Zhang 11359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite)); 11369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite)); 11379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite)); 11389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeType_C",MatCompositeSetMergeType_Composite)); 11399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMatStructure_C",MatCompositeSetMatStructure_Composite)); 11409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMatStructure_C",MatCompositeGetMatStructure_Composite)); 11419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite)); 11429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNumberMat_C",MatCompositeGetNumberMat_Composite)); 11439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite)); 11449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetScalings_C",MatCompositeSetScalings_Composite)); 114541cd0125SJakub Kruzik 11469566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE)); 114741cd0125SJakub Kruzik PetscFunctionReturn(0); 114841cd0125SJakub Kruzik } 1149