1793850ffSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3793850ffSBarry Smith 4793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink; 5793850ffSBarry Smith struct _Mat_CompositeLink { 6793850ffSBarry Smith Mat mat; 76d7c1e57SBarry Smith Vec work; 86d7c1e57SBarry Smith Mat_CompositeLink next,prev; 9793850ffSBarry Smith }; 10793850ffSBarry Smith 11793850ffSBarry Smith typedef struct { 126d7c1e57SBarry Smith MatCompositeType type; 136d7c1e57SBarry Smith Mat_CompositeLink head,tail; 14793850ffSBarry Smith Vec work; 15e4fc5df0SBarry Smith PetscScalar scale; /* scale factor supplied with MatScale() */ 16e4fc5df0SBarry Smith Vec left,right; /* left and right diagonal scaling provided with MatDiagonalScale() */ 17e4fc5df0SBarry Smith Vec leftwork,rightwork; 186a7ede75SJakub Kruzik PetscInt nmat; 19793850ffSBarry Smith } Mat_Composite; 20793850ffSBarry Smith 21793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat) 22793850ffSBarry Smith { 23793850ffSBarry Smith PetscErrorCode ierr; 242c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite*)mat->data; 256d7c1e57SBarry Smith Mat_CompositeLink next = shell->head,oldnext; 26793850ffSBarry Smith 27793850ffSBarry Smith PetscFunctionBegin; 28793850ffSBarry Smith while (next) { 296bf464f9SBarry Smith ierr = MatDestroy(&next->mat);CHKERRQ(ierr); 306d7c1e57SBarry Smith if (next->work && (!next->next || next->work != next->next->work)) { 316bf464f9SBarry Smith ierr = VecDestroy(&next->work);CHKERRQ(ierr); 326d7c1e57SBarry Smith } 336d7c1e57SBarry Smith oldnext = next; 34793850ffSBarry Smith next = next->next; 356d7c1e57SBarry Smith ierr = PetscFree(oldnext);CHKERRQ(ierr); 36793850ffSBarry Smith } 376bf464f9SBarry Smith ierr = VecDestroy(&shell->work);CHKERRQ(ierr); 386bf464f9SBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 396bf464f9SBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 406bf464f9SBarry Smith ierr = VecDestroy(&shell->leftwork);CHKERRQ(ierr); 416bf464f9SBarry Smith ierr = VecDestroy(&shell->rightwork);CHKERRQ(ierr); 426bf464f9SBarry Smith ierr = PetscFree(mat->data);CHKERRQ(ierr); 43793850ffSBarry Smith PetscFunctionReturn(0); 44793850ffSBarry Smith } 45793850ffSBarry Smith 466d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 476d7c1e57SBarry Smith { 486d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 496d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 506d7c1e57SBarry Smith PetscErrorCode ierr; 516d7c1e57SBarry Smith Vec in,out; 526d7c1e57SBarry Smith 536d7c1e57SBarry Smith PetscFunctionBegin; 54e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 556d7c1e57SBarry Smith in = x; 56e4fc5df0SBarry Smith if (shell->right) { 57e4fc5df0SBarry Smith if (!shell->rightwork) { 58e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 59e4fc5df0SBarry Smith } 60e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 61e4fc5df0SBarry Smith in = shell->rightwork; 62e4fc5df0SBarry Smith } 636d7c1e57SBarry Smith while (next->next) { 646d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 652a7a6963SBarry Smith ierr = MatCreateVecs(next->mat,NULL,&next->work);CHKERRQ(ierr); 666d7c1e57SBarry Smith } 676d7c1e57SBarry Smith out = next->work; 686d7c1e57SBarry Smith ierr = MatMult(next->mat,in,out);CHKERRQ(ierr); 696d7c1e57SBarry Smith in = out; 706d7c1e57SBarry Smith next = next->next; 716d7c1e57SBarry Smith } 726d7c1e57SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 73e4fc5df0SBarry Smith if (shell->left) { 74e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 75e4fc5df0SBarry Smith } 76e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 776d7c1e57SBarry Smith PetscFunctionReturn(0); 786d7c1e57SBarry Smith } 796d7c1e57SBarry Smith 806d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 816d7c1e57SBarry Smith { 826d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 836d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 846d7c1e57SBarry Smith PetscErrorCode ierr; 856d7c1e57SBarry Smith Vec in,out; 866d7c1e57SBarry Smith 876d7c1e57SBarry Smith PetscFunctionBegin; 88e32f2f54SBarry Smith if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 896d7c1e57SBarry Smith in = x; 90e4fc5df0SBarry Smith if (shell->left) { 91e4fc5df0SBarry Smith if (!shell->leftwork) { 92e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 93e4fc5df0SBarry Smith } 94e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 95e4fc5df0SBarry Smith in = shell->leftwork; 96e4fc5df0SBarry Smith } 976d7c1e57SBarry Smith while (tail->prev) { 986d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 992a7a6963SBarry Smith ierr = MatCreateVecs(tail->mat,NULL,&tail->prev->work);CHKERRQ(ierr); 1006d7c1e57SBarry Smith } 1016d7c1e57SBarry Smith out = tail->prev->work; 1026d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1036d7c1e57SBarry Smith in = out; 1046d7c1e57SBarry Smith tail = tail->prev; 1056d7c1e57SBarry Smith } 1066d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 107e4fc5df0SBarry Smith if (shell->right) { 108e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 109e4fc5df0SBarry Smith } 110e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1116d7c1e57SBarry Smith PetscFunctionReturn(0); 1126d7c1e57SBarry Smith } 1136d7c1e57SBarry Smith 114793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 115793850ffSBarry Smith { 116793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 117793850ffSBarry Smith Mat_CompositeLink next = shell->head; 118793850ffSBarry Smith PetscErrorCode ierr; 119e4fc5df0SBarry Smith Vec in; 120793850ffSBarry Smith 121793850ffSBarry Smith PetscFunctionBegin; 122e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 123e4fc5df0SBarry Smith in = x; 124e4fc5df0SBarry Smith if (shell->right) { 125e4fc5df0SBarry Smith if (!shell->rightwork) { 126e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 127793850ffSBarry Smith } 128e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 129e4fc5df0SBarry Smith in = shell->rightwork; 130e4fc5df0SBarry Smith } 131e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 132e4fc5df0SBarry Smith while ((next = next->next)) { 133e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 134e4fc5df0SBarry Smith } 135e4fc5df0SBarry Smith if (shell->left) { 136e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 137e4fc5df0SBarry Smith } 138e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 139793850ffSBarry Smith PetscFunctionReturn(0); 140793850ffSBarry Smith } 141793850ffSBarry Smith 142793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 143793850ffSBarry Smith { 144793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 145793850ffSBarry Smith Mat_CompositeLink next = shell->head; 146793850ffSBarry Smith PetscErrorCode ierr; 147e4fc5df0SBarry Smith Vec in; 148793850ffSBarry Smith 149793850ffSBarry Smith PetscFunctionBegin; 150e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 151e4fc5df0SBarry Smith in = x; 152e4fc5df0SBarry Smith if (shell->left) { 153e4fc5df0SBarry Smith if (!shell->leftwork) { 154e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 155793850ffSBarry Smith } 156e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 157e4fc5df0SBarry Smith in = shell->leftwork; 158e4fc5df0SBarry Smith } 159e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 160e4fc5df0SBarry Smith while ((next = next->next)) { 161e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 162e4fc5df0SBarry Smith } 163e4fc5df0SBarry Smith if (shell->right) { 164e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 165e4fc5df0SBarry Smith } 166e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 167793850ffSBarry Smith PetscFunctionReturn(0); 168793850ffSBarry Smith } 169793850ffSBarry Smith 170793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 171793850ffSBarry Smith { 172793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 173793850ffSBarry Smith Mat_CompositeLink next = shell->head; 174793850ffSBarry Smith PetscErrorCode ierr; 175793850ffSBarry Smith 176793850ffSBarry Smith PetscFunctionBegin; 177e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 178e32f2f54SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 179e4fc5df0SBarry Smith 180793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 181793850ffSBarry Smith if (next->next && !shell->work) { 182793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 183793850ffSBarry Smith } 184793850ffSBarry Smith while ((next = next->next)) { 185793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 186793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 187793850ffSBarry Smith } 188e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 189793850ffSBarry Smith PetscFunctionReturn(0); 190793850ffSBarry Smith } 191793850ffSBarry Smith 192793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 193793850ffSBarry Smith { 194b52f573bSBarry Smith PetscErrorCode ierr; 195ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 196b52f573bSBarry Smith 197793850ffSBarry Smith PetscFunctionBegin; 198c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Y)->options,((PetscObject)Y)->prefix,"-mat_composite_merge",&flg,NULL);CHKERRQ(ierr); 199b52f573bSBarry Smith if (flg) { 200b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 201b52f573bSBarry Smith } 202793850ffSBarry Smith PetscFunctionReturn(0); 203793850ffSBarry Smith } 204793850ffSBarry Smith 205e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 206e4fc5df0SBarry Smith { 207e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 208e4fc5df0SBarry Smith 209e4fc5df0SBarry Smith PetscFunctionBegin; 210321429dbSBarry Smith a->scale *= alpha; 211e4fc5df0SBarry Smith PetscFunctionReturn(0); 212e4fc5df0SBarry Smith } 213e4fc5df0SBarry Smith 214e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 215e4fc5df0SBarry Smith { 216e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 217e4fc5df0SBarry Smith PetscErrorCode ierr; 218e4fc5df0SBarry Smith 219e4fc5df0SBarry Smith PetscFunctionBegin; 220e4fc5df0SBarry Smith if (left) { 221321429dbSBarry Smith if (!a->left) { 222e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 223e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 224321429dbSBarry Smith } else { 225321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 226321429dbSBarry Smith } 227e4fc5df0SBarry Smith } 228e4fc5df0SBarry Smith if (right) { 229321429dbSBarry Smith if (!a->right) { 230e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 231e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 232321429dbSBarry Smith } else { 233321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 234321429dbSBarry Smith } 235e4fc5df0SBarry Smith } 236e4fc5df0SBarry Smith PetscFunctionReturn(0); 237e4fc5df0SBarry Smith } 238793850ffSBarry Smith 239793850ffSBarry Smith static struct _MatOps MatOps_Values = {0, 240793850ffSBarry Smith 0, 241793850ffSBarry Smith 0, 242793850ffSBarry Smith MatMult_Composite, 243793850ffSBarry Smith 0, 244793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite, 245793850ffSBarry Smith 0, 246793850ffSBarry Smith 0, 247793850ffSBarry Smith 0, 248793850ffSBarry Smith 0, 249793850ffSBarry Smith /* 10*/ 0, 250793850ffSBarry Smith 0, 251793850ffSBarry Smith 0, 252793850ffSBarry Smith 0, 253793850ffSBarry Smith 0, 254793850ffSBarry Smith /* 15*/ 0, 255793850ffSBarry Smith 0, 256793850ffSBarry Smith MatGetDiagonal_Composite, 257e4fc5df0SBarry Smith MatDiagonalScale_Composite, 258793850ffSBarry Smith 0, 259793850ffSBarry Smith /* 20*/ 0, 260793850ffSBarry Smith MatAssemblyEnd_Composite, 261793850ffSBarry Smith 0, 262793850ffSBarry Smith 0, 263d519adbfSMatthew Knepley /* 24*/ 0, 264793850ffSBarry Smith 0, 265793850ffSBarry Smith 0, 266793850ffSBarry Smith 0, 267793850ffSBarry Smith 0, 268d519adbfSMatthew Knepley /* 29*/ 0, 269793850ffSBarry Smith 0, 270793850ffSBarry Smith 0, 271793850ffSBarry Smith 0, 272793850ffSBarry Smith 0, 273d519adbfSMatthew Knepley /* 34*/ 0, 274793850ffSBarry Smith 0, 275793850ffSBarry Smith 0, 276793850ffSBarry Smith 0, 277793850ffSBarry Smith 0, 278d519adbfSMatthew Knepley /* 39*/ 0, 279793850ffSBarry Smith 0, 280793850ffSBarry Smith 0, 281793850ffSBarry Smith 0, 282793850ffSBarry Smith 0, 283d519adbfSMatthew Knepley /* 44*/ 0, 284e4fc5df0SBarry Smith MatScale_Composite, 2857d68702bSBarry Smith MatShift_Basic, 286793850ffSBarry Smith 0, 287793850ffSBarry Smith 0, 288d519adbfSMatthew Knepley /* 49*/ 0, 289793850ffSBarry Smith 0, 290793850ffSBarry Smith 0, 291793850ffSBarry Smith 0, 292793850ffSBarry Smith 0, 293d519adbfSMatthew Knepley /* 54*/ 0, 294793850ffSBarry Smith 0, 295793850ffSBarry Smith 0, 296793850ffSBarry Smith 0, 297793850ffSBarry Smith 0, 298d519adbfSMatthew Knepley /* 59*/ 0, 299793850ffSBarry Smith MatDestroy_Composite, 300793850ffSBarry Smith 0, 301793850ffSBarry Smith 0, 302793850ffSBarry Smith 0, 303d519adbfSMatthew Knepley /* 64*/ 0, 304793850ffSBarry Smith 0, 305793850ffSBarry Smith 0, 306793850ffSBarry Smith 0, 307793850ffSBarry Smith 0, 308d519adbfSMatthew Knepley /* 69*/ 0, 309793850ffSBarry Smith 0, 310793850ffSBarry Smith 0, 311793850ffSBarry Smith 0, 312793850ffSBarry Smith 0, 313d519adbfSMatthew Knepley /* 74*/ 0, 314793850ffSBarry Smith 0, 315793850ffSBarry Smith 0, 316793850ffSBarry Smith 0, 317793850ffSBarry Smith 0, 318d519adbfSMatthew Knepley /* 79*/ 0, 319793850ffSBarry Smith 0, 320793850ffSBarry Smith 0, 321793850ffSBarry Smith 0, 322793850ffSBarry Smith 0, 323d519adbfSMatthew Knepley /* 84*/ 0, 324793850ffSBarry Smith 0, 325793850ffSBarry Smith 0, 326793850ffSBarry Smith 0, 327793850ffSBarry Smith 0, 328d519adbfSMatthew Knepley /* 89*/ 0, 329793850ffSBarry Smith 0, 330793850ffSBarry Smith 0, 331793850ffSBarry Smith 0, 332793850ffSBarry Smith 0, 333d519adbfSMatthew Knepley /* 94*/ 0, 334793850ffSBarry Smith 0, 335793850ffSBarry Smith 0, 3363964eb88SJed Brown 0, 3373964eb88SJed Brown 0, 3383964eb88SJed Brown /*99*/ 0, 3393964eb88SJed Brown 0, 3403964eb88SJed Brown 0, 3413964eb88SJed Brown 0, 3423964eb88SJed Brown 0, 3433964eb88SJed Brown /*104*/ 0, 3443964eb88SJed Brown 0, 3453964eb88SJed Brown 0, 3463964eb88SJed Brown 0, 3473964eb88SJed Brown 0, 3483964eb88SJed Brown /*109*/ 0, 3493964eb88SJed Brown 0, 3503964eb88SJed Brown 0, 3513964eb88SJed Brown 0, 3523964eb88SJed Brown 0, 3533964eb88SJed Brown /*114*/ 0, 3543964eb88SJed Brown 0, 3553964eb88SJed Brown 0, 3563964eb88SJed Brown 0, 3573964eb88SJed Brown 0, 3583964eb88SJed Brown /*119*/ 0, 3593964eb88SJed Brown 0, 3603964eb88SJed Brown 0, 3613964eb88SJed Brown 0, 3623964eb88SJed Brown 0, 3633964eb88SJed Brown /*124*/ 0, 3643964eb88SJed Brown 0, 3653964eb88SJed Brown 0, 3663964eb88SJed Brown 0, 3673964eb88SJed Brown 0, 3683964eb88SJed Brown /*129*/ 0, 3693964eb88SJed Brown 0, 3703964eb88SJed Brown 0, 3713964eb88SJed Brown 0, 3723964eb88SJed Brown 0, 3733964eb88SJed Brown /*134*/ 0, 3743964eb88SJed Brown 0, 3753964eb88SJed Brown 0, 3763964eb88SJed Brown 0, 3773964eb88SJed Brown 0, 3783964eb88SJed Brown /*139*/ 0, 379f9426fe0SMark Adams 0, 3803964eb88SJed Brown 0 3813964eb88SJed Brown }; 382793850ffSBarry Smith 383793850ffSBarry Smith /*MC 384c4f13f43SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 385c4f13f43SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 3866d7c1e57SBarry Smith 38795452b02SPatrick Sanan Notes: 38895452b02SPatrick Sanan to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 389793850ffSBarry Smith 390793850ffSBarry Smith Level: advanced 391793850ffSBarry Smith 3926d7c1e57SBarry Smith .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 393793850ffSBarry Smith M*/ 394793850ffSBarry Smith 3958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 396793850ffSBarry Smith { 397793850ffSBarry Smith Mat_Composite *b; 398793850ffSBarry Smith PetscErrorCode ierr; 399793850ffSBarry Smith 400793850ffSBarry Smith PetscFunctionBegin; 401b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 402793850ffSBarry Smith A->data = (void*)b; 40338f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 404793850ffSBarry Smith 40526283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 40626283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 407793850ffSBarry Smith 408793850ffSBarry Smith A->assembled = PETSC_TRUE; 4093db03f37SJed Brown A->preallocated = PETSC_TRUE; 4106d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 411e4fc5df0SBarry Smith b->scale = 1.0; 4126a7ede75SJakub Kruzik b->nmat = 0; 413793850ffSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 414*d7f81bf2SJakub Kruzik 415*d7f81bf2SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite);CHKERRQ(ierr); 416*d7f81bf2SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite);CHKERRQ(ierr); 417*d7f81bf2SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite);CHKERRQ(ierr); 418*d7f81bf2SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNMat_C",MatCompositeGetNMat_Composite);CHKERRQ(ierr); 419*d7f81bf2SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite);CHKERRQ(ierr); 420*d7f81bf2SJakub Kruzik 421*d7f81bf2SJakub Kruzik ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 422793850ffSBarry Smith PetscFunctionReturn(0); 423793850ffSBarry Smith } 424793850ffSBarry Smith 4252c0821f3SBarry Smith /*@ 426793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 427793850ffSBarry Smith 428793850ffSBarry Smith Collective on MPI_Comm 429793850ffSBarry Smith 430793850ffSBarry Smith Input Parameters: 431793850ffSBarry Smith + comm - MPI communicator 432793850ffSBarry Smith . nmat - number of matrices to put in 433793850ffSBarry Smith - mats - the matrices 434793850ffSBarry Smith 435793850ffSBarry Smith Output Parameter: 436db36af27SMatthew G. Knepley . mat - the matrix 437793850ffSBarry Smith 438793850ffSBarry Smith Level: advanced 439793850ffSBarry Smith 440793850ffSBarry Smith Notes: 441793850ffSBarry Smith Alternative construction 442793850ffSBarry Smith $ MatCreate(comm,&mat); 443793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 444793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 445793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 446793850ffSBarry Smith $ .... 447793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 448b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 449b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 450793850ffSBarry Smith 4516d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4526d7c1e57SBarry Smith 4536d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 454793850ffSBarry Smith 455793850ffSBarry Smith @*/ 4567087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 457793850ffSBarry Smith { 458793850ffSBarry Smith PetscErrorCode ierr; 459793850ffSBarry Smith PetscInt m,n,M,N,i; 460793850ffSBarry Smith 461793850ffSBarry Smith PetscFunctionBegin; 462e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 463f3f84630SBarry Smith PetscValidPointer(mat,3); 464793850ffSBarry Smith 465c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[0],PETSC_IGNORE,&n);CHKERRQ(ierr); 466c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE);CHKERRQ(ierr); 467c4f13f43SJakub Kruzik ierr = MatGetSize(mats[0],PETSC_IGNORE,&N);CHKERRQ(ierr); 468c4f13f43SJakub Kruzik ierr = MatGetSize(mats[nmat-1],&M,PETSC_IGNORE);CHKERRQ(ierr); 469793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 470793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 471793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 472793850ffSBarry Smith for (i=0; i<nmat; i++) { 473793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 474793850ffSBarry Smith } 475b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 476b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 477793850ffSBarry Smith PetscFunctionReturn(0); 478793850ffSBarry Smith } 479793850ffSBarry Smith 480*d7f81bf2SJakub Kruzik 481*d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 482*d7f81bf2SJakub Kruzik { 483*d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 484*d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 485*d7f81bf2SJakub Kruzik PetscErrorCode ierr; 486*d7f81bf2SJakub Kruzik 487*d7f81bf2SJakub Kruzik PetscFunctionBegin; 488*d7f81bf2SJakub Kruzik ierr = PetscNewLog(mat,&ilink);CHKERRQ(ierr); 489*d7f81bf2SJakub Kruzik ilink->next = 0; 490*d7f81bf2SJakub Kruzik ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 491*d7f81bf2SJakub Kruzik ilink->mat = smat; 492*d7f81bf2SJakub Kruzik 493*d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 494*d7f81bf2SJakub Kruzik else { 495*d7f81bf2SJakub Kruzik while (next->next) { 496*d7f81bf2SJakub Kruzik next = next->next; 497*d7f81bf2SJakub Kruzik } 498*d7f81bf2SJakub Kruzik next->next = ilink; 499*d7f81bf2SJakub Kruzik ilink->prev = next; 500*d7f81bf2SJakub Kruzik } 501*d7f81bf2SJakub Kruzik shell->tail = ilink; 502*d7f81bf2SJakub Kruzik shell->nmat += 1; 503*d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 504*d7f81bf2SJakub Kruzik } 505*d7f81bf2SJakub Kruzik 506793850ffSBarry Smith /*@ 507793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 508793850ffSBarry Smith 509793850ffSBarry Smith Collective on Mat 510793850ffSBarry Smith 511793850ffSBarry Smith Input Parameters: 512793850ffSBarry Smith + mat - the composite matrix 513793850ffSBarry Smith - smat - the partial matrix 514793850ffSBarry Smith 515793850ffSBarry Smith Level: advanced 516793850ffSBarry Smith 517793850ffSBarry Smith .seealso: MatCreateComposite() 518793850ffSBarry Smith @*/ 5197087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 520793850ffSBarry Smith { 521793850ffSBarry Smith PetscErrorCode ierr; 522793850ffSBarry Smith 523793850ffSBarry Smith PetscFunctionBegin; 5240700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5250700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 526*d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));CHKERRQ(ierr); 527*d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 528*d7f81bf2SJakub Kruzik } 529793850ffSBarry Smith 530*d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 531*d7f81bf2SJakub Kruzik { 532*d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 533*d7f81bf2SJakub Kruzik 534*d7f81bf2SJakub Kruzik PetscFunctionBegin; 535*d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 536*d7f81bf2SJakub Kruzik mat->ops->getdiagonal = 0; 537*d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 538*d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 539*d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_MULTIPLICATIVE; 540*d7f81bf2SJakub Kruzik } else { 541*d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 542*d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 543*d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 544*d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 545793850ffSBarry Smith } 5466d7c1e57SBarry Smith PetscFunctionReturn(0); 5476d7c1e57SBarry Smith } 5486d7c1e57SBarry Smith 5492c0821f3SBarry Smith /*@ 5506d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 5516d7c1e57SBarry Smith 5526d7c1e57SBarry Smith Collective on MPI_Comm 5536d7c1e57SBarry Smith 5546d7c1e57SBarry Smith Input Parameters: 5556d7c1e57SBarry Smith . mat - the composite matrix 5566d7c1e57SBarry Smith 5576d7c1e57SBarry Smith 5586d7c1e57SBarry Smith Level: advanced 5596d7c1e57SBarry Smith 5606d7c1e57SBarry Smith Notes: 5616d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 5626d7c1e57SBarry Smith matrix in the composite matrix. 5636d7c1e57SBarry Smith 5646d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 5656d7c1e57SBarry Smith 5666d7c1e57SBarry Smith @*/ 5677087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 5686d7c1e57SBarry Smith { 5696d7c1e57SBarry Smith PetscErrorCode ierr; 5706d7c1e57SBarry Smith 5716d7c1e57SBarry Smith PetscFunctionBegin; 572*d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 573*d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));CHKERRQ(ierr); 574793850ffSBarry Smith PetscFunctionReturn(0); 575793850ffSBarry Smith } 576793850ffSBarry Smith 577*d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 578b52f573bSBarry Smith { 579b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 5806d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 581b52f573bSBarry Smith PetscErrorCode ierr; 5826d7c1e57SBarry Smith Mat tmat,newmat; 5831ba60692SJed Brown Vec left,right; 5841ba60692SJed Brown PetscScalar scale; 585b52f573bSBarry Smith 586b52f573bSBarry Smith PetscFunctionBegin; 587e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 588b52f573bSBarry Smith 589b52f573bSBarry Smith PetscFunctionBegin; 5906d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 591b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 592b52f573bSBarry Smith while ((next = next->next)) { 593b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 594b52f573bSBarry Smith } 5956d7c1e57SBarry Smith } else { 5966d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 597b6cfcaf5SJakub Kruzik while ((next = next->next)) { 598b6cfcaf5SJakub Kruzik ierr = MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 5996bf464f9SBarry Smith ierr = MatDestroy(&tmat);CHKERRQ(ierr); 6006d7c1e57SBarry Smith tmat = newmat; 6016d7c1e57SBarry Smith } 6026d7c1e57SBarry Smith } 6031ba60692SJed Brown 6041ba60692SJed Brown scale = shell->scale; 6051ba60692SJed Brown if ((left = shell->left)) {ierr = PetscObjectReference((PetscObject)left);CHKERRQ(ierr);} 6061ba60692SJed Brown if ((right = shell->right)) {ierr = PetscObjectReference((PetscObject)right);CHKERRQ(ierr);} 6071ba60692SJed Brown 60828be2f97SBarry Smith ierr = MatHeaderReplace(mat,&tmat);CHKERRQ(ierr); 6091ba60692SJed Brown 6101ba60692SJed Brown ierr = MatDiagonalScale(mat,left,right);CHKERRQ(ierr); 6111ba60692SJed Brown ierr = MatScale(mat,scale);CHKERRQ(ierr); 6121ba60692SJed Brown ierr = VecDestroy(&left);CHKERRQ(ierr); 6131ba60692SJed Brown ierr = VecDestroy(&right);CHKERRQ(ierr); 614b52f573bSBarry Smith PetscFunctionReturn(0); 615b52f573bSBarry Smith } 6166a7ede75SJakub Kruzik 6176a7ede75SJakub Kruzik /*@ 618*d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 619*d7f81bf2SJakub Kruzik by summing all the matrices inside the composite matrix. 620*d7f81bf2SJakub Kruzik 621*d7f81bf2SJakub Kruzik Collective on MPI_Comm 622*d7f81bf2SJakub Kruzik 623*d7f81bf2SJakub Kruzik Input Parameters: 624*d7f81bf2SJakub Kruzik . mat - the composite matrix 625*d7f81bf2SJakub Kruzik 626*d7f81bf2SJakub Kruzik 627*d7f81bf2SJakub Kruzik Options Database: 628*d7f81bf2SJakub Kruzik . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 629*d7f81bf2SJakub Kruzik 630*d7f81bf2SJakub Kruzik Level: advanced 631*d7f81bf2SJakub Kruzik 632*d7f81bf2SJakub Kruzik Notes: 633*d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 634*d7f81bf2SJakub Kruzik matrix in the composite matrix. 635*d7f81bf2SJakub Kruzik 636*d7f81bf2SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 637*d7f81bf2SJakub Kruzik 638*d7f81bf2SJakub Kruzik @*/ 639*d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 640*d7f81bf2SJakub Kruzik { 641*d7f81bf2SJakub Kruzik PetscErrorCode ierr; 642*d7f81bf2SJakub Kruzik 643*d7f81bf2SJakub Kruzik PetscFunctionBegin; 644*d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 645*d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));CHKERRQ(ierr); 646*d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 647*d7f81bf2SJakub Kruzik } 648*d7f81bf2SJakub Kruzik 649*d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetNMat_Composite(Mat mat,PetscInt *nmat) 650*d7f81bf2SJakub Kruzik { 651*d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 652*d7f81bf2SJakub Kruzik 653*d7f81bf2SJakub Kruzik PetscFunctionBegin; 654*d7f81bf2SJakub Kruzik *nmat = shell->nmat; 655*d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 656*d7f81bf2SJakub Kruzik } 657*d7f81bf2SJakub Kruzik 658*d7f81bf2SJakub Kruzik /*@ 6596a7ede75SJakub Kruzik MatCompositeGetNMat - Returns the number of matrices in composite. 6606a7ede75SJakub Kruzik 6616a7ede75SJakub Kruzik Not Collective 6626a7ede75SJakub Kruzik 6636a7ede75SJakub Kruzik Input Parameter: 664*d7f81bf2SJakub Kruzik . mat - the composite matrix 6656a7ede75SJakub Kruzik 6666a7ede75SJakub Kruzik Output Parameter: 6676a7ede75SJakub Kruzik . size - the local size 6686a7ede75SJakub Kruzik . nmat - number of matrices in composite 6696a7ede75SJakub Kruzik 6708b5c3584SJakub Kruzik Level: advanced 6716a7ede75SJakub Kruzik 6728b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat() 6736a7ede75SJakub Kruzik 6746a7ede75SJakub Kruzik @*/ 675*d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetNMat(Mat mat,PetscInt *nmat) 6766a7ede75SJakub Kruzik { 677*d7f81bf2SJakub Kruzik PetscErrorCode ierr; 6786a7ede75SJakub Kruzik 6796a7ede75SJakub Kruzik PetscFunctionBegin; 680*d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6816a7ede75SJakub Kruzik PetscValidPointer(nmat,2); 682*d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetNMat_C",(Mat,PetscInt*),(mat,nmat));CHKERRQ(ierr); 683*d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 684*d7f81bf2SJakub Kruzik } 685*d7f81bf2SJakub Kruzik 686*d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 687*d7f81bf2SJakub Kruzik { 688*d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 689*d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 690*d7f81bf2SJakub Kruzik PetscInt k; 691*d7f81bf2SJakub Kruzik 692*d7f81bf2SJakub Kruzik PetscFunctionBegin; 693*d7f81bf2SJakub Kruzik if (i >= shell->nmat) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %d >= %d",i,shell->nmat); 694*d7f81bf2SJakub Kruzik ilink = shell->head; 695*d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 696*d7f81bf2SJakub Kruzik ilink = ilink->next; 697*d7f81bf2SJakub Kruzik } 698*d7f81bf2SJakub Kruzik *Ai = ilink->mat; 6996a7ede75SJakub Kruzik PetscFunctionReturn(0); 7006a7ede75SJakub Kruzik } 7016a7ede75SJakub Kruzik 7028b5c3584SJakub Kruzik /*@ 7038b5c3584SJakub Kruzik MatCompositeGetMat - Returns the ith matrix from composite. 7048b5c3584SJakub Kruzik 7058b5c3584SJakub Kruzik Logically Collective on Mat 7068b5c3584SJakub Kruzik 7078b5c3584SJakub Kruzik Input Parameter: 708*d7f81bf2SJakub Kruzik + mat - the composite matrix 7098b5c3584SJakub Kruzik - i - the number of requested matrix 7108b5c3584SJakub Kruzik 7118b5c3584SJakub Kruzik Output Parameter: 7128b5c3584SJakub Kruzik . Ai - ith matrix in composite 7138b5c3584SJakub Kruzik 7148b5c3584SJakub Kruzik Level: advanced 7158b5c3584SJakub Kruzik 7168b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNMat() 7178b5c3584SJakub Kruzik 7188b5c3584SJakub Kruzik @*/ 719*d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 7208b5c3584SJakub Kruzik { 721*d7f81bf2SJakub Kruzik PetscErrorCode ierr; 7228b5c3584SJakub Kruzik 7238b5c3584SJakub Kruzik PetscFunctionBegin; 724*d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 725*d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 7268b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 727*d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));CHKERRQ(ierr); 7288b5c3584SJakub Kruzik PetscFunctionReturn(0); 7298b5c3584SJakub Kruzik } 7308b5c3584SJakub Kruzik 731