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 2392c0821f3SBarry Smith /*@ 240793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 241793850ffSBarry Smith 242793850ffSBarry Smith Collective on MPI_Comm 243793850ffSBarry Smith 244793850ffSBarry Smith Input Parameters: 245793850ffSBarry Smith + comm - MPI communicator 246793850ffSBarry Smith . nmat - number of matrices to put in 247793850ffSBarry Smith - mats - the matrices 248793850ffSBarry Smith 249793850ffSBarry Smith Output Parameter: 250db36af27SMatthew G. Knepley . mat - the matrix 251793850ffSBarry Smith 252793850ffSBarry Smith Level: advanced 253793850ffSBarry Smith 254793850ffSBarry Smith Notes: 255793850ffSBarry Smith Alternative construction 256793850ffSBarry Smith $ MatCreate(comm,&mat); 257793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 258793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 259793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 260793850ffSBarry Smith $ .... 261793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 262b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 263b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 264793850ffSBarry Smith 2656d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 2666d7c1e57SBarry Smith 2676d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 268793850ffSBarry Smith 269793850ffSBarry Smith @*/ 2707087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 271793850ffSBarry Smith { 272793850ffSBarry Smith PetscErrorCode ierr; 273793850ffSBarry Smith PetscInt m,n,M,N,i; 274793850ffSBarry Smith 275793850ffSBarry Smith PetscFunctionBegin; 276e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 277f3f84630SBarry Smith PetscValidPointer(mat,3); 278793850ffSBarry Smith 279c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[0],PETSC_IGNORE,&n);CHKERRQ(ierr); 280c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE);CHKERRQ(ierr); 281c4f13f43SJakub Kruzik ierr = MatGetSize(mats[0],PETSC_IGNORE,&N);CHKERRQ(ierr); 282c4f13f43SJakub Kruzik ierr = MatGetSize(mats[nmat-1],&M,PETSC_IGNORE);CHKERRQ(ierr); 283793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 284793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 285793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 286793850ffSBarry Smith for (i=0; i<nmat; i++) { 287793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 288793850ffSBarry Smith } 289b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 290b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 291793850ffSBarry Smith PetscFunctionReturn(0); 292793850ffSBarry Smith } 293793850ffSBarry Smith 294d7f81bf2SJakub Kruzik 295d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 296d7f81bf2SJakub Kruzik { 297d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 298d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 299d7f81bf2SJakub Kruzik PetscErrorCode ierr; 300d7f81bf2SJakub Kruzik 301d7f81bf2SJakub Kruzik PetscFunctionBegin; 302d7f81bf2SJakub Kruzik ierr = PetscNewLog(mat,&ilink);CHKERRQ(ierr); 303d7f81bf2SJakub Kruzik ilink->next = 0; 304d7f81bf2SJakub Kruzik ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 305d7f81bf2SJakub Kruzik ilink->mat = smat; 306d7f81bf2SJakub Kruzik 307d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 308d7f81bf2SJakub Kruzik else { 309d7f81bf2SJakub Kruzik while (next->next) { 310d7f81bf2SJakub Kruzik next = next->next; 311d7f81bf2SJakub Kruzik } 312d7f81bf2SJakub Kruzik next->next = ilink; 313d7f81bf2SJakub Kruzik ilink->prev = next; 314d7f81bf2SJakub Kruzik } 315d7f81bf2SJakub Kruzik shell->tail = ilink; 316d7f81bf2SJakub Kruzik shell->nmat += 1; 317d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 318d7f81bf2SJakub Kruzik } 319d7f81bf2SJakub Kruzik 320793850ffSBarry Smith /*@ 321793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 322793850ffSBarry Smith 323793850ffSBarry Smith Collective on Mat 324793850ffSBarry Smith 325793850ffSBarry Smith Input Parameters: 326793850ffSBarry Smith + mat - the composite matrix 327793850ffSBarry Smith - smat - the partial matrix 328793850ffSBarry Smith 329793850ffSBarry Smith Level: advanced 330793850ffSBarry Smith 331793850ffSBarry Smith .seealso: MatCreateComposite() 332793850ffSBarry Smith @*/ 3337087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 334793850ffSBarry Smith { 335793850ffSBarry Smith PetscErrorCode ierr; 336793850ffSBarry Smith 337793850ffSBarry Smith PetscFunctionBegin; 3380700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 3390700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 340d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));CHKERRQ(ierr); 341d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 342d7f81bf2SJakub Kruzik } 343793850ffSBarry Smith 344d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 345d7f81bf2SJakub Kruzik { 346d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 347d7f81bf2SJakub Kruzik 348d7f81bf2SJakub Kruzik PetscFunctionBegin; 349d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 350d7f81bf2SJakub Kruzik mat->ops->getdiagonal = 0; 351d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 352d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 353d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_MULTIPLICATIVE; 354d7f81bf2SJakub Kruzik } else { 355d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 356d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 357d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 358d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 359793850ffSBarry Smith } 3606d7c1e57SBarry Smith PetscFunctionReturn(0); 3616d7c1e57SBarry Smith } 3626d7c1e57SBarry Smith 3632c0821f3SBarry Smith /*@ 3646d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 3656d7c1e57SBarry Smith 3666d7c1e57SBarry Smith Collective on MPI_Comm 3676d7c1e57SBarry Smith 3686d7c1e57SBarry Smith Input Parameters: 3696d7c1e57SBarry Smith . mat - the composite matrix 3706d7c1e57SBarry Smith 3716d7c1e57SBarry Smith 3726d7c1e57SBarry Smith Level: advanced 3736d7c1e57SBarry Smith 3746d7c1e57SBarry Smith Notes: 3756d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 3766d7c1e57SBarry Smith matrix in the composite matrix. 3776d7c1e57SBarry Smith 3786d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 3796d7c1e57SBarry Smith 3806d7c1e57SBarry Smith @*/ 3817087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 3826d7c1e57SBarry Smith { 3836d7c1e57SBarry Smith PetscErrorCode ierr; 3846d7c1e57SBarry Smith 3856d7c1e57SBarry Smith PetscFunctionBegin; 386d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 387d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));CHKERRQ(ierr); 388793850ffSBarry Smith PetscFunctionReturn(0); 389793850ffSBarry Smith } 390793850ffSBarry Smith 391d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 392b52f573bSBarry Smith { 393b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 3946d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 395b52f573bSBarry Smith PetscErrorCode ierr; 3966d7c1e57SBarry Smith Mat tmat,newmat; 3971ba60692SJed Brown Vec left,right; 3981ba60692SJed Brown PetscScalar scale; 399b52f573bSBarry Smith 400b52f573bSBarry Smith PetscFunctionBegin; 401e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 402b52f573bSBarry Smith 403b52f573bSBarry Smith PetscFunctionBegin; 4046d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 405b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 406b52f573bSBarry Smith while ((next = next->next)) { 407b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 408b52f573bSBarry Smith } 4096d7c1e57SBarry Smith } else { 4106d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 411b6cfcaf5SJakub Kruzik while ((next = next->next)) { 412b6cfcaf5SJakub Kruzik ierr = MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 4136bf464f9SBarry Smith ierr = MatDestroy(&tmat);CHKERRQ(ierr); 4146d7c1e57SBarry Smith tmat = newmat; 4156d7c1e57SBarry Smith } 4166d7c1e57SBarry Smith } 4171ba60692SJed Brown 4181ba60692SJed Brown scale = shell->scale; 4191ba60692SJed Brown if ((left = shell->left)) {ierr = PetscObjectReference((PetscObject)left);CHKERRQ(ierr);} 4201ba60692SJed Brown if ((right = shell->right)) {ierr = PetscObjectReference((PetscObject)right);CHKERRQ(ierr);} 4211ba60692SJed Brown 42228be2f97SBarry Smith ierr = MatHeaderReplace(mat,&tmat);CHKERRQ(ierr); 4231ba60692SJed Brown 4241ba60692SJed Brown ierr = MatDiagonalScale(mat,left,right);CHKERRQ(ierr); 4251ba60692SJed Brown ierr = MatScale(mat,scale);CHKERRQ(ierr); 4261ba60692SJed Brown ierr = VecDestroy(&left);CHKERRQ(ierr); 4271ba60692SJed Brown ierr = VecDestroy(&right);CHKERRQ(ierr); 428b52f573bSBarry Smith PetscFunctionReturn(0); 429b52f573bSBarry Smith } 4306a7ede75SJakub Kruzik 4316a7ede75SJakub Kruzik /*@ 432d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 433d7f81bf2SJakub Kruzik by summing all the matrices inside the composite matrix. 434d7f81bf2SJakub Kruzik 435d7f81bf2SJakub Kruzik Collective on MPI_Comm 436d7f81bf2SJakub Kruzik 437d7f81bf2SJakub Kruzik Input Parameters: 438d7f81bf2SJakub Kruzik . mat - the composite matrix 439d7f81bf2SJakub Kruzik 440d7f81bf2SJakub Kruzik 441d7f81bf2SJakub Kruzik Options Database: 442d7f81bf2SJakub Kruzik . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 443d7f81bf2SJakub Kruzik 444d7f81bf2SJakub Kruzik Level: advanced 445d7f81bf2SJakub Kruzik 446d7f81bf2SJakub Kruzik Notes: 447d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 448d7f81bf2SJakub Kruzik matrix in the composite matrix. 449d7f81bf2SJakub Kruzik 450d7f81bf2SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 451d7f81bf2SJakub Kruzik 452d7f81bf2SJakub Kruzik @*/ 453d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 454d7f81bf2SJakub Kruzik { 455d7f81bf2SJakub Kruzik PetscErrorCode ierr; 456d7f81bf2SJakub Kruzik 457d7f81bf2SJakub Kruzik PetscFunctionBegin; 458d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 459d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));CHKERRQ(ierr); 460d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 461d7f81bf2SJakub Kruzik } 462d7f81bf2SJakub Kruzik 463d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetNMat_Composite(Mat mat,PetscInt *nmat) 464d7f81bf2SJakub Kruzik { 465d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 466d7f81bf2SJakub Kruzik 467d7f81bf2SJakub Kruzik PetscFunctionBegin; 468d7f81bf2SJakub Kruzik *nmat = shell->nmat; 469d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 470d7f81bf2SJakub Kruzik } 471d7f81bf2SJakub Kruzik 472d7f81bf2SJakub Kruzik /*@ 4736a7ede75SJakub Kruzik MatCompositeGetNMat - Returns the number of matrices in composite. 4746a7ede75SJakub Kruzik 4756a7ede75SJakub Kruzik Not Collective 4766a7ede75SJakub Kruzik 4776a7ede75SJakub Kruzik Input Parameter: 478d7f81bf2SJakub Kruzik . mat - the composite matrix 4796a7ede75SJakub Kruzik 4806a7ede75SJakub Kruzik Output Parameter: 4816a7ede75SJakub Kruzik . size - the local size 4826a7ede75SJakub Kruzik . nmat - number of matrices in composite 4836a7ede75SJakub Kruzik 4848b5c3584SJakub Kruzik Level: advanced 4856a7ede75SJakub Kruzik 4868b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat() 4876a7ede75SJakub Kruzik 4886a7ede75SJakub Kruzik @*/ 489d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetNMat(Mat mat,PetscInt *nmat) 4906a7ede75SJakub Kruzik { 491d7f81bf2SJakub Kruzik PetscErrorCode ierr; 4926a7ede75SJakub Kruzik 4936a7ede75SJakub Kruzik PetscFunctionBegin; 494d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4956a7ede75SJakub Kruzik PetscValidPointer(nmat,2); 496d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetNMat_C",(Mat,PetscInt*),(mat,nmat));CHKERRQ(ierr); 497d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 498d7f81bf2SJakub Kruzik } 499d7f81bf2SJakub Kruzik 500d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 501d7f81bf2SJakub Kruzik { 502d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 503d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 504d7f81bf2SJakub Kruzik PetscInt k; 505d7f81bf2SJakub Kruzik 506d7f81bf2SJakub Kruzik PetscFunctionBegin; 507d7f81bf2SJakub Kruzik if (i >= shell->nmat) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %d >= %d",i,shell->nmat); 508d7f81bf2SJakub Kruzik ilink = shell->head; 509d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 510d7f81bf2SJakub Kruzik ilink = ilink->next; 511d7f81bf2SJakub Kruzik } 512d7f81bf2SJakub Kruzik *Ai = ilink->mat; 5136a7ede75SJakub Kruzik PetscFunctionReturn(0); 5146a7ede75SJakub Kruzik } 5156a7ede75SJakub Kruzik 5168b5c3584SJakub Kruzik /*@ 5178b5c3584SJakub Kruzik MatCompositeGetMat - Returns the ith matrix from composite. 5188b5c3584SJakub Kruzik 5198b5c3584SJakub Kruzik Logically Collective on Mat 5208b5c3584SJakub Kruzik 5218b5c3584SJakub Kruzik Input Parameter: 522d7f81bf2SJakub Kruzik + mat - the composite matrix 5238b5c3584SJakub Kruzik - i - the number of requested matrix 5248b5c3584SJakub Kruzik 5258b5c3584SJakub Kruzik Output Parameter: 5268b5c3584SJakub Kruzik . Ai - ith matrix in composite 5278b5c3584SJakub Kruzik 5288b5c3584SJakub Kruzik Level: advanced 5298b5c3584SJakub Kruzik 5308b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNMat() 5318b5c3584SJakub Kruzik 5328b5c3584SJakub Kruzik @*/ 533d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 5348b5c3584SJakub Kruzik { 535d7f81bf2SJakub Kruzik PetscErrorCode ierr; 5368b5c3584SJakub Kruzik 5378b5c3584SJakub Kruzik PetscFunctionBegin; 538d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 539d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 5408b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 541d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));CHKERRQ(ierr); 5428b5c3584SJakub Kruzik PetscFunctionReturn(0); 5438b5c3584SJakub Kruzik } 5448b5c3584SJakub Kruzik 545*41cd0125SJakub Kruzik static struct _MatOps MatOps_Values = {0, 546*41cd0125SJakub Kruzik 0, 547*41cd0125SJakub Kruzik 0, 548*41cd0125SJakub Kruzik MatMult_Composite, 549*41cd0125SJakub Kruzik 0, 550*41cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 551*41cd0125SJakub Kruzik 0, 552*41cd0125SJakub Kruzik 0, 553*41cd0125SJakub Kruzik 0, 554*41cd0125SJakub Kruzik 0, 555*41cd0125SJakub Kruzik /* 10*/ 0, 556*41cd0125SJakub Kruzik 0, 557*41cd0125SJakub Kruzik 0, 558*41cd0125SJakub Kruzik 0, 559*41cd0125SJakub Kruzik 0, 560*41cd0125SJakub Kruzik /* 15*/ 0, 561*41cd0125SJakub Kruzik 0, 562*41cd0125SJakub Kruzik MatGetDiagonal_Composite, 563*41cd0125SJakub Kruzik MatDiagonalScale_Composite, 564*41cd0125SJakub Kruzik 0, 565*41cd0125SJakub Kruzik /* 20*/ 0, 566*41cd0125SJakub Kruzik MatAssemblyEnd_Composite, 567*41cd0125SJakub Kruzik 0, 568*41cd0125SJakub Kruzik 0, 569*41cd0125SJakub Kruzik /* 24*/ 0, 570*41cd0125SJakub Kruzik 0, 571*41cd0125SJakub Kruzik 0, 572*41cd0125SJakub Kruzik 0, 573*41cd0125SJakub Kruzik 0, 574*41cd0125SJakub Kruzik /* 29*/ 0, 575*41cd0125SJakub Kruzik 0, 576*41cd0125SJakub Kruzik 0, 577*41cd0125SJakub Kruzik 0, 578*41cd0125SJakub Kruzik 0, 579*41cd0125SJakub Kruzik /* 34*/ 0, 580*41cd0125SJakub Kruzik 0, 581*41cd0125SJakub Kruzik 0, 582*41cd0125SJakub Kruzik 0, 583*41cd0125SJakub Kruzik 0, 584*41cd0125SJakub Kruzik /* 39*/ 0, 585*41cd0125SJakub Kruzik 0, 586*41cd0125SJakub Kruzik 0, 587*41cd0125SJakub Kruzik 0, 588*41cd0125SJakub Kruzik 0, 589*41cd0125SJakub Kruzik /* 44*/ 0, 590*41cd0125SJakub Kruzik MatScale_Composite, 591*41cd0125SJakub Kruzik MatShift_Basic, 592*41cd0125SJakub Kruzik 0, 593*41cd0125SJakub Kruzik 0, 594*41cd0125SJakub Kruzik /* 49*/ 0, 595*41cd0125SJakub Kruzik 0, 596*41cd0125SJakub Kruzik 0, 597*41cd0125SJakub Kruzik 0, 598*41cd0125SJakub Kruzik 0, 599*41cd0125SJakub Kruzik /* 54*/ 0, 600*41cd0125SJakub Kruzik 0, 601*41cd0125SJakub Kruzik 0, 602*41cd0125SJakub Kruzik 0, 603*41cd0125SJakub Kruzik 0, 604*41cd0125SJakub Kruzik /* 59*/ 0, 605*41cd0125SJakub Kruzik MatDestroy_Composite, 606*41cd0125SJakub Kruzik 0, 607*41cd0125SJakub Kruzik 0, 608*41cd0125SJakub Kruzik 0, 609*41cd0125SJakub Kruzik /* 64*/ 0, 610*41cd0125SJakub Kruzik 0, 611*41cd0125SJakub Kruzik 0, 612*41cd0125SJakub Kruzik 0, 613*41cd0125SJakub Kruzik 0, 614*41cd0125SJakub Kruzik /* 69*/ 0, 615*41cd0125SJakub Kruzik 0, 616*41cd0125SJakub Kruzik 0, 617*41cd0125SJakub Kruzik 0, 618*41cd0125SJakub Kruzik 0, 619*41cd0125SJakub Kruzik /* 74*/ 0, 620*41cd0125SJakub Kruzik 0, 621*41cd0125SJakub Kruzik 0, 622*41cd0125SJakub Kruzik 0, 623*41cd0125SJakub Kruzik 0, 624*41cd0125SJakub Kruzik /* 79*/ 0, 625*41cd0125SJakub Kruzik 0, 626*41cd0125SJakub Kruzik 0, 627*41cd0125SJakub Kruzik 0, 628*41cd0125SJakub Kruzik 0, 629*41cd0125SJakub Kruzik /* 84*/ 0, 630*41cd0125SJakub Kruzik 0, 631*41cd0125SJakub Kruzik 0, 632*41cd0125SJakub Kruzik 0, 633*41cd0125SJakub Kruzik 0, 634*41cd0125SJakub Kruzik /* 89*/ 0, 635*41cd0125SJakub Kruzik 0, 636*41cd0125SJakub Kruzik 0, 637*41cd0125SJakub Kruzik 0, 638*41cd0125SJakub Kruzik 0, 639*41cd0125SJakub Kruzik /* 94*/ 0, 640*41cd0125SJakub Kruzik 0, 641*41cd0125SJakub Kruzik 0, 642*41cd0125SJakub Kruzik 0, 643*41cd0125SJakub Kruzik 0, 644*41cd0125SJakub Kruzik /*99*/ 0, 645*41cd0125SJakub Kruzik 0, 646*41cd0125SJakub Kruzik 0, 647*41cd0125SJakub Kruzik 0, 648*41cd0125SJakub Kruzik 0, 649*41cd0125SJakub Kruzik /*104*/ 0, 650*41cd0125SJakub Kruzik 0, 651*41cd0125SJakub Kruzik 0, 652*41cd0125SJakub Kruzik 0, 653*41cd0125SJakub Kruzik 0, 654*41cd0125SJakub Kruzik /*109*/ 0, 655*41cd0125SJakub Kruzik 0, 656*41cd0125SJakub Kruzik 0, 657*41cd0125SJakub Kruzik 0, 658*41cd0125SJakub Kruzik 0, 659*41cd0125SJakub Kruzik /*114*/ 0, 660*41cd0125SJakub Kruzik 0, 661*41cd0125SJakub Kruzik 0, 662*41cd0125SJakub Kruzik 0, 663*41cd0125SJakub Kruzik 0, 664*41cd0125SJakub Kruzik /*119*/ 0, 665*41cd0125SJakub Kruzik 0, 666*41cd0125SJakub Kruzik 0, 667*41cd0125SJakub Kruzik 0, 668*41cd0125SJakub Kruzik 0, 669*41cd0125SJakub Kruzik /*124*/ 0, 670*41cd0125SJakub Kruzik 0, 671*41cd0125SJakub Kruzik 0, 672*41cd0125SJakub Kruzik 0, 673*41cd0125SJakub Kruzik 0, 674*41cd0125SJakub Kruzik /*129*/ 0, 675*41cd0125SJakub Kruzik 0, 676*41cd0125SJakub Kruzik 0, 677*41cd0125SJakub Kruzik 0, 678*41cd0125SJakub Kruzik 0, 679*41cd0125SJakub Kruzik /*134*/ 0, 680*41cd0125SJakub Kruzik 0, 681*41cd0125SJakub Kruzik 0, 682*41cd0125SJakub Kruzik 0, 683*41cd0125SJakub Kruzik 0, 684*41cd0125SJakub Kruzik /*139*/ 0, 685*41cd0125SJakub Kruzik 0, 686*41cd0125SJakub Kruzik 0 687*41cd0125SJakub Kruzik }; 688*41cd0125SJakub Kruzik 689*41cd0125SJakub Kruzik /*MC 690*41cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 691*41cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 692*41cd0125SJakub Kruzik 693*41cd0125SJakub Kruzik Notes: 694*41cd0125SJakub Kruzik to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 695*41cd0125SJakub Kruzik 696*41cd0125SJakub Kruzik Level: advanced 697*41cd0125SJakub Kruzik 698*41cd0125SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 699*41cd0125SJakub Kruzik M*/ 700*41cd0125SJakub Kruzik 701*41cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 702*41cd0125SJakub Kruzik { 703*41cd0125SJakub Kruzik Mat_Composite *b; 704*41cd0125SJakub Kruzik PetscErrorCode ierr; 705*41cd0125SJakub Kruzik 706*41cd0125SJakub Kruzik PetscFunctionBegin; 707*41cd0125SJakub Kruzik ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 708*41cd0125SJakub Kruzik A->data = (void*)b; 709*41cd0125SJakub Kruzik ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 710*41cd0125SJakub Kruzik 711*41cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 712*41cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 713*41cd0125SJakub Kruzik 714*41cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 715*41cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 716*41cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 717*41cd0125SJakub Kruzik b->scale = 1.0; 718*41cd0125SJakub Kruzik b->nmat = 0; 719*41cd0125SJakub Kruzik 720*41cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite);CHKERRQ(ierr); 721*41cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite);CHKERRQ(ierr); 722*41cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite);CHKERRQ(ierr); 723*41cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNMat_C",MatCompositeGetNMat_Composite);CHKERRQ(ierr); 724*41cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite);CHKERRQ(ierr); 725*41cd0125SJakub Kruzik 726*41cd0125SJakub Kruzik ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 727*41cd0125SJakub Kruzik PetscFunctionReturn(0); 728*41cd0125SJakub Kruzik } 729*41cd0125SJakub Kruzik 730