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 359371c9d4SSatish Balay PetscErrorCode MatDestroy_Composite(Mat mat) { 362c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite *)mat->data; 376d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, oldnext; 3803049c21SJunchao Zhang PetscInt i; 39793850ffSBarry Smith 40793850ffSBarry Smith PetscFunctionBegin; 41793850ffSBarry Smith while (next) { 429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&next->mat)); 4348a46eb9SPierre Jolivet if (next->work && (!next->next || next->work != next->next->work)) PetscCall(VecDestroy(&next->work)); 446d7c1e57SBarry Smith oldnext = next; 45793850ffSBarry Smith next = next->next; 469566063dSJacob Faibussowitsch PetscCall(PetscFree(oldnext)); 47793850ffSBarry Smith } 489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->work)); 499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left)); 509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right)); 519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->leftwork)); 529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->rightwork)); 539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->leftwork2)); 549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->rightwork2)); 5503049c21SJunchao Zhang 5603049c21SJunchao Zhang if (shell->Mvctx) { 579566063dSJacob Faibussowitsch for (i = 0; i < shell->nmat; i++) PetscCall(VecDestroy(&shell->lvecs[i])); 589566063dSJacob Faibussowitsch PetscCall(PetscFree3(shell->location, shell->larray, shell->lvecs)); 599566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->larray)); 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->gvec)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->Mvctx)); 6203049c21SJunchao Zhang } 6303049c21SJunchao Zhang 649566063dSJacob Faibussowitsch PetscCall(PetscFree(shell->scalings)); 652e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeAddMat_C", NULL)); 662e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetType_C", NULL)); 672e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetType_C", NULL)); 682e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMergeType_C", NULL)); 692e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetMatStructure_C", NULL)); 702e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMatStructure_C", NULL)); 712e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeMerge_C", NULL)); 722e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetNumberMat_C", NULL)); 732e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeGetMat_C", NULL)); 742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatCompositeSetScalings_C", NULL)); 759566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 76793850ffSBarry Smith PetscFunctionReturn(0); 77793850ffSBarry Smith } 78793850ffSBarry Smith 799371c9d4SSatish Balay PetscErrorCode MatMult_Composite_Multiplicative(Mat A, Vec x, Vec y) { 806d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite *)A->data; 816d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 826d7c1e57SBarry Smith Vec in, out; 8303049c21SJunchao Zhang PetscScalar scale; 8403049c21SJunchao Zhang PetscInt i; 856d7c1e57SBarry Smith 866d7c1e57SBarry Smith PetscFunctionBegin; 8728b400f6SJacob Faibussowitsch PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()"); 886d7c1e57SBarry Smith in = x; 89e4fc5df0SBarry Smith if (shell->right) { 9048a46eb9SPierre Jolivet if (!shell->rightwork) PetscCall(VecDuplicate(shell->right, &shell->rightwork)); 919566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->rightwork, shell->right, in)); 92e4fc5df0SBarry Smith in = shell->rightwork; 93e4fc5df0SBarry Smith } 946d7c1e57SBarry Smith while (next->next) { 956d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 969566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(next->mat, NULL, &next->work)); 976d7c1e57SBarry Smith } 986d7c1e57SBarry Smith out = next->work; 999566063dSJacob Faibussowitsch PetscCall(MatMult(next->mat, in, out)); 1006d7c1e57SBarry Smith in = out; 1016d7c1e57SBarry Smith next = next->next; 1026d7c1e57SBarry Smith } 1039566063dSJacob Faibussowitsch PetscCall(MatMult(next->mat, in, y)); 1041baa6e33SBarry Smith if (shell->left) PetscCall(VecPointwiseMult(y, shell->left, y)); 10503049c21SJunchao Zhang scale = shell->scale; 1069371c9d4SSatish Balay if (shell->scalings) { 1079371c9d4SSatish Balay for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i]; 1089371c9d4SSatish Balay } 1099566063dSJacob Faibussowitsch PetscCall(VecScale(y, scale)); 1106d7c1e57SBarry Smith PetscFunctionReturn(0); 1116d7c1e57SBarry Smith } 1126d7c1e57SBarry Smith 1139371c9d4SSatish Balay PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A, Vec x, Vec y) { 1146d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite *)A->data; 1156d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 1166d7c1e57SBarry Smith Vec in, out; 11703049c21SJunchao Zhang PetscScalar scale; 11803049c21SJunchao Zhang PetscInt i; 1196d7c1e57SBarry Smith 1206d7c1e57SBarry Smith PetscFunctionBegin; 12128b400f6SJacob Faibussowitsch PetscCheck(tail, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()"); 1226d7c1e57SBarry Smith in = x; 123e4fc5df0SBarry Smith if (shell->left) { 12448a46eb9SPierre Jolivet if (!shell->leftwork) PetscCall(VecDuplicate(shell->left, &shell->leftwork)); 1259566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork, shell->left, in)); 126e4fc5df0SBarry Smith in = shell->leftwork; 127e4fc5df0SBarry Smith } 1286d7c1e57SBarry Smith while (tail->prev) { 1296d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(tail->mat, NULL, &tail->prev->work)); 1316d7c1e57SBarry Smith } 1326d7c1e57SBarry Smith out = tail->prev->work; 1339566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat, in, out)); 1346d7c1e57SBarry Smith in = out; 1356d7c1e57SBarry Smith tail = tail->prev; 1366d7c1e57SBarry Smith } 1379566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tail->mat, in, y)); 1381baa6e33SBarry Smith if (shell->right) PetscCall(VecPointwiseMult(y, shell->right, y)); 13903049c21SJunchao Zhang 14003049c21SJunchao Zhang scale = shell->scale; 1419371c9d4SSatish Balay if (shell->scalings) { 1429371c9d4SSatish Balay for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i]; 1439371c9d4SSatish Balay } 1449566063dSJacob Faibussowitsch PetscCall(VecScale(y, scale)); 1456d7c1e57SBarry Smith PetscFunctionReturn(0); 1466d7c1e57SBarry Smith } 1476d7c1e57SBarry Smith 1489371c9d4SSatish Balay PetscErrorCode MatMult_Composite(Mat mat, Vec x, Vec y) { 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) { 16348a46eb9SPierre Jolivet if (!shell->rightwork) PetscCall(VecDuplicate(shell->right, &shell->rightwork)); 1649566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->rightwork, shell->right, in)); 165e4fc5df0SBarry Smith in = shell->rightwork; 166e4fc5df0SBarry Smith } 16703049c21SJunchao Zhang 16803049c21SJunchao Zhang /* Try to merge Mvctx when instructed but not yet done. We did not do it in MatAssemblyEnd() since at that time 16903049c21SJunchao Zhang we did not know whether mat is ADDITIVE or MULTIPLICATIVE. Only now we are assured mat is ADDITIVE and 17003049c21SJunchao Zhang it is legal to merge Mvctx, because all component matrices have the same size. 17103049c21SJunchao Zhang */ 17203049c21SJunchao Zhang if (shell->merge_mvctx && !shell->Mvctx) { 17303049c21SJunchao Zhang /* Currently only implemented for MATMPIAIJ */ 17403049c21SJunchao Zhang for (cur = shell->head; cur; cur = cur->next) { 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)cur->mat, MATMPIAIJ, &match)); 17603049c21SJunchao Zhang if (!match) { 17703049c21SJunchao Zhang shell->merge_mvctx = PETSC_FALSE; 17803049c21SJunchao Zhang goto skip_merge_mvctx; 179e4fc5df0SBarry Smith } 180e4fc5df0SBarry Smith } 18103049c21SJunchao Zhang 18203049c21SJunchao Zhang /* Go through matrices first time to count total number of nonzero off-diag columns (may have dups) */ 18303049c21SJunchao Zhang tot = 0; 18403049c21SJunchao Zhang for (cur = shell->head; cur; cur = cur->next) { 1859566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, NULL)); 1869566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 18703049c21SJunchao Zhang tot += n; 18803049c21SJunchao Zhang } 1899566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tot, &shell->location, tot, &shell->larray, shell->nmat, &shell->lvecs)); 19003049c21SJunchao Zhang shell->len = tot; 19103049c21SJunchao Zhang 19203049c21SJunchao Zhang /* Go through matrices second time to sort off-diag columns and remove dups */ 1939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot, &gindices)); /* No Malloc2() since we will give one to petsc and free the other */ 1949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tot, &buf)); 19503049c21SJunchao Zhang nuniq = 0; /* Number of unique nonzero columns */ 19603049c21SJunchao Zhang for (cur = shell->head; cur; cur = cur->next) { 1979566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray)); 1989566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 19903049c21SJunchao Zhang /* Merge pre-sorted garray[0,n) and gindices[0,nuniq) to buf[] */ 20003049c21SJunchao Zhang i = j = k = 0; 20103049c21SJunchao Zhang while (i < n && j < nuniq) { 20203049c21SJunchao Zhang if (garray[i] < gindices[j]) buf[k++] = garray[i++]; 20303049c21SJunchao Zhang else if (garray[i] > gindices[j]) buf[k++] = gindices[j++]; 2049371c9d4SSatish Balay else { 2059371c9d4SSatish Balay buf[k++] = garray[i++]; 2069371c9d4SSatish Balay j++; 2079371c9d4SSatish Balay } 20803049c21SJunchao Zhang } 20903049c21SJunchao Zhang /* Copy leftover in garray[] or gindices[] */ 21003049c21SJunchao Zhang if (i < n) { 2119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf + k, garray + i, n - i)); 21203049c21SJunchao Zhang nuniq = k + n - i; 21303049c21SJunchao Zhang } else if (j < nuniq) { 2149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buf + k, gindices + j, nuniq - j)); 21503049c21SJunchao Zhang nuniq = k + nuniq - j; 21603049c21SJunchao Zhang } else nuniq = k; 21703049c21SJunchao Zhang /* Swap gindices and buf to merge garray of the next matrix */ 21803049c21SJunchao Zhang tmp = gindices; 21903049c21SJunchao Zhang gindices = buf; 22003049c21SJunchao Zhang buf = tmp; 22103049c21SJunchao Zhang } 2229566063dSJacob Faibussowitsch PetscCall(PetscFree(buf)); 22303049c21SJunchao Zhang 22403049c21SJunchao Zhang /* Go through matrices third time to build a map from gindices[] to garray[] */ 22503049c21SJunchao Zhang tot = 0; 22603049c21SJunchao Zhang for (cur = shell->head, j = 0; cur; cur = cur->next, j++) { /* j-th matrix */ 2279566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, NULL, &B, &garray)); 2289566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 2299566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, n, NULL, &shell->lvecs[j])); 23003049c21SJunchao Zhang /* This is an optimized PetscFindInt(garray[i],nuniq,gindices,&shell->location[tot+i]), using the fact that garray[] is also sorted */ 23103049c21SJunchao Zhang lo = 0; 23203049c21SJunchao Zhang for (i = 0; i < n; i++) { 23303049c21SJunchao Zhang hi = nuniq; 23403049c21SJunchao Zhang while (hi - lo > 1) { 23503049c21SJunchao Zhang mid = lo + (hi - lo) / 2; 23603049c21SJunchao Zhang if (garray[i] < gindices[mid]) hi = mid; 23703049c21SJunchao Zhang else lo = mid; 23803049c21SJunchao Zhang } 23903049c21SJunchao Zhang shell->location[tot + i] = lo; /* gindices[lo] = garray[i] */ 24003049c21SJunchao Zhang lo++; /* Since garray[i+1] > garray[i], we can safely advance lo */ 24103049c21SJunchao Zhang } 24203049c21SJunchao Zhang tot += n; 24303049c21SJunchao Zhang } 24403049c21SJunchao Zhang 24503049c21SJunchao Zhang /* Build merged Mvctx */ 2469566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nuniq, gindices, PETSC_OWN_POINTER, &ix)); 2479566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, nuniq, 0, 1, &iy)); 2489566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &xin)); 2499566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, nuniq, &shell->gvec)); 2509566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xin, ix, shell->gvec, iy, &shell->Mvctx)); 2519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xin)); 2529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ix)); 2539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iy)); 25403049c21SJunchao Zhang } 25503049c21SJunchao Zhang 25603049c21SJunchao Zhang skip_merge_mvctx: 2579566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0)); 2589566063dSJacob Faibussowitsch if (!shell->leftwork2) PetscCall(VecDuplicate(y, &shell->leftwork2)); 25903049c21SJunchao Zhang y2 = shell->leftwork2; 26003049c21SJunchao Zhang 26103049c21SJunchao Zhang if (shell->Mvctx) { /* Have a merged Mvctx */ 26203049c21SJunchao 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 26303049c21SJunchao 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 26403049c21SJunchao Zhang overlap communication/computation since we can do sAx1 while communicating x2. Here, we use the former approach. 26503049c21SJunchao Zhang */ 2669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->Mvctx, in, shell->gvec, INSERT_VALUES, SCATTER_FORWARD)); 2679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->Mvctx, in, shell->gvec, INSERT_VALUES, SCATTER_FORWARD)); 26803049c21SJunchao Zhang 2699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(shell->gvec, &vals)); 27003049c21SJunchao Zhang for (i = 0; i < shell->len; i++) shell->larray[i] = vals[shell->location[i]]; 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(shell->gvec, &vals)); 27203049c21SJunchao Zhang 27303049c21SJunchao Zhang for (cur = shell->head, tot = i = 0; cur; cur = cur->next, i++) { /* i-th matrix */ 2749566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(cur->mat, &A, &B, NULL)); 275dbbe0bcdSBarry Smith PetscUseTypeMethod(A, mult, in, y2); 2769566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, NULL, &n)); 2779566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(shell->lvecs[i], &shell->larray[tot])); 2789566063dSJacob Faibussowitsch PetscCall((*B->ops->multadd)(B, shell->lvecs[i], y2, y2)); 2799566063dSJacob Faibussowitsch PetscCall(VecResetArray(shell->lvecs[i])); 2809566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, (shell->scalings ? shell->scalings[i] : 1.0), y2)); 28103049c21SJunchao Zhang tot += n; 28203049c21SJunchao Zhang } 28303049c21SJunchao Zhang } else { 28403049c21SJunchao Zhang if (shell->scalings) { 28503049c21SJunchao Zhang for (cur = shell->head, i = 0; cur; cur = cur->next, i++) { 2869566063dSJacob Faibussowitsch PetscCall(MatMult(cur->mat, in, y2)); 2879566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, shell->scalings[i], y2)); 28803049c21SJunchao Zhang } 28903049c21SJunchao Zhang } else { 2909566063dSJacob Faibussowitsch for (cur = shell->head; cur; cur = cur->next) PetscCall(MatMultAdd(cur->mat, in, y, y)); 29103049c21SJunchao Zhang } 29203049c21SJunchao Zhang } 29303049c21SJunchao Zhang 2949566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(y, shell->left, y)); 2959566063dSJacob Faibussowitsch PetscCall(VecScale(y, shell->scale)); 296793850ffSBarry Smith PetscFunctionReturn(0); 297793850ffSBarry Smith } 298793850ffSBarry Smith 2999371c9d4SSatish Balay PetscErrorCode MatMultTranspose_Composite(Mat A, Vec x, Vec y) { 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) { 30948a46eb9SPierre Jolivet if (!shell->leftwork) PetscCall(VecDuplicate(shell->left, &shell->leftwork)); 3109566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->leftwork, shell->left, in)); 311e4fc5df0SBarry Smith in = shell->leftwork; 312e4fc5df0SBarry Smith } 31303049c21SJunchao Zhang 3149566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat, in, y)); 31503049c21SJunchao Zhang if (shell->scalings) { 3169566063dSJacob Faibussowitsch PetscCall(VecScale(y, shell->scalings[0])); 3179566063dSJacob Faibussowitsch if (!shell->rightwork2) PetscCall(VecDuplicate(y, &shell->rightwork2)); 31803049c21SJunchao Zhang y2 = shell->rightwork2; 31903049c21SJunchao Zhang } 32003049c21SJunchao Zhang i = 1; 321e4fc5df0SBarry Smith while ((next = next->next)) { 3229566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(MatMultTransposeAdd(next->mat, in, y, y)); 32303049c21SJunchao Zhang else { 3249566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(next->mat, in, y2)); 3259566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, shell->scalings[i++], y2)); 32603049c21SJunchao Zhang } 327e4fc5df0SBarry Smith } 3281baa6e33SBarry Smith if (shell->right) PetscCall(VecPointwiseMult(y, shell->right, y)); 3299566063dSJacob Faibussowitsch PetscCall(VecScale(y, shell->scale)); 330793850ffSBarry Smith PetscFunctionReturn(0); 331793850ffSBarry Smith } 332793850ffSBarry Smith 3339371c9d4SSatish Balay PetscErrorCode MatMultAdd_Composite(Mat A, Vec x, Vec y, Vec z) { 3347bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite *)A->data; 3357bf3a71bSJakub Kruzik 3367bf3a71bSJakub Kruzik PetscFunctionBegin; 3377bf3a71bSJakub Kruzik if (y != z) { 3389566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, z)); 3399566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, y)); 3407bf3a71bSJakub Kruzik } else { 34148a46eb9SPierre Jolivet if (!shell->leftwork) PetscCall(VecDuplicate(z, &shell->leftwork)); 3429566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, shell->leftwork)); 3439566063dSJacob Faibussowitsch PetscCall(VecCopy(y, z)); 3449566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, shell->leftwork)); 3457bf3a71bSJakub Kruzik } 3467bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3477bf3a71bSJakub Kruzik } 3487bf3a71bSJakub Kruzik 3499371c9d4SSatish Balay PetscErrorCode MatMultTransposeAdd_Composite(Mat A, Vec x, Vec y, Vec z) { 3507bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite *)A->data; 3517bf3a71bSJakub Kruzik 3527bf3a71bSJakub Kruzik PetscFunctionBegin; 3537bf3a71bSJakub Kruzik if (y != z) { 3549566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, z)); 3559566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, y)); 3567bf3a71bSJakub Kruzik } else { 35748a46eb9SPierre Jolivet if (!shell->rightwork) PetscCall(VecDuplicate(z, &shell->rightwork)); 3589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, shell->rightwork)); 3599566063dSJacob Faibussowitsch PetscCall(VecCopy(y, z)); 3609566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, shell->rightwork)); 3617bf3a71bSJakub Kruzik } 3627bf3a71bSJakub Kruzik PetscFunctionReturn(0); 3637bf3a71bSJakub Kruzik } 3647bf3a71bSJakub Kruzik 3659371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_Composite(Mat A, Vec v) { 366793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite *)A->data; 367793850ffSBarry Smith Mat_CompositeLink next = shell->head; 36803049c21SJunchao Zhang PetscInt i; 369793850ffSBarry Smith 370793850ffSBarry Smith PetscFunctionBegin; 37128b400f6SJacob Faibussowitsch PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()"); 37208401ef6SPierre Jolivet PetscCheck(!shell->right && !shell->left, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get diagonal if left or right scaling"); 373e4fc5df0SBarry Smith 3749566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat, v)); 3759566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(VecScale(v, shell->scalings[0])); 37603049c21SJunchao Zhang 37748a46eb9SPierre Jolivet if (next->next && !shell->work) PetscCall(VecDuplicate(v, &shell->work)); 37803049c21SJunchao Zhang i = 1; 379793850ffSBarry Smith while ((next = next->next)) { 3809566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(next->mat, shell->work)); 3819566063dSJacob Faibussowitsch PetscCall(VecAXPY(v, (shell->scalings ? shell->scalings[i++] : 1.0), shell->work)); 382793850ffSBarry Smith } 3839566063dSJacob Faibussowitsch PetscCall(VecScale(v, shell->scale)); 384793850ffSBarry Smith PetscFunctionReturn(0); 385793850ffSBarry Smith } 386793850ffSBarry Smith 3879371c9d4SSatish Balay PetscErrorCode MatAssemblyEnd_Composite(Mat Y, MatAssemblyType t) { 3884b2558d6SJakub Kruzik Mat_Composite *shell = (Mat_Composite *)Y->data; 389b52f573bSBarry Smith 390793850ffSBarry Smith PetscFunctionBegin; 3911baa6e33SBarry Smith if (shell->merge) PetscCall(MatCompositeMerge(Y)); 392793850ffSBarry Smith PetscFunctionReturn(0); 393793850ffSBarry Smith } 394793850ffSBarry Smith 3959371c9d4SSatish Balay PetscErrorCode MatScale_Composite(Mat inA, PetscScalar alpha) { 396e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite *)inA->data; 397e4fc5df0SBarry Smith 398e4fc5df0SBarry Smith PetscFunctionBegin; 399321429dbSBarry Smith a->scale *= alpha; 400e4fc5df0SBarry Smith PetscFunctionReturn(0); 401e4fc5df0SBarry Smith } 402e4fc5df0SBarry Smith 4039371c9d4SSatish Balay PetscErrorCode MatDiagonalScale_Composite(Mat inA, Vec left, Vec right) { 404e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite *)inA->data; 405e4fc5df0SBarry Smith 406e4fc5df0SBarry Smith PetscFunctionBegin; 407e4fc5df0SBarry Smith if (left) { 408321429dbSBarry Smith if (!a->left) { 4099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &a->left)); 4109566063dSJacob Faibussowitsch PetscCall(VecCopy(left, a->left)); 411321429dbSBarry Smith } else { 4129566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left, left, a->left)); 413321429dbSBarry Smith } 414e4fc5df0SBarry Smith } 415e4fc5df0SBarry Smith if (right) { 416321429dbSBarry Smith if (!a->right) { 4179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &a->right)); 4189566063dSJacob Faibussowitsch PetscCall(VecCopy(right, a->right)); 419321429dbSBarry Smith } else { 4209566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right, right, a->right)); 421321429dbSBarry Smith } 422e4fc5df0SBarry Smith } 423e4fc5df0SBarry Smith PetscFunctionReturn(0); 424e4fc5df0SBarry Smith } 425793850ffSBarry Smith 4269371c9d4SSatish Balay PetscErrorCode MatSetFromOptions_Composite(Mat A, PetscOptionItems *PetscOptionsObject) { 4274b2558d6SJakub Kruzik Mat_Composite *a = (Mat_Composite *)A->data; 4284b2558d6SJakub Kruzik 4294b2558d6SJakub Kruzik PetscFunctionBegin; 430d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "MATCOMPOSITE options"); 4319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge", "Merge at MatAssemblyEnd", "MatCompositeMerge", a->merge, &a->merge, NULL)); 4329566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_composite_merge_type", "Set composite merge direction", "MatCompositeSetMergeType", MatCompositeMergeTypes, (PetscEnum)a->mergetype, (PetscEnum *)&a->mergetype, NULL)); 4339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_composite_merge_mvctx", "Merge MatMult() vecscat contexts", "MatCreateComposite", a->merge_mvctx, &a->merge_mvctx, NULL)); 434d0609cedSBarry Smith PetscOptionsHeadEnd(); 4354b2558d6SJakub Kruzik PetscFunctionReturn(0); 4364b2558d6SJakub Kruzik } 4374b2558d6SJakub Kruzik 4382c0821f3SBarry Smith /*@ 4398bd739bdSJakub Kruzik MatCreateComposite - Creates a matrix as the sum or product of one or more matrices 440793850ffSBarry Smith 441d083f849SBarry Smith Collective 442793850ffSBarry Smith 443793850ffSBarry Smith Input Parameters: 444793850ffSBarry Smith + comm - MPI communicator 445793850ffSBarry Smith . nmat - number of matrices to put in 446793850ffSBarry Smith - mats - the matrices 447793850ffSBarry Smith 448793850ffSBarry Smith Output Parameter: 449db36af27SMatthew G. Knepley . mat - the matrix 450793850ffSBarry Smith 4514b2558d6SJakub Kruzik Options Database Keys: 45211a5261eSBarry Smith + -mat_composite_merge - merge in `MatAssemblyEnd()` 45311a5261eSBarry Smith . -mat_composite_merge_mvctx - merge Mvctx of component matrices to optimize communication in `MatMult()` for ADDITIVE matrices 454b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 4554b2558d6SJakub Kruzik 456793850ffSBarry Smith Level: advanced 457793850ffSBarry Smith 45811a5261eSBarry Smith Note: 459793850ffSBarry Smith Alternative construction 460793850ffSBarry Smith $ MatCreate(comm,&mat); 461793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 462793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 463793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 464793850ffSBarry Smith $ .... 465793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 466b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 467b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 468793850ffSBarry Smith 4696d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4706d7c1e57SBarry Smith 471db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCompositeGetMat()`, `MatCompositeMerge()`, `MatCompositeSetType()`, `MATCOMPOSITE` 472793850ffSBarry Smith @*/ 4739371c9d4SSatish Balay PetscErrorCode MatCreateComposite(MPI_Comm comm, PetscInt nmat, const Mat *mats, Mat *mat) { 474793850ffSBarry Smith PetscInt m, n, M, N, i; 475793850ffSBarry Smith 476793850ffSBarry Smith PetscFunctionBegin; 47708401ef6SPierre Jolivet PetscCheck(nmat >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must pass in at least one matrix"); 478064a246eSJacob Faibussowitsch PetscValidPointer(mat, 4); 479793850ffSBarry Smith 4809566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[0], PETSC_IGNORE, &n)); 4819566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mats[nmat - 1], &m, PETSC_IGNORE)); 4829566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[0], PETSC_IGNORE, &N)); 4839566063dSJacob Faibussowitsch PetscCall(MatGetSize(mats[nmat - 1], &M, PETSC_IGNORE)); 4849566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, mat)); 4859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*mat, m, n, M, N)); 4869566063dSJacob Faibussowitsch PetscCall(MatSetType(*mat, MATCOMPOSITE)); 48748a46eb9SPierre Jolivet for (i = 0; i < nmat; i++) PetscCall(MatCompositeAddMat(*mat, mats[i])); 4889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY)); 4899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY)); 490793850ffSBarry Smith PetscFunctionReturn(0); 491793850ffSBarry Smith } 492793850ffSBarry Smith 4939371c9d4SSatish Balay static PetscErrorCode MatCompositeAddMat_Composite(Mat mat, Mat smat) { 494d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite *)mat->data; 495d7f81bf2SJakub Kruzik Mat_CompositeLink ilink, next = shell->head; 496d7f81bf2SJakub Kruzik 497d7f81bf2SJakub Kruzik PetscFunctionBegin; 498*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilink)); 499f4259b30SLisandro Dalcin ilink->next = NULL; 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)smat)); 501d7f81bf2SJakub Kruzik ilink->mat = smat; 502d7f81bf2SJakub Kruzik 503d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 504d7f81bf2SJakub Kruzik else { 505ad540459SPierre Jolivet while (next->next) next = next->next; 506d7f81bf2SJakub Kruzik next->next = ilink; 507d7f81bf2SJakub Kruzik ilink->prev = next; 508d7f81bf2SJakub Kruzik } 509d7f81bf2SJakub Kruzik shell->tail = ilink; 510d7f81bf2SJakub Kruzik shell->nmat += 1; 51103049c21SJunchao Zhang 51203049c21SJunchao Zhang /* Retain the old scalings (if any) and expand it with a 1.0 for the newly added matrix */ 51303049c21SJunchao Zhang if (shell->scalings) { 5149566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(PetscScalar) * shell->nmat, &shell->scalings)); 51503049c21SJunchao Zhang shell->scalings[shell->nmat - 1] = 1.0; 51603049c21SJunchao Zhang } 517d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 518d7f81bf2SJakub Kruzik } 519d7f81bf2SJakub Kruzik 520793850ffSBarry Smith /*@ 5218bd739bdSJakub Kruzik MatCompositeAddMat - Add another matrix to a composite matrix. 522793850ffSBarry Smith 52311a5261eSBarry Smith Collective on mat 524793850ffSBarry Smith 525793850ffSBarry Smith Input Parameters: 526793850ffSBarry Smith + mat - the composite matrix 527793850ffSBarry Smith - smat - the partial matrix 528793850ffSBarry Smith 529793850ffSBarry Smith Level: advanced 530793850ffSBarry Smith 531db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE` 532793850ffSBarry Smith @*/ 5339371c9d4SSatish Balay PetscErrorCode MatCompositeAddMat(Mat mat, Mat smat) { 534793850ffSBarry Smith PetscFunctionBegin; 5350700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5360700a824SBarry Smith PetscValidHeaderSpecific(smat, MAT_CLASSID, 2); 537cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeAddMat_C", (Mat, Mat), (mat, smat)); 538d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 539d7f81bf2SJakub Kruzik } 540793850ffSBarry Smith 5419371c9d4SSatish Balay static PetscErrorCode MatCompositeSetType_Composite(Mat mat, MatCompositeType type) { 542d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite *)mat->data; 543d7f81bf2SJakub Kruzik 544d7f81bf2SJakub Kruzik PetscFunctionBegin; 545ced1ac25SJakub Kruzik b->type = type; 546d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 547f4259b30SLisandro Dalcin mat->ops->getdiagonal = NULL; 548d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 549d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 55003049c21SJunchao Zhang b->merge_mvctx = PETSC_FALSE; 551d7f81bf2SJakub Kruzik } else { 552d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 553d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 554d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 555793850ffSBarry Smith } 5566d7c1e57SBarry Smith PetscFunctionReturn(0); 5576d7c1e57SBarry Smith } 5586d7c1e57SBarry Smith 5592c0821f3SBarry Smith /*@ 5608bd739bdSJakub Kruzik MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product. 5616d7c1e57SBarry Smith 56211a5261eSBarry Smith Logically Collective on mat 5636d7c1e57SBarry Smith 5646d7c1e57SBarry Smith Input Parameters: 5656d7c1e57SBarry Smith . mat - the composite matrix 5666d7c1e57SBarry Smith 5676d7c1e57SBarry Smith Level: advanced 5686d7c1e57SBarry Smith 569db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeGetType()`, `MATCOMPOSITE` 5706d7c1e57SBarry Smith @*/ 5719371c9d4SSatish Balay PetscErrorCode MatCompositeSetType(Mat mat, MatCompositeType type) { 5726d7c1e57SBarry Smith PetscFunctionBegin; 573d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 574b2b245abSJakub Kruzik PetscValidLogicalCollectiveEnum(mat, type, 2); 575cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeSetType_C", (Mat, MatCompositeType), (mat, type)); 576793850ffSBarry Smith PetscFunctionReturn(0); 577793850ffSBarry Smith } 578793850ffSBarry Smith 5799371c9d4SSatish Balay static PetscErrorCode MatCompositeGetType_Composite(Mat mat, MatCompositeType *type) { 5806dbc55e5SJakub Kruzik Mat_Composite *b = (Mat_Composite *)mat->data; 5816dbc55e5SJakub Kruzik 5826dbc55e5SJakub Kruzik PetscFunctionBegin; 5836dbc55e5SJakub Kruzik *type = b->type; 5846dbc55e5SJakub Kruzik PetscFunctionReturn(0); 5856dbc55e5SJakub Kruzik } 5866dbc55e5SJakub Kruzik 5876dbc55e5SJakub Kruzik /*@ 5886dbc55e5SJakub Kruzik MatCompositeGetType - Returns type of composite. 5896dbc55e5SJakub Kruzik 5906dbc55e5SJakub Kruzik Not Collective 5916dbc55e5SJakub Kruzik 5926dbc55e5SJakub Kruzik Input Parameter: 5936dbc55e5SJakub Kruzik . mat - the composite matrix 5946dbc55e5SJakub Kruzik 5956dbc55e5SJakub Kruzik Output Parameter: 5966dbc55e5SJakub Kruzik . type - type of composite 5976dbc55e5SJakub Kruzik 5986dbc55e5SJakub Kruzik Level: advanced 5996dbc55e5SJakub Kruzik 60011a5261eSBarry Smith .seealso: `MatCreateComposite()`, `MatCompositeSetType()`, `MATCOMPOSITE`, `MatCompositeType` 6016dbc55e5SJakub Kruzik @*/ 6029371c9d4SSatish Balay PetscErrorCode MatCompositeGetType(Mat mat, MatCompositeType *type) { 6036dbc55e5SJakub Kruzik PetscFunctionBegin; 6046dbc55e5SJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6056dbc55e5SJakub Kruzik PetscValidPointer(type, 2); 606cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeGetType_C", (Mat, MatCompositeType *), (mat, type)); 6076dbc55e5SJakub Kruzik PetscFunctionReturn(0); 6086dbc55e5SJakub Kruzik } 6096dbc55e5SJakub Kruzik 6109371c9d4SSatish Balay static PetscErrorCode MatCompositeSetMatStructure_Composite(Mat mat, MatStructure str) { 6113b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite *)mat->data; 6123b35acafSJakub Kruzik 6133b35acafSJakub Kruzik PetscFunctionBegin; 6143b35acafSJakub Kruzik b->structure = str; 6153b35acafSJakub Kruzik PetscFunctionReturn(0); 6163b35acafSJakub Kruzik } 6173b35acafSJakub Kruzik 6183b35acafSJakub Kruzik /*@ 6193b35acafSJakub Kruzik MatCompositeSetMatStructure - Indicates structure of matrices in the composite matrix. 6203b35acafSJakub Kruzik 6213b35acafSJakub Kruzik Not Collective 6223b35acafSJakub Kruzik 6233b35acafSJakub Kruzik Input Parameters: 6243b35acafSJakub Kruzik + mat - the composite matrix 62511a5261eSBarry Smith - str - either `SAME_NONZERO_PATTERN`, `DIFFERENT_NONZERO_PATTERN` (default) or `SUBSET_NONZERO_PATTERN` 6263b35acafSJakub Kruzik 6273b35acafSJakub Kruzik Level: advanced 6283b35acafSJakub Kruzik 62911a5261eSBarry Smith Note: 63011a5261eSBarry Smith Information about the matrices structure is used in `MatCompositeMerge()` for additive composite matrix. 6313b35acafSJakub Kruzik 632db781477SPatrick Sanan .seealso: `MatAXPY()`, `MatCreateComposite()`, `MatCompositeMerge()` `MatCompositeGetMatStructure()`, `MATCOMPOSITE` 6333b35acafSJakub Kruzik @*/ 6349371c9d4SSatish Balay PetscErrorCode MatCompositeSetMatStructure(Mat mat, MatStructure str) { 6353b35acafSJakub Kruzik PetscFunctionBegin; 6363b35acafSJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 637cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeSetMatStructure_C", (Mat, MatStructure), (mat, str)); 6383b35acafSJakub Kruzik PetscFunctionReturn(0); 6393b35acafSJakub Kruzik } 6403b35acafSJakub Kruzik 6419371c9d4SSatish Balay static PetscErrorCode MatCompositeGetMatStructure_Composite(Mat mat, MatStructure *str) { 6423b35acafSJakub Kruzik Mat_Composite *b = (Mat_Composite *)mat->data; 6433b35acafSJakub Kruzik 6443b35acafSJakub Kruzik PetscFunctionBegin; 6453b35acafSJakub Kruzik *str = b->structure; 6463b35acafSJakub Kruzik PetscFunctionReturn(0); 6473b35acafSJakub Kruzik } 6483b35acafSJakub Kruzik 6493b35acafSJakub Kruzik /*@ 6503b35acafSJakub Kruzik MatCompositeGetMatStructure - Returns the structure of matrices in the composite matrix. 6513b35acafSJakub Kruzik 6523b35acafSJakub Kruzik Not Collective 6533b35acafSJakub Kruzik 6543b35acafSJakub Kruzik Input Parameter: 6553b35acafSJakub Kruzik . mat - the composite matrix 6563b35acafSJakub Kruzik 6573b35acafSJakub Kruzik Output Parameter: 6583b35acafSJakub Kruzik . str - structure of the matrices 6593b35acafSJakub Kruzik 6603b35acafSJakub Kruzik Level: advanced 6613b35acafSJakub Kruzik 662db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MATCOMPOSITE` 6633b35acafSJakub Kruzik @*/ 6649371c9d4SSatish Balay PetscErrorCode MatCompositeGetMatStructure(Mat mat, MatStructure *str) { 6653b35acafSJakub Kruzik PetscFunctionBegin; 6663b35acafSJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6673b35acafSJakub Kruzik PetscValidPointer(str, 2); 668cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeGetMatStructure_C", (Mat, MatStructure *), (mat, str)); 6693b35acafSJakub Kruzik PetscFunctionReturn(0); 6703b35acafSJakub Kruzik } 6713b35acafSJakub Kruzik 6729371c9d4SSatish Balay static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat, MatCompositeMergeType type) { 6738c8613bfSJakub Kruzik Mat_Composite *shell = (Mat_Composite *)mat->data; 6748c8613bfSJakub Kruzik 6758c8613bfSJakub Kruzik PetscFunctionBegin; 6768c8613bfSJakub Kruzik shell->mergetype = type; 6778c8613bfSJakub Kruzik PetscFunctionReturn(0); 6788c8613bfSJakub Kruzik } 6798c8613bfSJakub Kruzik 6808c8613bfSJakub Kruzik /*@ 68111a5261eSBarry Smith MatCompositeSetMergeType - Sets order of `MatCompositeMerge()`. 6828c8613bfSJakub Kruzik 68311a5261eSBarry Smith Logically Collective on mat 6848c8613bfSJakub Kruzik 685d8d19677SJose E. Roman Input Parameters: 6868c8613bfSJakub Kruzik + mat - the composite matrix 68711a5261eSBarry Smith - type - `MAT_COMPOSITE_MERGE RIGHT` (default) to start merge from right with the first added matrix (mat[0]), 68811a5261eSBarry Smith `MAT_COMPOSITE_MERGE_LEFT` to start merge from left with the last added matrix (mat[nmat-1]) 6898c8613bfSJakub Kruzik 6908c8613bfSJakub Kruzik Level: advanced 6918c8613bfSJakub Kruzik 69211a5261eSBarry Smith Note: 69311a5261eSBarry Smith The resulting matrix is the same regardles of the `MatCompositeMergeType`. Only the order of operation is changed. 69411a5261eSBarry Smith If set to `MAT_COMPOSITE_MERGE_RIGHT` the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0]))) 6958c8613bfSJakub Kruzik otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0]. 6968c8613bfSJakub Kruzik 697db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeMerge()`, `MATCOMPOSITE` 6988c8613bfSJakub Kruzik @*/ 6999371c9d4SSatish Balay PetscErrorCode MatCompositeSetMergeType(Mat mat, MatCompositeMergeType type) { 7008c8613bfSJakub Kruzik PetscFunctionBegin; 7018c8613bfSJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7028c8613bfSJakub Kruzik PetscValidLogicalCollectiveEnum(mat, type, 2); 703cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeSetMergeType_C", (Mat, MatCompositeMergeType), (mat, type)); 7048c8613bfSJakub Kruzik PetscFunctionReturn(0); 7058c8613bfSJakub Kruzik } 7068c8613bfSJakub Kruzik 7079371c9d4SSatish Balay static PetscErrorCode MatCompositeMerge_Composite(Mat mat) { 708b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite *)mat->data; 7096d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 7106d7c1e57SBarry Smith Mat tmat, newmat; 7111ba60692SJed Brown Vec left, right; 7121ba60692SJed Brown PetscScalar scale; 71303049c21SJunchao Zhang PetscInt i; 714b52f573bSBarry Smith 715b52f573bSBarry Smith PetscFunctionBegin; 71628b400f6SJacob Faibussowitsch PetscCheck(next, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least one matrix with MatCompositeAddMat()"); 71703049c21SJunchao Zhang scale = shell->scale; 7186d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 7198c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 72003049c21SJunchao Zhang i = 0; 7219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat)); 7229566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i++])); 72348a46eb9SPierre Jolivet while ((next = next->next)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i++] : 1.0), next->mat, shell->structure)); 7246d7c1e57SBarry Smith } else { 72503049c21SJunchao Zhang i = shell->nmat - 1; 7269566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat)); 7279566063dSJacob Faibussowitsch if (shell->scalings) PetscCall(MatScale(tmat, shell->scalings[i--])); 72848a46eb9SPierre Jolivet while ((prev = prev->prev)) PetscCall(MatAXPY(tmat, (shell->scalings ? shell->scalings[i--] : 1.0), prev->mat, shell->structure)); 7298c8613bfSJakub Kruzik } 7308c8613bfSJakub Kruzik } else { 7318c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 7329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(next->mat, MAT_COPY_VALUES, &tmat)); 733b6cfcaf5SJakub Kruzik while ((next = next->next)) { 7349566063dSJacob Faibussowitsch PetscCall(MatMatMult(next->mat, tmat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat)); 7359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 7366d7c1e57SBarry Smith tmat = newmat; 7376d7c1e57SBarry Smith } 73804d534ceSJakub Kruzik } else { 7399566063dSJacob Faibussowitsch PetscCall(MatDuplicate(prev->mat, MAT_COPY_VALUES, &tmat)); 74004d534ceSJakub Kruzik while ((prev = prev->prev)) { 7419566063dSJacob Faibussowitsch PetscCall(MatMatMult(tmat, prev->mat, MAT_INITIAL_MATRIX, PETSC_DECIDE, &newmat)); 7429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 74304d534ceSJakub Kruzik tmat = newmat; 74404d534ceSJakub Kruzik } 74504d534ceSJakub Kruzik } 7469371c9d4SSatish Balay if (shell->scalings) { 7479371c9d4SSatish Balay for (i = 0; i < shell->nmat; i++) scale *= shell->scalings[i]; 7489371c9d4SSatish Balay } 7496d7c1e57SBarry Smith } 7501ba60692SJed Brown 7519566063dSJacob Faibussowitsch if ((left = shell->left)) PetscCall(PetscObjectReference((PetscObject)left)); 7529566063dSJacob Faibussowitsch if ((right = shell->right)) PetscCall(PetscObjectReference((PetscObject)right)); 7531ba60692SJed Brown 7549566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(mat, &tmat)); 7551ba60692SJed Brown 7569566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mat, left, right)); 7579566063dSJacob Faibussowitsch PetscCall(MatScale(mat, scale)); 7589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 7599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 760b52f573bSBarry Smith PetscFunctionReturn(0); 761b52f573bSBarry Smith } 7626a7ede75SJakub Kruzik 7636a7ede75SJakub Kruzik /*@ 764d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 7658bd739bdSJakub Kruzik by summing or computing the product of all the matrices inside the composite matrix. 766d7f81bf2SJakub Kruzik 767b2b245abSJakub Kruzik Collective 768d7f81bf2SJakub Kruzik 769f899ff85SJose E. Roman Input Parameter: 770d7f81bf2SJakub Kruzik . mat - the composite matrix 771d7f81bf2SJakub Kruzik 7724b2558d6SJakub Kruzik Options Database Keys: 77311a5261eSBarry Smith + -mat_composite_merge - merge in `MatAssemblyEnd()` 774b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 775d7f81bf2SJakub Kruzik 776d7f81bf2SJakub Kruzik Level: advanced 777d7f81bf2SJakub Kruzik 77811a5261eSBarry Smith Note: 77911a5261eSBarry Smith The `MatType` of the resulting matrix will be the same as the `MatType` of the FIRST matrix in the composite matrix. 780d7f81bf2SJakub Kruzik 781db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMult()`, `MatCompositeAddMat()`, `MatCreateComposite()`, `MatCompositeSetMatStructure()`, `MatCompositeSetMergeType()`, `MATCOMPOSITE` 782d7f81bf2SJakub Kruzik @*/ 7839371c9d4SSatish Balay PetscErrorCode MatCompositeMerge(Mat mat) { 784d7f81bf2SJakub Kruzik PetscFunctionBegin; 785d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 786cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeMerge_C", (Mat), (mat)); 787d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 788d7f81bf2SJakub Kruzik } 789d7f81bf2SJakub Kruzik 7909371c9d4SSatish Balay static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat, PetscInt *nmat) { 791d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite *)mat->data; 792d7f81bf2SJakub Kruzik 793d7f81bf2SJakub Kruzik PetscFunctionBegin; 794d7f81bf2SJakub Kruzik *nmat = shell->nmat; 795d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 796d7f81bf2SJakub Kruzik } 797d7f81bf2SJakub Kruzik 798d7f81bf2SJakub Kruzik /*@ 7996d0add67SJakub Kruzik MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix. 8006a7ede75SJakub Kruzik 8016a7ede75SJakub Kruzik Not Collective 8026a7ede75SJakub Kruzik 8036a7ede75SJakub Kruzik Input Parameter: 804d7f81bf2SJakub Kruzik . mat - the composite matrix 8056a7ede75SJakub Kruzik 8066a7ede75SJakub Kruzik Output Parameter: 8076d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix 8086a7ede75SJakub Kruzik 8098b5c3584SJakub Kruzik Level: advanced 8106a7ede75SJakub Kruzik 811db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetMat()`, `MATCOMPOSITE` 8126a7ede75SJakub Kruzik @*/ 8139371c9d4SSatish Balay PetscErrorCode MatCompositeGetNumberMat(Mat mat, PetscInt *nmat) { 8146a7ede75SJakub Kruzik PetscFunctionBegin; 815d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 816dadcf809SJacob Faibussowitsch PetscValidIntPointer(nmat, 2); 817cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeGetNumberMat_C", (Mat, PetscInt *), (mat, nmat)); 818d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 819d7f81bf2SJakub Kruzik } 820d7f81bf2SJakub Kruzik 8219371c9d4SSatish Balay static PetscErrorCode MatCompositeGetMat_Composite(Mat mat, PetscInt i, Mat *Ai) { 822d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite *)mat->data; 823d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 824d7f81bf2SJakub Kruzik PetscInt k; 825d7f81bf2SJakub Kruzik 826d7f81bf2SJakub Kruzik PetscFunctionBegin; 82708401ef6SPierre Jolivet PetscCheck(i < shell->nmat, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "index out of range: %" PetscInt_FMT " >= %" PetscInt_FMT, i, shell->nmat); 828d7f81bf2SJakub Kruzik ilink = shell->head; 829ad540459SPierre Jolivet for (k = 0; k < i; k++) ilink = ilink->next; 830d7f81bf2SJakub Kruzik *Ai = ilink->mat; 8316a7ede75SJakub Kruzik PetscFunctionReturn(0); 8326a7ede75SJakub Kruzik } 8336a7ede75SJakub Kruzik 8348b5c3584SJakub Kruzik /*@ 8358bd739bdSJakub Kruzik MatCompositeGetMat - Returns the ith matrix from the composite matrix. 8368b5c3584SJakub Kruzik 83711a5261eSBarry Smith Logically Collective on mat 8388b5c3584SJakub Kruzik 839d8d19677SJose E. Roman Input Parameters: 840d7f81bf2SJakub Kruzik + mat - the composite matrix 8418b5c3584SJakub Kruzik - i - the number of requested matrix 8428b5c3584SJakub Kruzik 8438b5c3584SJakub Kruzik Output Parameter: 8448b5c3584SJakub Kruzik . Ai - ith matrix in composite 8458b5c3584SJakub Kruzik 8468b5c3584SJakub Kruzik Level: advanced 8478b5c3584SJakub Kruzik 848db781477SPatrick Sanan .seealso: `MatCreateComposite()`, `MatCompositeGetNumberMat()`, `MatCompositeAddMat()`, `MATCOMPOSITE` 8498b5c3584SJakub Kruzik @*/ 8509371c9d4SSatish Balay PetscErrorCode MatCompositeGetMat(Mat mat, PetscInt i, Mat *Ai) { 8518b5c3584SJakub Kruzik PetscFunctionBegin; 852d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 853d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat, i, 2); 8548b5c3584SJakub Kruzik PetscValidPointer(Ai, 3); 855cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeGetMat_C", (Mat, PetscInt, Mat *), (mat, i, Ai)); 8568b5c3584SJakub Kruzik PetscFunctionReturn(0); 8578b5c3584SJakub Kruzik } 8588b5c3584SJakub Kruzik 8599371c9d4SSatish Balay PetscErrorCode MatCompositeSetScalings_Composite(Mat mat, const PetscScalar *scalings) { 86003049c21SJunchao Zhang Mat_Composite *shell = (Mat_Composite *)mat->data; 86103049c21SJunchao Zhang PetscInt nmat; 86203049c21SJunchao Zhang 86303049c21SJunchao Zhang PetscFunctionBegin; 8649566063dSJacob Faibussowitsch PetscCall(MatCompositeGetNumberMat(mat, &nmat)); 8659566063dSJacob Faibussowitsch if (!shell->scalings) PetscCall(PetscMalloc1(nmat, &shell->scalings)); 8669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(shell->scalings, scalings, nmat)); 86703049c21SJunchao Zhang PetscFunctionReturn(0); 86803049c21SJunchao Zhang } 86903049c21SJunchao Zhang 87003049c21SJunchao Zhang /*@ 87103049c21SJunchao Zhang MatCompositeSetScalings - Sets separate scaling factors for component matrices. 87203049c21SJunchao Zhang 87311a5261eSBarry Smith Logically Collective on mat 87403049c21SJunchao Zhang 875d8d19677SJose E. Roman Input Parameters: 87603049c21SJunchao Zhang + mat - the composite matrix 87703049c21SJunchao Zhang - scalings - array of scaling factors with scalings[i] being factor of i-th matrix, for i in [0, nmat) 87803049c21SJunchao Zhang 87903049c21SJunchao Zhang Level: advanced 88003049c21SJunchao Zhang 881db781477SPatrick Sanan .seealso: `MatScale()`, `MatDiagonalScale()`, `MATCOMPOSITE` 88203049c21SJunchao Zhang @*/ 8839371c9d4SSatish Balay PetscErrorCode MatCompositeSetScalings(Mat mat, const PetscScalar *scalings) { 88403049c21SJunchao Zhang PetscFunctionBegin; 88503049c21SJunchao Zhang PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 886dadcf809SJacob Faibussowitsch PetscValidScalarPointer(scalings, 2); 88703049c21SJunchao Zhang PetscValidLogicalCollectiveScalar(mat, *scalings, 2); 888cac4c232SBarry Smith PetscUseMethod(mat, "MatCompositeSetScalings_C", (Mat, const PetscScalar *), (mat, scalings)); 88903049c21SJunchao Zhang PetscFunctionReturn(0); 89003049c21SJunchao Zhang } 89103049c21SJunchao Zhang 892f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 893f4259b30SLisandro Dalcin NULL, 894f4259b30SLisandro Dalcin NULL, 89541cd0125SJakub Kruzik MatMult_Composite, 8967bf3a71bSJakub Kruzik MatMultAdd_Composite, 89741cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 8987bf3a71bSJakub Kruzik MatMultTransposeAdd_Composite, 899f4259b30SLisandro Dalcin NULL, 900f4259b30SLisandro Dalcin NULL, 901f4259b30SLisandro Dalcin NULL, 902f4259b30SLisandro Dalcin /* 10*/ NULL, 903f4259b30SLisandro Dalcin NULL, 904f4259b30SLisandro Dalcin NULL, 905f4259b30SLisandro Dalcin NULL, 906f4259b30SLisandro Dalcin NULL, 907f4259b30SLisandro Dalcin /* 15*/ NULL, 908f4259b30SLisandro Dalcin NULL, 90941cd0125SJakub Kruzik MatGetDiagonal_Composite, 91041cd0125SJakub Kruzik MatDiagonalScale_Composite, 911f4259b30SLisandro Dalcin NULL, 912f4259b30SLisandro Dalcin /* 20*/ NULL, 91341cd0125SJakub Kruzik MatAssemblyEnd_Composite, 914f4259b30SLisandro Dalcin NULL, 915f4259b30SLisandro Dalcin NULL, 916f4259b30SLisandro Dalcin /* 24*/ NULL, 917f4259b30SLisandro Dalcin NULL, 918f4259b30SLisandro Dalcin NULL, 919f4259b30SLisandro Dalcin NULL, 920f4259b30SLisandro Dalcin NULL, 921f4259b30SLisandro Dalcin /* 29*/ NULL, 922f4259b30SLisandro Dalcin NULL, 923f4259b30SLisandro Dalcin NULL, 924f4259b30SLisandro Dalcin NULL, 925f4259b30SLisandro Dalcin NULL, 926f4259b30SLisandro Dalcin /* 34*/ NULL, 927f4259b30SLisandro Dalcin NULL, 928f4259b30SLisandro Dalcin NULL, 929f4259b30SLisandro Dalcin NULL, 930f4259b30SLisandro Dalcin NULL, 931f4259b30SLisandro Dalcin /* 39*/ NULL, 932f4259b30SLisandro Dalcin NULL, 933f4259b30SLisandro Dalcin NULL, 934f4259b30SLisandro Dalcin NULL, 935f4259b30SLisandro Dalcin NULL, 936f4259b30SLisandro Dalcin /* 44*/ NULL, 93741cd0125SJakub Kruzik MatScale_Composite, 93841cd0125SJakub Kruzik MatShift_Basic, 939f4259b30SLisandro Dalcin NULL, 940f4259b30SLisandro Dalcin NULL, 941f4259b30SLisandro Dalcin /* 49*/ NULL, 942f4259b30SLisandro Dalcin NULL, 943f4259b30SLisandro Dalcin NULL, 944f4259b30SLisandro Dalcin NULL, 945f4259b30SLisandro Dalcin NULL, 946f4259b30SLisandro Dalcin /* 54*/ NULL, 947f4259b30SLisandro Dalcin NULL, 948f4259b30SLisandro Dalcin NULL, 949f4259b30SLisandro Dalcin NULL, 950f4259b30SLisandro Dalcin NULL, 951f4259b30SLisandro Dalcin /* 59*/ NULL, 95241cd0125SJakub Kruzik MatDestroy_Composite, 953f4259b30SLisandro Dalcin NULL, 954f4259b30SLisandro Dalcin NULL, 955f4259b30SLisandro Dalcin NULL, 956f4259b30SLisandro Dalcin /* 64*/ NULL, 957f4259b30SLisandro Dalcin NULL, 958f4259b30SLisandro Dalcin NULL, 959f4259b30SLisandro Dalcin NULL, 960f4259b30SLisandro Dalcin NULL, 961f4259b30SLisandro Dalcin /* 69*/ NULL, 962f4259b30SLisandro Dalcin NULL, 963f4259b30SLisandro Dalcin NULL, 964f4259b30SLisandro Dalcin NULL, 965f4259b30SLisandro Dalcin NULL, 966f4259b30SLisandro Dalcin /* 74*/ NULL, 967f4259b30SLisandro Dalcin NULL, 9684b2558d6SJakub Kruzik MatSetFromOptions_Composite, 969f4259b30SLisandro Dalcin NULL, 970f4259b30SLisandro Dalcin NULL, 971f4259b30SLisandro Dalcin /* 79*/ NULL, 972f4259b30SLisandro Dalcin NULL, 973f4259b30SLisandro Dalcin NULL, 974f4259b30SLisandro Dalcin NULL, 975f4259b30SLisandro Dalcin NULL, 976f4259b30SLisandro Dalcin /* 84*/ NULL, 977f4259b30SLisandro Dalcin NULL, 978f4259b30SLisandro Dalcin NULL, 979f4259b30SLisandro Dalcin NULL, 980f4259b30SLisandro Dalcin NULL, 981f4259b30SLisandro Dalcin /* 89*/ NULL, 982f4259b30SLisandro Dalcin NULL, 983f4259b30SLisandro Dalcin NULL, 984f4259b30SLisandro Dalcin NULL, 985f4259b30SLisandro Dalcin NULL, 986f4259b30SLisandro Dalcin /* 94*/ NULL, 987f4259b30SLisandro Dalcin NULL, 988f4259b30SLisandro Dalcin NULL, 989f4259b30SLisandro Dalcin NULL, 990f4259b30SLisandro Dalcin NULL, 991f4259b30SLisandro Dalcin /*99*/ NULL, 992f4259b30SLisandro Dalcin NULL, 993f4259b30SLisandro Dalcin NULL, 994f4259b30SLisandro Dalcin NULL, 995f4259b30SLisandro Dalcin NULL, 996f4259b30SLisandro Dalcin /*104*/ NULL, 997f4259b30SLisandro Dalcin NULL, 998f4259b30SLisandro Dalcin NULL, 999f4259b30SLisandro Dalcin NULL, 1000f4259b30SLisandro Dalcin NULL, 1001f4259b30SLisandro Dalcin /*109*/ NULL, 1002f4259b30SLisandro Dalcin NULL, 1003f4259b30SLisandro Dalcin NULL, 1004f4259b30SLisandro Dalcin NULL, 1005f4259b30SLisandro Dalcin NULL, 1006f4259b30SLisandro Dalcin /*114*/ NULL, 1007f4259b30SLisandro Dalcin NULL, 1008f4259b30SLisandro Dalcin NULL, 1009f4259b30SLisandro Dalcin NULL, 1010f4259b30SLisandro Dalcin NULL, 1011f4259b30SLisandro Dalcin /*119*/ NULL, 1012f4259b30SLisandro Dalcin NULL, 1013f4259b30SLisandro Dalcin NULL, 1014f4259b30SLisandro Dalcin NULL, 1015f4259b30SLisandro Dalcin NULL, 1016f4259b30SLisandro Dalcin /*124*/ NULL, 1017f4259b30SLisandro Dalcin NULL, 1018f4259b30SLisandro Dalcin NULL, 1019f4259b30SLisandro Dalcin NULL, 1020f4259b30SLisandro Dalcin NULL, 1021f4259b30SLisandro Dalcin /*129*/ NULL, 1022f4259b30SLisandro Dalcin NULL, 1023f4259b30SLisandro Dalcin NULL, 1024f4259b30SLisandro Dalcin NULL, 1025f4259b30SLisandro Dalcin NULL, 1026f4259b30SLisandro Dalcin /*134*/ NULL, 1027f4259b30SLisandro Dalcin NULL, 1028f4259b30SLisandro Dalcin NULL, 1029f4259b30SLisandro Dalcin NULL, 1030f4259b30SLisandro Dalcin NULL, 1031f4259b30SLisandro Dalcin /*139*/ NULL, 1032f4259b30SLisandro Dalcin NULL, 1033d70f29a3SPierre Jolivet NULL, 1034d70f29a3SPierre Jolivet NULL, 1035d70f29a3SPierre Jolivet NULL, 1036d70f29a3SPierre Jolivet /*144*/ NULL, 1037d70f29a3SPierre Jolivet NULL, 1038d70f29a3SPierre Jolivet NULL, 103999a7f59eSMark Adams NULL, 104099a7f59eSMark Adams NULL, 10417fb60732SBarry Smith NULL, 10429371c9d4SSatish Balay /*150*/ NULL}; 104341cd0125SJakub Kruzik 104441cd0125SJakub Kruzik /*MC 104541cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 104641cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 104741cd0125SJakub Kruzik 104811a5261eSBarry Smith Note: 104911a5261eSBarry Smith To use the product of the matrices call `MatCompositeSetType`(mat,`MAT_COMPOSITE_MULTIPLICATIVE`); 105041cd0125SJakub Kruzik 105141cd0125SJakub Kruzik Level: advanced 105241cd0125SJakub Kruzik 105311a5261eSBarry Smith .seealso: `MatCreateComposite()`, `MatCompositeSetScalings()`, `MatCompositeAddMat()`, `MatSetType()`, `MatCompositeSetType()`, `MatCompositeGetType()`, 105411a5261eSBarry Smith `MatCompositeSetMatStructure()`, `MatCompositeGetMatStructure()`, `MatCompositeMerge()`, `MatCompositeSetMergeType()`, `MatCompositeGetNumberMat()`, `MatCompositeGetMat()` 105541cd0125SJakub Kruzik M*/ 105641cd0125SJakub Kruzik 10579371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) { 105841cd0125SJakub Kruzik Mat_Composite *b; 105941cd0125SJakub Kruzik 106041cd0125SJakub Kruzik PetscFunctionBegin; 1061*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 106241cd0125SJakub Kruzik A->data = (void *)b; 10639566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops, &MatOps_Values, sizeof(struct _MatOps))); 106441cd0125SJakub Kruzik 10659566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 10669566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 106741cd0125SJakub Kruzik 106841cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 106941cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 107041cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 107141cd0125SJakub Kruzik b->scale = 1.0; 107241cd0125SJakub Kruzik b->nmat = 0; 10734b2558d6SJakub Kruzik b->merge = PETSC_FALSE; 10748c8613bfSJakub Kruzik b->mergetype = MAT_COMPOSITE_MERGE_RIGHT; 10753b35acafSJakub Kruzik b->structure = DIFFERENT_NONZERO_PATTERN; 107603049c21SJunchao Zhang b->merge_mvctx = PETSC_TRUE; 107703049c21SJunchao Zhang 10789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeAddMat_C", MatCompositeAddMat_Composite)); 10799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetType_C", MatCompositeSetType_Composite)); 10809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetType_C", MatCompositeGetType_Composite)); 10819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMergeType_C", MatCompositeSetMergeType_Composite)); 10829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetMatStructure_C", MatCompositeSetMatStructure_Composite)); 10839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMatStructure_C", MatCompositeGetMatStructure_Composite)); 10849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeMerge_C", MatCompositeMerge_Composite)); 10859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetNumberMat_C", MatCompositeGetNumberMat_Composite)); 10869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeGetMat_C", MatCompositeGetMat_Composite)); 10879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatCompositeSetScalings_C", MatCompositeSetScalings_Composite)); 108841cd0125SJakub Kruzik 10899566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCOMPOSITE)); 109041cd0125SJakub Kruzik PetscFunctionReturn(0); 109141cd0125SJakub Kruzik } 1092