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; 1904d534ceSJakub Kruzik PetscBool mergefromright; 20793850ffSBarry Smith } Mat_Composite; 21793850ffSBarry Smith 22793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat) 23793850ffSBarry Smith { 24793850ffSBarry Smith PetscErrorCode ierr; 252c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite*)mat->data; 266d7c1e57SBarry Smith Mat_CompositeLink next = shell->head,oldnext; 27793850ffSBarry Smith 28793850ffSBarry Smith PetscFunctionBegin; 29793850ffSBarry Smith while (next) { 306bf464f9SBarry Smith ierr = MatDestroy(&next->mat);CHKERRQ(ierr); 316d7c1e57SBarry Smith if (next->work && (!next->next || next->work != next->next->work)) { 326bf464f9SBarry Smith ierr = VecDestroy(&next->work);CHKERRQ(ierr); 336d7c1e57SBarry Smith } 346d7c1e57SBarry Smith oldnext = next; 35793850ffSBarry Smith next = next->next; 366d7c1e57SBarry Smith ierr = PetscFree(oldnext);CHKERRQ(ierr); 37793850ffSBarry Smith } 386bf464f9SBarry Smith ierr = VecDestroy(&shell->work);CHKERRQ(ierr); 396bf464f9SBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 406bf464f9SBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 416bf464f9SBarry Smith ierr = VecDestroy(&shell->leftwork);CHKERRQ(ierr); 426bf464f9SBarry Smith ierr = VecDestroy(&shell->rightwork);CHKERRQ(ierr); 436bf464f9SBarry Smith ierr = PetscFree(mat->data);CHKERRQ(ierr); 44793850ffSBarry Smith PetscFunctionReturn(0); 45793850ffSBarry Smith } 46793850ffSBarry Smith 476d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 486d7c1e57SBarry Smith { 496d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 506d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 516d7c1e57SBarry Smith PetscErrorCode ierr; 526d7c1e57SBarry Smith Vec in,out; 536d7c1e57SBarry Smith 546d7c1e57SBarry Smith PetscFunctionBegin; 55e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 566d7c1e57SBarry Smith in = x; 57e4fc5df0SBarry Smith if (shell->right) { 58e4fc5df0SBarry Smith if (!shell->rightwork) { 59e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 60e4fc5df0SBarry Smith } 61e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 62e4fc5df0SBarry Smith in = shell->rightwork; 63e4fc5df0SBarry Smith } 646d7c1e57SBarry Smith while (next->next) { 656d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 662a7a6963SBarry Smith ierr = MatCreateVecs(next->mat,NULL,&next->work);CHKERRQ(ierr); 676d7c1e57SBarry Smith } 686d7c1e57SBarry Smith out = next->work; 696d7c1e57SBarry Smith ierr = MatMult(next->mat,in,out);CHKERRQ(ierr); 706d7c1e57SBarry Smith in = out; 716d7c1e57SBarry Smith next = next->next; 726d7c1e57SBarry Smith } 736d7c1e57SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 74e4fc5df0SBarry Smith if (shell->left) { 75e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 76e4fc5df0SBarry Smith } 77e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 786d7c1e57SBarry Smith PetscFunctionReturn(0); 796d7c1e57SBarry Smith } 806d7c1e57SBarry Smith 816d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 826d7c1e57SBarry Smith { 836d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 846d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 856d7c1e57SBarry Smith PetscErrorCode ierr; 866d7c1e57SBarry Smith Vec in,out; 876d7c1e57SBarry Smith 886d7c1e57SBarry Smith PetscFunctionBegin; 89e32f2f54SBarry Smith if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 906d7c1e57SBarry Smith in = x; 91e4fc5df0SBarry Smith if (shell->left) { 92e4fc5df0SBarry Smith if (!shell->leftwork) { 93e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 94e4fc5df0SBarry Smith } 95e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 96e4fc5df0SBarry Smith in = shell->leftwork; 97e4fc5df0SBarry Smith } 986d7c1e57SBarry Smith while (tail->prev) { 996d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1002a7a6963SBarry Smith ierr = MatCreateVecs(tail->mat,NULL,&tail->prev->work);CHKERRQ(ierr); 1016d7c1e57SBarry Smith } 1026d7c1e57SBarry Smith out = tail->prev->work; 1036d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1046d7c1e57SBarry Smith in = out; 1056d7c1e57SBarry Smith tail = tail->prev; 1066d7c1e57SBarry Smith } 1076d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 108e4fc5df0SBarry Smith if (shell->right) { 109e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 110e4fc5df0SBarry Smith } 111e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1126d7c1e57SBarry Smith PetscFunctionReturn(0); 1136d7c1e57SBarry Smith } 1146d7c1e57SBarry Smith 115793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 116793850ffSBarry Smith { 117793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 118793850ffSBarry Smith Mat_CompositeLink next = shell->head; 119793850ffSBarry Smith PetscErrorCode ierr; 120e4fc5df0SBarry Smith Vec in; 121793850ffSBarry Smith 122793850ffSBarry Smith PetscFunctionBegin; 123e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 124e4fc5df0SBarry Smith in = x; 125e4fc5df0SBarry Smith if (shell->right) { 126e4fc5df0SBarry Smith if (!shell->rightwork) { 127e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 128793850ffSBarry Smith } 129e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 130e4fc5df0SBarry Smith in = shell->rightwork; 131e4fc5df0SBarry Smith } 132e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 133e4fc5df0SBarry Smith while ((next = next->next)) { 134e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 135e4fc5df0SBarry Smith } 136e4fc5df0SBarry Smith if (shell->left) { 137e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 138e4fc5df0SBarry Smith } 139e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 140793850ffSBarry Smith PetscFunctionReturn(0); 141793850ffSBarry Smith } 142793850ffSBarry Smith 143793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 144793850ffSBarry Smith { 145793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 146793850ffSBarry Smith Mat_CompositeLink next = shell->head; 147793850ffSBarry Smith PetscErrorCode ierr; 148e4fc5df0SBarry Smith Vec in; 149793850ffSBarry Smith 150793850ffSBarry Smith PetscFunctionBegin; 151e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 152e4fc5df0SBarry Smith in = x; 153e4fc5df0SBarry Smith if (shell->left) { 154e4fc5df0SBarry Smith if (!shell->leftwork) { 155e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 156793850ffSBarry Smith } 157e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 158e4fc5df0SBarry Smith in = shell->leftwork; 159e4fc5df0SBarry Smith } 160e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 161e4fc5df0SBarry Smith while ((next = next->next)) { 162e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 163e4fc5df0SBarry Smith } 164e4fc5df0SBarry Smith if (shell->right) { 165e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 166e4fc5df0SBarry Smith } 167e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 168793850ffSBarry Smith PetscFunctionReturn(0); 169793850ffSBarry Smith } 170793850ffSBarry Smith 171*7bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z) 172*7bf3a71bSJakub Kruzik { 173*7bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 174*7bf3a71bSJakub Kruzik PetscErrorCode ierr; 175*7bf3a71bSJakub Kruzik 176*7bf3a71bSJakub Kruzik PetscFunctionBegin; 177*7bf3a71bSJakub Kruzik if (y != z) { 178*7bf3a71bSJakub Kruzik ierr = MatMult(A,x,z);CHKERRQ(ierr); 179*7bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 180*7bf3a71bSJakub Kruzik } else { 181*7bf3a71bSJakub Kruzik if (!shell->leftwork) { 182*7bf3a71bSJakub Kruzik ierr = VecDuplicate(z,&shell->leftwork);CHKERRQ(ierr); 183*7bf3a71bSJakub Kruzik } 184*7bf3a71bSJakub Kruzik ierr = MatMult(A,x,shell->leftwork);CHKERRQ(ierr); 185*7bf3a71bSJakub Kruzik ierr = VecCopy(y,z);CHKERRQ(ierr); 186*7bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,shell->leftwork);CHKERRQ(ierr); 187*7bf3a71bSJakub Kruzik } 188*7bf3a71bSJakub Kruzik PetscFunctionReturn(0); 189*7bf3a71bSJakub Kruzik } 190*7bf3a71bSJakub Kruzik 191*7bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z) 192*7bf3a71bSJakub Kruzik { 193*7bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 194*7bf3a71bSJakub Kruzik PetscErrorCode ierr; 195*7bf3a71bSJakub Kruzik 196*7bf3a71bSJakub Kruzik PetscFunctionBegin; 197*7bf3a71bSJakub Kruzik if (y != z) { 198*7bf3a71bSJakub Kruzik ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 199*7bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 200*7bf3a71bSJakub Kruzik } else { 201*7bf3a71bSJakub Kruzik if (!shell->rightwork) { 202*7bf3a71bSJakub Kruzik ierr = VecDuplicate(z,&shell->rightwork);CHKERRQ(ierr); 203*7bf3a71bSJakub Kruzik } 204*7bf3a71bSJakub Kruzik ierr = MatMultTranspose(A,x,shell->rightwork);CHKERRQ(ierr); 205*7bf3a71bSJakub Kruzik ierr = VecCopy(y,z);CHKERRQ(ierr); 206*7bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,shell->rightwork);CHKERRQ(ierr); 207*7bf3a71bSJakub Kruzik } 208*7bf3a71bSJakub Kruzik PetscFunctionReturn(0); 209*7bf3a71bSJakub Kruzik } 210*7bf3a71bSJakub Kruzik 211793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 212793850ffSBarry Smith { 213793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 214793850ffSBarry Smith Mat_CompositeLink next = shell->head; 215793850ffSBarry Smith PetscErrorCode ierr; 216793850ffSBarry Smith 217793850ffSBarry Smith PetscFunctionBegin; 218e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 219e32f2f54SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 220e4fc5df0SBarry Smith 221793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 222793850ffSBarry Smith if (next->next && !shell->work) { 223793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 224793850ffSBarry Smith } 225793850ffSBarry Smith while ((next = next->next)) { 226793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 227793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 228793850ffSBarry Smith } 229e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 230793850ffSBarry Smith PetscFunctionReturn(0); 231793850ffSBarry Smith } 232793850ffSBarry Smith 233793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 234793850ffSBarry Smith { 235b52f573bSBarry Smith PetscErrorCode ierr; 236ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 237b52f573bSBarry Smith 238793850ffSBarry Smith PetscFunctionBegin; 239c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Y)->options,((PetscObject)Y)->prefix,"-mat_composite_merge",&flg,NULL);CHKERRQ(ierr); 240b52f573bSBarry Smith if (flg) { 241b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 242b52f573bSBarry Smith } 243793850ffSBarry Smith PetscFunctionReturn(0); 244793850ffSBarry Smith } 245793850ffSBarry Smith 246e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 247e4fc5df0SBarry Smith { 248e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 249e4fc5df0SBarry Smith 250e4fc5df0SBarry Smith PetscFunctionBegin; 251321429dbSBarry Smith a->scale *= alpha; 252e4fc5df0SBarry Smith PetscFunctionReturn(0); 253e4fc5df0SBarry Smith } 254e4fc5df0SBarry Smith 255e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 256e4fc5df0SBarry Smith { 257e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 258e4fc5df0SBarry Smith PetscErrorCode ierr; 259e4fc5df0SBarry Smith 260e4fc5df0SBarry Smith PetscFunctionBegin; 261e4fc5df0SBarry Smith if (left) { 262321429dbSBarry Smith if (!a->left) { 263e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 264e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 265321429dbSBarry Smith } else { 266321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 267321429dbSBarry Smith } 268e4fc5df0SBarry Smith } 269e4fc5df0SBarry Smith if (right) { 270321429dbSBarry Smith if (!a->right) { 271e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 272e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 273321429dbSBarry Smith } else { 274321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 275321429dbSBarry Smith } 276e4fc5df0SBarry Smith } 277e4fc5df0SBarry Smith PetscFunctionReturn(0); 278e4fc5df0SBarry Smith } 279793850ffSBarry Smith 2802c0821f3SBarry Smith /*@ 281793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 282793850ffSBarry Smith 283793850ffSBarry Smith Collective on MPI_Comm 284793850ffSBarry Smith 285793850ffSBarry Smith Input Parameters: 286793850ffSBarry Smith + comm - MPI communicator 287793850ffSBarry Smith . nmat - number of matrices to put in 288793850ffSBarry Smith - mats - the matrices 289793850ffSBarry Smith 290793850ffSBarry Smith Output Parameter: 291db36af27SMatthew G. Knepley . mat - the matrix 292793850ffSBarry Smith 293793850ffSBarry Smith Level: advanced 294793850ffSBarry Smith 295793850ffSBarry Smith Notes: 296793850ffSBarry Smith Alternative construction 297793850ffSBarry Smith $ MatCreate(comm,&mat); 298793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 299793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 300793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 301793850ffSBarry Smith $ .... 302793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 303b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 304b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 305793850ffSBarry Smith 3066d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 3076d7c1e57SBarry Smith 3086dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeGetMat(), MatCompositeMerge(), MatCompositeSetType(), MATCOMPOSITE 309793850ffSBarry Smith 310793850ffSBarry Smith @*/ 3117087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 312793850ffSBarry Smith { 313793850ffSBarry Smith PetscErrorCode ierr; 314793850ffSBarry Smith PetscInt m,n,M,N,i; 315793850ffSBarry Smith 316793850ffSBarry Smith PetscFunctionBegin; 317e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 318f3f84630SBarry Smith PetscValidPointer(mat,3); 319793850ffSBarry Smith 320c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[0],PETSC_IGNORE,&n);CHKERRQ(ierr); 321c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE);CHKERRQ(ierr); 322c4f13f43SJakub Kruzik ierr = MatGetSize(mats[0],PETSC_IGNORE,&N);CHKERRQ(ierr); 323c4f13f43SJakub Kruzik ierr = MatGetSize(mats[nmat-1],&M,PETSC_IGNORE);CHKERRQ(ierr); 324793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 325793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 326793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 327793850ffSBarry Smith for (i=0; i<nmat; i++) { 328793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 329793850ffSBarry Smith } 330b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 331b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 332793850ffSBarry Smith PetscFunctionReturn(0); 333793850ffSBarry Smith } 334793850ffSBarry Smith 335d7f81bf2SJakub Kruzik 336d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 337d7f81bf2SJakub Kruzik { 338d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 339d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 340d7f81bf2SJakub Kruzik PetscErrorCode ierr; 341d7f81bf2SJakub Kruzik 342d7f81bf2SJakub Kruzik PetscFunctionBegin; 343d7f81bf2SJakub Kruzik ierr = PetscNewLog(mat,&ilink);CHKERRQ(ierr); 344d7f81bf2SJakub Kruzik ilink->next = 0; 345d7f81bf2SJakub Kruzik ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 346d7f81bf2SJakub Kruzik ilink->mat = smat; 347d7f81bf2SJakub Kruzik 348d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 349d7f81bf2SJakub Kruzik else { 350d7f81bf2SJakub Kruzik while (next->next) { 351d7f81bf2SJakub Kruzik next = next->next; 352d7f81bf2SJakub Kruzik } 353d7f81bf2SJakub Kruzik next->next = ilink; 354d7f81bf2SJakub Kruzik ilink->prev = next; 355d7f81bf2SJakub Kruzik } 356d7f81bf2SJakub Kruzik shell->tail = ilink; 357d7f81bf2SJakub Kruzik shell->nmat += 1; 358d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 359d7f81bf2SJakub Kruzik } 360d7f81bf2SJakub Kruzik 361793850ffSBarry Smith /*@ 362793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 363793850ffSBarry Smith 364793850ffSBarry Smith Collective on Mat 365793850ffSBarry Smith 366793850ffSBarry Smith Input Parameters: 367793850ffSBarry Smith + mat - the composite matrix 368793850ffSBarry Smith - smat - the partial matrix 369793850ffSBarry Smith 370793850ffSBarry Smith Level: advanced 371793850ffSBarry Smith 372793850ffSBarry Smith .seealso: MatCreateComposite() 373793850ffSBarry Smith @*/ 3747087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 375793850ffSBarry Smith { 376793850ffSBarry Smith PetscErrorCode ierr; 377793850ffSBarry Smith 378793850ffSBarry Smith PetscFunctionBegin; 3790700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 3800700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 381d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));CHKERRQ(ierr); 382d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 383d7f81bf2SJakub Kruzik } 384793850ffSBarry Smith 385d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 386d7f81bf2SJakub Kruzik { 387d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 388d7f81bf2SJakub Kruzik 389d7f81bf2SJakub Kruzik PetscFunctionBegin; 390d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 391d7f81bf2SJakub Kruzik mat->ops->getdiagonal = 0; 392d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 393d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 394d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_MULTIPLICATIVE; 395d7f81bf2SJakub Kruzik } else { 396d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 397d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 398d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 399d7f81bf2SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 400793850ffSBarry Smith } 4016d7c1e57SBarry Smith PetscFunctionReturn(0); 4026d7c1e57SBarry Smith } 4036d7c1e57SBarry Smith 4042c0821f3SBarry Smith /*@ 4056d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 4066d7c1e57SBarry Smith 4076d7c1e57SBarry Smith Collective on MPI_Comm 4086d7c1e57SBarry Smith 4096d7c1e57SBarry Smith Input Parameters: 4106d7c1e57SBarry Smith . mat - the composite matrix 4116d7c1e57SBarry Smith 4126d7c1e57SBarry Smith 4136d7c1e57SBarry Smith Level: advanced 4146d7c1e57SBarry Smith 4156d7c1e57SBarry Smith Notes: 4166d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 4176d7c1e57SBarry Smith matrix in the composite matrix. 4186d7c1e57SBarry Smith 4196dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeGetType(), MATCOMPOSITE 4206d7c1e57SBarry Smith 4216d7c1e57SBarry Smith @*/ 4227087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 4236d7c1e57SBarry Smith { 4246d7c1e57SBarry Smith PetscErrorCode ierr; 4256d7c1e57SBarry Smith 4266d7c1e57SBarry Smith PetscFunctionBegin; 427d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 428d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));CHKERRQ(ierr); 429793850ffSBarry Smith PetscFunctionReturn(0); 430793850ffSBarry Smith } 431793850ffSBarry Smith 4326dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type) 4336dbc55e5SJakub Kruzik { 4346dbc55e5SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 4356dbc55e5SJakub Kruzik 4366dbc55e5SJakub Kruzik PetscFunctionBegin; 4376dbc55e5SJakub Kruzik *type = b->type; 4386dbc55e5SJakub Kruzik PetscFunctionReturn(0); 4396dbc55e5SJakub Kruzik } 4406dbc55e5SJakub Kruzik 4416dbc55e5SJakub Kruzik /*@ 4426dbc55e5SJakub Kruzik MatCompositeGetType - Returns type of composite. 4436dbc55e5SJakub Kruzik 4446dbc55e5SJakub Kruzik Not Collective 4456dbc55e5SJakub Kruzik 4466dbc55e5SJakub Kruzik Input Parameter: 4476dbc55e5SJakub Kruzik . mat - the composite matrix 4486dbc55e5SJakub Kruzik 4496dbc55e5SJakub Kruzik Output Parameter: 4506dbc55e5SJakub Kruzik . type - type of composite 4516dbc55e5SJakub Kruzik 4526dbc55e5SJakub Kruzik Level: advanced 4536dbc55e5SJakub Kruzik 4546dbc55e5SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetType(), MATCOMPOSITE 4556dbc55e5SJakub Kruzik 4566dbc55e5SJakub Kruzik @*/ 4576dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type) 4586dbc55e5SJakub Kruzik { 4596dbc55e5SJakub Kruzik PetscErrorCode ierr; 4606dbc55e5SJakub Kruzik 4616dbc55e5SJakub Kruzik PetscFunctionBegin; 4626dbc55e5SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4636dbc55e5SJakub Kruzik PetscValidPointer(type,2); 4646dbc55e5SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type));CHKERRQ(ierr); 4656dbc55e5SJakub Kruzik PetscFunctionReturn(0); 4666dbc55e5SJakub Kruzik } 4676dbc55e5SJakub Kruzik 468d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 469b52f573bSBarry Smith { 470b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 4716d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 472b52f573bSBarry Smith PetscErrorCode ierr; 4736d7c1e57SBarry Smith Mat tmat,newmat; 4741ba60692SJed Brown Vec left,right; 4751ba60692SJed Brown PetscScalar scale; 476b52f573bSBarry Smith 477b52f573bSBarry Smith PetscFunctionBegin; 478e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 479b52f573bSBarry Smith 480b52f573bSBarry Smith PetscFunctionBegin; 4816d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 482b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 483b52f573bSBarry Smith while ((next = next->next)) { 484b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 485b52f573bSBarry Smith } 4866d7c1e57SBarry Smith } else { 48704d534ceSJakub Kruzik if (shell->mergefromright) { 4886d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 489b6cfcaf5SJakub Kruzik while ((next = next->next)) { 490b6cfcaf5SJakub Kruzik ierr = MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 4916bf464f9SBarry Smith ierr = MatDestroy(&tmat);CHKERRQ(ierr); 4926d7c1e57SBarry Smith tmat = newmat; 4936d7c1e57SBarry Smith } 49404d534ceSJakub Kruzik } else { 49504d534ceSJakub Kruzik ierr = MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 49604d534ceSJakub Kruzik while ((prev = prev->prev)) { 49704d534ceSJakub Kruzik ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 49804d534ceSJakub Kruzik ierr = MatDestroy(&tmat);CHKERRQ(ierr); 49904d534ceSJakub Kruzik tmat = newmat; 50004d534ceSJakub Kruzik } 50104d534ceSJakub Kruzik } 5026d7c1e57SBarry Smith } 5031ba60692SJed Brown 5041ba60692SJed Brown scale = shell->scale; 5051ba60692SJed Brown if ((left = shell->left)) {ierr = PetscObjectReference((PetscObject)left);CHKERRQ(ierr);} 5061ba60692SJed Brown if ((right = shell->right)) {ierr = PetscObjectReference((PetscObject)right);CHKERRQ(ierr);} 5071ba60692SJed Brown 50828be2f97SBarry Smith ierr = MatHeaderReplace(mat,&tmat);CHKERRQ(ierr); 5091ba60692SJed Brown 5101ba60692SJed Brown ierr = MatDiagonalScale(mat,left,right);CHKERRQ(ierr); 5111ba60692SJed Brown ierr = MatScale(mat,scale);CHKERRQ(ierr); 5121ba60692SJed Brown ierr = VecDestroy(&left);CHKERRQ(ierr); 5131ba60692SJed Brown ierr = VecDestroy(&right);CHKERRQ(ierr); 514b52f573bSBarry Smith PetscFunctionReturn(0); 515b52f573bSBarry Smith } 5166a7ede75SJakub Kruzik 5176a7ede75SJakub Kruzik /*@ 518d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 519d7f81bf2SJakub Kruzik by summing all the matrices inside the composite matrix. 520d7f81bf2SJakub Kruzik 521d7f81bf2SJakub Kruzik Collective on MPI_Comm 522d7f81bf2SJakub Kruzik 523d7f81bf2SJakub Kruzik Input Parameters: 524d7f81bf2SJakub Kruzik . mat - the composite matrix 525d7f81bf2SJakub Kruzik 526d7f81bf2SJakub Kruzik 527d7f81bf2SJakub Kruzik Options Database: 528d7f81bf2SJakub Kruzik . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 529d7f81bf2SJakub Kruzik 530d7f81bf2SJakub Kruzik Level: advanced 531d7f81bf2SJakub Kruzik 532d7f81bf2SJakub Kruzik Notes: 533d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 534d7f81bf2SJakub Kruzik matrix in the composite matrix. 535d7f81bf2SJakub Kruzik 53604d534ceSJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeSetMergeFromRight(), MATCOMPOSITE 537d7f81bf2SJakub Kruzik 538d7f81bf2SJakub Kruzik @*/ 539d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 540d7f81bf2SJakub Kruzik { 541d7f81bf2SJakub Kruzik PetscErrorCode ierr; 542d7f81bf2SJakub Kruzik 543d7f81bf2SJakub Kruzik PetscFunctionBegin; 544d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 545d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));CHKERRQ(ierr); 546d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 547d7f81bf2SJakub Kruzik } 548d7f81bf2SJakub Kruzik 549d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetNMat_Composite(Mat mat,PetscInt *nmat) 550d7f81bf2SJakub Kruzik { 551d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 552d7f81bf2SJakub Kruzik 553d7f81bf2SJakub Kruzik PetscFunctionBegin; 554d7f81bf2SJakub Kruzik *nmat = shell->nmat; 555d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 556d7f81bf2SJakub Kruzik } 557d7f81bf2SJakub Kruzik 558d7f81bf2SJakub Kruzik /*@ 5596a7ede75SJakub Kruzik MatCompositeGetNMat - Returns the number of matrices in composite. 5606a7ede75SJakub Kruzik 5616a7ede75SJakub Kruzik Not Collective 5626a7ede75SJakub Kruzik 5636a7ede75SJakub Kruzik Input Parameter: 564d7f81bf2SJakub Kruzik . mat - the composite matrix 5656a7ede75SJakub Kruzik 5666a7ede75SJakub Kruzik Output Parameter: 5676dbc55e5SJakub Kruzik + size - the local size 5686dbc55e5SJakub Kruzik - nmat - number of matrices in composite 5696a7ede75SJakub Kruzik 5708b5c3584SJakub Kruzik Level: advanced 5716a7ede75SJakub Kruzik 5728b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat() 5736a7ede75SJakub Kruzik 5746a7ede75SJakub Kruzik @*/ 575d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetNMat(Mat mat,PetscInt *nmat) 5766a7ede75SJakub Kruzik { 577d7f81bf2SJakub Kruzik PetscErrorCode ierr; 5786a7ede75SJakub Kruzik 5796a7ede75SJakub Kruzik PetscFunctionBegin; 580d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5816a7ede75SJakub Kruzik PetscValidPointer(nmat,2); 582d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetNMat_C",(Mat,PetscInt*),(mat,nmat));CHKERRQ(ierr); 583d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 584d7f81bf2SJakub Kruzik } 585d7f81bf2SJakub Kruzik 586d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 587d7f81bf2SJakub Kruzik { 588d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 589d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 590d7f81bf2SJakub Kruzik PetscInt k; 591d7f81bf2SJakub Kruzik 592d7f81bf2SJakub Kruzik PetscFunctionBegin; 593d7f81bf2SJakub Kruzik if (i >= shell->nmat) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %d >= %d",i,shell->nmat); 594d7f81bf2SJakub Kruzik ilink = shell->head; 595d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 596d7f81bf2SJakub Kruzik ilink = ilink->next; 597d7f81bf2SJakub Kruzik } 598d7f81bf2SJakub Kruzik *Ai = ilink->mat; 5996a7ede75SJakub Kruzik PetscFunctionReturn(0); 6006a7ede75SJakub Kruzik } 6016a7ede75SJakub Kruzik 6028b5c3584SJakub Kruzik /*@ 6038b5c3584SJakub Kruzik MatCompositeGetMat - Returns the ith matrix from composite. 6048b5c3584SJakub Kruzik 6058b5c3584SJakub Kruzik Logically Collective on Mat 6068b5c3584SJakub Kruzik 6078b5c3584SJakub Kruzik Input Parameter: 608d7f81bf2SJakub Kruzik + mat - the composite matrix 6098b5c3584SJakub Kruzik - i - the number of requested matrix 6108b5c3584SJakub Kruzik 6118b5c3584SJakub Kruzik Output Parameter: 6128b5c3584SJakub Kruzik . Ai - ith matrix in composite 6138b5c3584SJakub Kruzik 6148b5c3584SJakub Kruzik Level: advanced 6158b5c3584SJakub Kruzik 6168b5c3584SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNMat() 6178b5c3584SJakub Kruzik 6188b5c3584SJakub Kruzik @*/ 619d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 6208b5c3584SJakub Kruzik { 621d7f81bf2SJakub Kruzik PetscErrorCode ierr; 6228b5c3584SJakub Kruzik 6238b5c3584SJakub Kruzik PetscFunctionBegin; 624d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 625d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 6268b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 627d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));CHKERRQ(ierr); 6288b5c3584SJakub Kruzik PetscFunctionReturn(0); 6298b5c3584SJakub Kruzik } 6308b5c3584SJakub Kruzik 63104d534ceSJakub Kruzik static PetscErrorCode MatCompositeSetMergeFromRight_Composite(Mat mat,PetscBool flg) 63204d534ceSJakub Kruzik { 63304d534ceSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 63404d534ceSJakub Kruzik 63504d534ceSJakub Kruzik PetscFunctionBegin; 63604d534ceSJakub Kruzik shell->mergefromright = flg; 63704d534ceSJakub Kruzik PetscFunctionReturn(0); 63804d534ceSJakub Kruzik } 63904d534ceSJakub Kruzik 64004d534ceSJakub Kruzik /*@ 64104d534ceSJakub Kruzik MatCompositeSetMergeFromRight - Sets MatCompositeMerge() to start from right. 64204d534ceSJakub Kruzik 64304d534ceSJakub Kruzik Logically Collective on Mat 64404d534ceSJakub Kruzik 64504d534ceSJakub Kruzik Input Parameter: 64604d534ceSJakub Kruzik + mat - the composite matrix 64704d534ceSJakub Kruzik - flg - if true (default) the merge starts with the first added matrix (mat[0]) 64804d534ceSJakub Kruzik 64904d534ceSJakub Kruzik Level: advanced 65004d534ceSJakub Kruzik 65104d534ceSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeMerge() 65204d534ceSJakub Kruzik 65304d534ceSJakub Kruzik @*/ 65404d534ceSJakub Kruzik PetscErrorCode MatCompositeSetMergeFromRight(Mat mat,PetscBool flg) 65504d534ceSJakub Kruzik { 65604d534ceSJakub Kruzik PetscErrorCode ierr; 65704d534ceSJakub Kruzik 65804d534ceSJakub Kruzik PetscFunctionBegin; 65904d534ceSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 66004d534ceSJakub Kruzik PetscValidLogicalCollectiveBool(mat,flg,2); 66104d534ceSJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetMergeFromRight_C",(Mat,PetscBool),(mat,flg));CHKERRQ(ierr); 66204d534ceSJakub Kruzik PetscFunctionReturn(0); 66304d534ceSJakub Kruzik } 66404d534ceSJakub Kruzik 66541cd0125SJakub Kruzik static struct _MatOps MatOps_Values = {0, 66641cd0125SJakub Kruzik 0, 66741cd0125SJakub Kruzik 0, 66841cd0125SJakub Kruzik MatMult_Composite, 669*7bf3a71bSJakub Kruzik MatMultAdd_Composite, 67041cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 671*7bf3a71bSJakub Kruzik MatMultTransposeAdd_Composite, 67241cd0125SJakub Kruzik 0, 67341cd0125SJakub Kruzik 0, 67441cd0125SJakub Kruzik 0, 67541cd0125SJakub Kruzik /* 10*/ 0, 67641cd0125SJakub Kruzik 0, 67741cd0125SJakub Kruzik 0, 67841cd0125SJakub Kruzik 0, 67941cd0125SJakub Kruzik 0, 68041cd0125SJakub Kruzik /* 15*/ 0, 68141cd0125SJakub Kruzik 0, 68241cd0125SJakub Kruzik MatGetDiagonal_Composite, 68341cd0125SJakub Kruzik MatDiagonalScale_Composite, 68441cd0125SJakub Kruzik 0, 68541cd0125SJakub Kruzik /* 20*/ 0, 68641cd0125SJakub Kruzik MatAssemblyEnd_Composite, 68741cd0125SJakub Kruzik 0, 68841cd0125SJakub Kruzik 0, 68941cd0125SJakub Kruzik /* 24*/ 0, 69041cd0125SJakub Kruzik 0, 69141cd0125SJakub Kruzik 0, 69241cd0125SJakub Kruzik 0, 69341cd0125SJakub Kruzik 0, 69441cd0125SJakub Kruzik /* 29*/ 0, 69541cd0125SJakub Kruzik 0, 69641cd0125SJakub Kruzik 0, 69741cd0125SJakub Kruzik 0, 69841cd0125SJakub Kruzik 0, 69941cd0125SJakub Kruzik /* 34*/ 0, 70041cd0125SJakub Kruzik 0, 70141cd0125SJakub Kruzik 0, 70241cd0125SJakub Kruzik 0, 70341cd0125SJakub Kruzik 0, 70441cd0125SJakub Kruzik /* 39*/ 0, 70541cd0125SJakub Kruzik 0, 70641cd0125SJakub Kruzik 0, 70741cd0125SJakub Kruzik 0, 70841cd0125SJakub Kruzik 0, 70941cd0125SJakub Kruzik /* 44*/ 0, 71041cd0125SJakub Kruzik MatScale_Composite, 71141cd0125SJakub Kruzik MatShift_Basic, 71241cd0125SJakub Kruzik 0, 71341cd0125SJakub Kruzik 0, 71441cd0125SJakub Kruzik /* 49*/ 0, 71541cd0125SJakub Kruzik 0, 71641cd0125SJakub Kruzik 0, 71741cd0125SJakub Kruzik 0, 71841cd0125SJakub Kruzik 0, 71941cd0125SJakub Kruzik /* 54*/ 0, 72041cd0125SJakub Kruzik 0, 72141cd0125SJakub Kruzik 0, 72241cd0125SJakub Kruzik 0, 72341cd0125SJakub Kruzik 0, 72441cd0125SJakub Kruzik /* 59*/ 0, 72541cd0125SJakub Kruzik MatDestroy_Composite, 72641cd0125SJakub Kruzik 0, 72741cd0125SJakub Kruzik 0, 72841cd0125SJakub Kruzik 0, 72941cd0125SJakub Kruzik /* 64*/ 0, 73041cd0125SJakub Kruzik 0, 73141cd0125SJakub Kruzik 0, 73241cd0125SJakub Kruzik 0, 73341cd0125SJakub Kruzik 0, 73441cd0125SJakub Kruzik /* 69*/ 0, 73541cd0125SJakub Kruzik 0, 73641cd0125SJakub Kruzik 0, 73741cd0125SJakub Kruzik 0, 73841cd0125SJakub Kruzik 0, 73941cd0125SJakub Kruzik /* 74*/ 0, 74041cd0125SJakub Kruzik 0, 74141cd0125SJakub Kruzik 0, 74241cd0125SJakub Kruzik 0, 74341cd0125SJakub Kruzik 0, 74441cd0125SJakub Kruzik /* 79*/ 0, 74541cd0125SJakub Kruzik 0, 74641cd0125SJakub Kruzik 0, 74741cd0125SJakub Kruzik 0, 74841cd0125SJakub Kruzik 0, 74941cd0125SJakub Kruzik /* 84*/ 0, 75041cd0125SJakub Kruzik 0, 75141cd0125SJakub Kruzik 0, 75241cd0125SJakub Kruzik 0, 75341cd0125SJakub Kruzik 0, 75441cd0125SJakub Kruzik /* 89*/ 0, 75541cd0125SJakub Kruzik 0, 75641cd0125SJakub Kruzik 0, 75741cd0125SJakub Kruzik 0, 75841cd0125SJakub Kruzik 0, 75941cd0125SJakub Kruzik /* 94*/ 0, 76041cd0125SJakub Kruzik 0, 76141cd0125SJakub Kruzik 0, 76241cd0125SJakub Kruzik 0, 76341cd0125SJakub Kruzik 0, 76441cd0125SJakub Kruzik /*99*/ 0, 76541cd0125SJakub Kruzik 0, 76641cd0125SJakub Kruzik 0, 76741cd0125SJakub Kruzik 0, 76841cd0125SJakub Kruzik 0, 76941cd0125SJakub Kruzik /*104*/ 0, 77041cd0125SJakub Kruzik 0, 77141cd0125SJakub Kruzik 0, 77241cd0125SJakub Kruzik 0, 77341cd0125SJakub Kruzik 0, 77441cd0125SJakub Kruzik /*109*/ 0, 77541cd0125SJakub Kruzik 0, 77641cd0125SJakub Kruzik 0, 77741cd0125SJakub Kruzik 0, 77841cd0125SJakub Kruzik 0, 77941cd0125SJakub Kruzik /*114*/ 0, 78041cd0125SJakub Kruzik 0, 78141cd0125SJakub Kruzik 0, 78241cd0125SJakub Kruzik 0, 78341cd0125SJakub Kruzik 0, 78441cd0125SJakub Kruzik /*119*/ 0, 78541cd0125SJakub Kruzik 0, 78641cd0125SJakub Kruzik 0, 78741cd0125SJakub Kruzik 0, 78841cd0125SJakub Kruzik 0, 78941cd0125SJakub Kruzik /*124*/ 0, 79041cd0125SJakub Kruzik 0, 79141cd0125SJakub Kruzik 0, 79241cd0125SJakub Kruzik 0, 79341cd0125SJakub Kruzik 0, 79441cd0125SJakub Kruzik /*129*/ 0, 79541cd0125SJakub Kruzik 0, 79641cd0125SJakub Kruzik 0, 79741cd0125SJakub Kruzik 0, 79841cd0125SJakub Kruzik 0, 79941cd0125SJakub Kruzik /*134*/ 0, 80041cd0125SJakub Kruzik 0, 80141cd0125SJakub Kruzik 0, 80241cd0125SJakub Kruzik 0, 80341cd0125SJakub Kruzik 0, 80441cd0125SJakub Kruzik /*139*/ 0, 80541cd0125SJakub Kruzik 0, 80641cd0125SJakub Kruzik 0 80741cd0125SJakub Kruzik }; 80841cd0125SJakub Kruzik 80941cd0125SJakub Kruzik /*MC 81041cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 81141cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 81241cd0125SJakub Kruzik 81341cd0125SJakub Kruzik Notes: 81441cd0125SJakub Kruzik to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 81541cd0125SJakub Kruzik 81641cd0125SJakub Kruzik Level: advanced 81741cd0125SJakub Kruzik 8186dbc55e5SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeSetType(), MatCompositeGetType(), MatCompositeMerge(), MatCompositeSetMergeFromRight(), MatCompositeGetNmat(), MatCompositeGetMat() 81941cd0125SJakub Kruzik M*/ 82041cd0125SJakub Kruzik 82141cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 82241cd0125SJakub Kruzik { 82341cd0125SJakub Kruzik Mat_Composite *b; 82441cd0125SJakub Kruzik PetscErrorCode ierr; 82541cd0125SJakub Kruzik 82641cd0125SJakub Kruzik PetscFunctionBegin; 82741cd0125SJakub Kruzik ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 82841cd0125SJakub Kruzik A->data = (void*)b; 82941cd0125SJakub Kruzik ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 83041cd0125SJakub Kruzik 83141cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 83241cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 83341cd0125SJakub Kruzik 83441cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 83541cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 83641cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 83741cd0125SJakub Kruzik b->scale = 1.0; 83841cd0125SJakub Kruzik b->nmat = 0; 83904d534ceSJakub Kruzik b->mergefromright = PETSC_TRUE; 84041cd0125SJakub Kruzik 84141cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite);CHKERRQ(ierr); 84241cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite);CHKERRQ(ierr); 8436dbc55e5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite);CHKERRQ(ierr); 84441cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite);CHKERRQ(ierr); 84541cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNMat_C",MatCompositeGetNMat_Composite);CHKERRQ(ierr); 84641cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite);CHKERRQ(ierr); 84704d534ceSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeFromRight_C",MatCompositeSetMergeFromRight_Composite);CHKERRQ(ierr); 84841cd0125SJakub Kruzik 84941cd0125SJakub Kruzik ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 85041cd0125SJakub Kruzik PetscFunctionReturn(0); 85141cd0125SJakub Kruzik } 85241cd0125SJakub Kruzik 853