1793850ffSBarry Smith 2b45d2f2cSJed Brown #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; 18793850ffSBarry Smith } Mat_Composite; 19793850ffSBarry Smith 20793850ffSBarry Smith #undef __FUNCT__ 21793850ffSBarry Smith #define __FUNCT__ "MatDestroy_Composite" 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 47793850ffSBarry Smith #undef __FUNCT__ 486d7c1e57SBarry Smith #define __FUNCT__ "MatMult_Composite_Multiplicative" 496d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 506d7c1e57SBarry Smith { 516d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 526d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 536d7c1e57SBarry Smith PetscErrorCode ierr; 546d7c1e57SBarry Smith Vec in,out; 556d7c1e57SBarry Smith 566d7c1e57SBarry Smith PetscFunctionBegin; 57e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 586d7c1e57SBarry Smith in = x; 59e4fc5df0SBarry Smith if (shell->right) { 60e4fc5df0SBarry Smith if (!shell->rightwork) { 61e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 62e4fc5df0SBarry Smith } 63e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 64e4fc5df0SBarry Smith in = shell->rightwork; 65e4fc5df0SBarry Smith } 666d7c1e57SBarry Smith while (next->next) { 676d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 686d7c1e57SBarry Smith ierr = MatGetVecs(next->mat,PETSC_NULL,&next->work);CHKERRQ(ierr); 696d7c1e57SBarry Smith } 706d7c1e57SBarry Smith out = next->work; 716d7c1e57SBarry Smith ierr = MatMult(next->mat,in,out);CHKERRQ(ierr); 726d7c1e57SBarry Smith in = out; 736d7c1e57SBarry Smith next = next->next; 746d7c1e57SBarry Smith } 756d7c1e57SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 76e4fc5df0SBarry Smith if (shell->left) { 77e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 78e4fc5df0SBarry Smith } 79e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 806d7c1e57SBarry Smith PetscFunctionReturn(0); 816d7c1e57SBarry Smith } 826d7c1e57SBarry Smith 836d7c1e57SBarry Smith #undef __FUNCT__ 846d7c1e57SBarry Smith #define __FUNCT__ "MatMultTranspose_Composite_Multiplicative" 856d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 866d7c1e57SBarry Smith { 876d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 886d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 896d7c1e57SBarry Smith PetscErrorCode ierr; 906d7c1e57SBarry Smith Vec in,out; 916d7c1e57SBarry Smith 926d7c1e57SBarry Smith PetscFunctionBegin; 93e32f2f54SBarry Smith if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 946d7c1e57SBarry Smith in = x; 95e4fc5df0SBarry Smith if (shell->left) { 96e4fc5df0SBarry Smith if (!shell->leftwork) { 97e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 98e4fc5df0SBarry Smith } 99e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 100e4fc5df0SBarry Smith in = shell->leftwork; 101e4fc5df0SBarry Smith } 1026d7c1e57SBarry Smith while (tail->prev) { 1036d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1046d7c1e57SBarry Smith ierr = MatGetVecs(tail->mat,PETSC_NULL,&tail->prev->work);CHKERRQ(ierr); 1056d7c1e57SBarry Smith } 1066d7c1e57SBarry Smith out = tail->prev->work; 1076d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1086d7c1e57SBarry Smith in = out; 1096d7c1e57SBarry Smith tail = tail->prev; 1106d7c1e57SBarry Smith } 1116d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 112e4fc5df0SBarry Smith if (shell->right) { 113e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 114e4fc5df0SBarry Smith } 115e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1166d7c1e57SBarry Smith PetscFunctionReturn(0); 1176d7c1e57SBarry Smith } 1186d7c1e57SBarry Smith 1196d7c1e57SBarry Smith #undef __FUNCT__ 120793850ffSBarry Smith #define __FUNCT__ "MatMult_Composite" 121793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 122793850ffSBarry Smith { 123793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 124793850ffSBarry Smith Mat_CompositeLink next = shell->head; 125793850ffSBarry Smith PetscErrorCode ierr; 126e4fc5df0SBarry Smith Vec in; 127793850ffSBarry Smith 128793850ffSBarry Smith PetscFunctionBegin; 129e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 130e4fc5df0SBarry Smith in = x; 131e4fc5df0SBarry Smith if (shell->right) { 132e4fc5df0SBarry Smith if (!shell->rightwork) { 133e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 134793850ffSBarry Smith } 135e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 136e4fc5df0SBarry Smith in = shell->rightwork; 137e4fc5df0SBarry Smith } 138e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 139e4fc5df0SBarry Smith while ((next = next->next)) { 140e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 141e4fc5df0SBarry Smith } 142e4fc5df0SBarry Smith if (shell->left) { 143e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 144e4fc5df0SBarry Smith } 145e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 146793850ffSBarry Smith PetscFunctionReturn(0); 147793850ffSBarry Smith } 148793850ffSBarry Smith 149793850ffSBarry Smith #undef __FUNCT__ 150793850ffSBarry Smith #define __FUNCT__ "MatMultTranspose_Composite" 151793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 152793850ffSBarry Smith { 153793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 154793850ffSBarry Smith Mat_CompositeLink next = shell->head; 155793850ffSBarry Smith PetscErrorCode ierr; 156e4fc5df0SBarry Smith Vec in; 157793850ffSBarry Smith 158793850ffSBarry Smith PetscFunctionBegin; 159e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 160e4fc5df0SBarry Smith in = x; 161e4fc5df0SBarry Smith if (shell->left) { 162e4fc5df0SBarry Smith if (!shell->leftwork) { 163e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 164793850ffSBarry Smith } 165e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 166e4fc5df0SBarry Smith in = shell->leftwork; 167e4fc5df0SBarry Smith } 168e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 169e4fc5df0SBarry Smith while ((next = next->next)) { 170e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 171e4fc5df0SBarry Smith } 172e4fc5df0SBarry Smith if (shell->right) { 173e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 174e4fc5df0SBarry Smith } 175e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 176793850ffSBarry Smith PetscFunctionReturn(0); 177793850ffSBarry Smith } 178793850ffSBarry Smith 179793850ffSBarry Smith #undef __FUNCT__ 180793850ffSBarry Smith #define __FUNCT__ "MatGetDiagonal_Composite" 181793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 182793850ffSBarry Smith { 183793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 184793850ffSBarry Smith Mat_CompositeLink next = shell->head; 185793850ffSBarry Smith PetscErrorCode ierr; 186793850ffSBarry Smith 187793850ffSBarry Smith PetscFunctionBegin; 188e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 189e32f2f54SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 190e4fc5df0SBarry Smith 191793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 192793850ffSBarry Smith if (next->next && !shell->work) { 193793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 194793850ffSBarry Smith } 195793850ffSBarry Smith while ((next = next->next)) { 196793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 197793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 198793850ffSBarry Smith } 199e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 200793850ffSBarry Smith PetscFunctionReturn(0); 201793850ffSBarry Smith } 202793850ffSBarry Smith 203793850ffSBarry Smith #undef __FUNCT__ 204793850ffSBarry Smith #define __FUNCT__ "MatAssemblyEnd_Composite" 205793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 206793850ffSBarry Smith { 207b52f573bSBarry Smith PetscErrorCode ierr; 208ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 209b52f573bSBarry Smith 210793850ffSBarry Smith PetscFunctionBegin; 211acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)Y)->prefix,"-mat_composite_merge",&flg,PETSC_NULL);CHKERRQ(ierr); 212b52f573bSBarry Smith if (flg) { 213b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 214b52f573bSBarry Smith } 215793850ffSBarry Smith PetscFunctionReturn(0); 216793850ffSBarry Smith } 217793850ffSBarry Smith 218e4fc5df0SBarry Smith #undef __FUNCT__ 219e4fc5df0SBarry Smith #define __FUNCT__ "MatScale_Composite" 220e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 221e4fc5df0SBarry Smith { 222e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 223e4fc5df0SBarry Smith 224e4fc5df0SBarry Smith PetscFunctionBegin; 225321429dbSBarry Smith a->scale *= alpha; 226e4fc5df0SBarry Smith PetscFunctionReturn(0); 227e4fc5df0SBarry Smith } 228e4fc5df0SBarry Smith 229e4fc5df0SBarry Smith #undef __FUNCT__ 230e4fc5df0SBarry Smith #define __FUNCT__ "MatDiagonalScale_Composite" 231e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 232e4fc5df0SBarry Smith { 233e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 234e4fc5df0SBarry Smith PetscErrorCode ierr; 235e4fc5df0SBarry Smith 236e4fc5df0SBarry Smith PetscFunctionBegin; 237e4fc5df0SBarry Smith if (left) { 238321429dbSBarry Smith if (!a->left) { 239e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 240e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 241321429dbSBarry Smith } else { 242321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 243321429dbSBarry Smith } 244e4fc5df0SBarry Smith } 245e4fc5df0SBarry Smith if (right) { 246321429dbSBarry Smith if (!a->right) { 247e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 248e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 249321429dbSBarry Smith } else { 250321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 251321429dbSBarry Smith } 252e4fc5df0SBarry Smith } 253e4fc5df0SBarry Smith PetscFunctionReturn(0); 254e4fc5df0SBarry Smith } 255793850ffSBarry Smith 256793850ffSBarry Smith static struct _MatOps MatOps_Values = {0, 257793850ffSBarry Smith 0, 258793850ffSBarry Smith 0, 259793850ffSBarry Smith MatMult_Composite, 260793850ffSBarry Smith 0, 261793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite, 262793850ffSBarry Smith 0, 263793850ffSBarry Smith 0, 264793850ffSBarry Smith 0, 265793850ffSBarry Smith 0, 266793850ffSBarry Smith /* 10*/ 0, 267793850ffSBarry Smith 0, 268793850ffSBarry Smith 0, 269793850ffSBarry Smith 0, 270793850ffSBarry Smith 0, 271793850ffSBarry Smith /* 15*/ 0, 272793850ffSBarry Smith 0, 273793850ffSBarry Smith MatGetDiagonal_Composite, 274e4fc5df0SBarry Smith MatDiagonalScale_Composite, 275793850ffSBarry Smith 0, 276793850ffSBarry Smith /* 20*/ 0, 277793850ffSBarry Smith MatAssemblyEnd_Composite, 278793850ffSBarry Smith 0, 279793850ffSBarry Smith 0, 280d519adbfSMatthew Knepley /* 24*/ 0, 281793850ffSBarry Smith 0, 282793850ffSBarry Smith 0, 283793850ffSBarry Smith 0, 284793850ffSBarry Smith 0, 285d519adbfSMatthew Knepley /* 29*/ 0, 286793850ffSBarry Smith 0, 287793850ffSBarry Smith 0, 288793850ffSBarry Smith 0, 289793850ffSBarry Smith 0, 290d519adbfSMatthew Knepley /* 34*/ 0, 291793850ffSBarry Smith 0, 292793850ffSBarry Smith 0, 293793850ffSBarry Smith 0, 294793850ffSBarry Smith 0, 295d519adbfSMatthew Knepley /* 39*/ 0, 296793850ffSBarry Smith 0, 297793850ffSBarry Smith 0, 298793850ffSBarry Smith 0, 299793850ffSBarry Smith 0, 300d519adbfSMatthew Knepley /* 44*/ 0, 301e4fc5df0SBarry Smith MatScale_Composite, 302793850ffSBarry Smith 0, 303793850ffSBarry Smith 0, 304793850ffSBarry Smith 0, 305d519adbfSMatthew Knepley /* 49*/ 0, 306793850ffSBarry Smith 0, 307793850ffSBarry Smith 0, 308793850ffSBarry Smith 0, 309793850ffSBarry Smith 0, 310d519adbfSMatthew Knepley /* 54*/ 0, 311793850ffSBarry Smith 0, 312793850ffSBarry Smith 0, 313793850ffSBarry Smith 0, 314793850ffSBarry Smith 0, 315d519adbfSMatthew Knepley /* 59*/ 0, 316793850ffSBarry Smith MatDestroy_Composite, 317793850ffSBarry Smith 0, 318793850ffSBarry Smith 0, 319793850ffSBarry Smith 0, 320d519adbfSMatthew Knepley /* 64*/ 0, 321793850ffSBarry Smith 0, 322793850ffSBarry Smith 0, 323793850ffSBarry Smith 0, 324793850ffSBarry Smith 0, 325d519adbfSMatthew Knepley /* 69*/ 0, 326793850ffSBarry Smith 0, 327793850ffSBarry Smith 0, 328793850ffSBarry Smith 0, 329793850ffSBarry Smith 0, 330d519adbfSMatthew Knepley /* 74*/ 0, 331793850ffSBarry Smith 0, 332793850ffSBarry Smith 0, 333793850ffSBarry Smith 0, 334793850ffSBarry Smith 0, 335d519adbfSMatthew Knepley /* 79*/ 0, 336793850ffSBarry Smith 0, 337793850ffSBarry Smith 0, 338793850ffSBarry Smith 0, 339793850ffSBarry Smith 0, 340d519adbfSMatthew Knepley /* 84*/ 0, 341793850ffSBarry Smith 0, 342793850ffSBarry Smith 0, 343793850ffSBarry Smith 0, 344793850ffSBarry Smith 0, 345d519adbfSMatthew Knepley /* 89*/ 0, 346793850ffSBarry Smith 0, 347793850ffSBarry Smith 0, 348793850ffSBarry Smith 0, 349793850ffSBarry Smith 0, 350d519adbfSMatthew Knepley /* 94*/ 0, 351793850ffSBarry Smith 0, 352793850ffSBarry Smith 0, 353793850ffSBarry Smith 0}; 354793850ffSBarry Smith 355793850ffSBarry Smith /*MC 3566d7c1e57SBarry Smith MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices (all matrices are of same size and parallel layout). 3576d7c1e57SBarry Smith 3586d7c1e57SBarry Smith Notes: to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 359793850ffSBarry Smith 360793850ffSBarry Smith Level: advanced 361793850ffSBarry Smith 3626d7c1e57SBarry Smith .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 363793850ffSBarry Smith M*/ 364793850ffSBarry Smith 365793850ffSBarry Smith EXTERN_C_BEGIN 366793850ffSBarry Smith #undef __FUNCT__ 367793850ffSBarry Smith #define __FUNCT__ "MatCreate_Composite" 3687087cfbeSBarry Smith PetscErrorCode MatCreate_Composite(Mat A) 369793850ffSBarry Smith { 370793850ffSBarry Smith Mat_Composite *b; 371793850ffSBarry Smith PetscErrorCode ierr; 372793850ffSBarry Smith 373793850ffSBarry Smith PetscFunctionBegin; 37438f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Composite,&b);CHKERRQ(ierr); 375793850ffSBarry Smith A->data = (void*)b; 37638f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 377793850ffSBarry Smith 37826283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 37926283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 380793850ffSBarry Smith 381793850ffSBarry Smith A->assembled = PETSC_TRUE; 3823db03f37SJed Brown A->preallocated = PETSC_TRUE; 3836d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 384e4fc5df0SBarry Smith b->scale = 1.0; 385793850ffSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 386793850ffSBarry Smith PetscFunctionReturn(0); 387793850ffSBarry Smith } 388793850ffSBarry Smith EXTERN_C_END 389793850ffSBarry Smith 390793850ffSBarry Smith #undef __FUNCT__ 391793850ffSBarry Smith #define __FUNCT__ "MatCreateComposite" 392793850ffSBarry Smith /*@C 393793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 394793850ffSBarry Smith 395793850ffSBarry Smith Collective on MPI_Comm 396793850ffSBarry Smith 397793850ffSBarry Smith Input Parameters: 398793850ffSBarry Smith + comm - MPI communicator 399793850ffSBarry Smith . nmat - number of matrices to put in 400793850ffSBarry Smith - mats - the matrices 401793850ffSBarry Smith 402793850ffSBarry Smith Output Parameter: 403793850ffSBarry Smith . A - the matrix 404793850ffSBarry Smith 405793850ffSBarry Smith Level: advanced 406793850ffSBarry Smith 407793850ffSBarry Smith Notes: 408793850ffSBarry Smith Alternative construction 409793850ffSBarry Smith $ MatCreate(comm,&mat); 410793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 411793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 412793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 413793850ffSBarry Smith $ .... 414793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 415b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 416b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 417793850ffSBarry Smith 4186d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4196d7c1e57SBarry Smith 4206d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 421793850ffSBarry Smith 422793850ffSBarry Smith @*/ 4237087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 424793850ffSBarry Smith { 425793850ffSBarry Smith PetscErrorCode ierr; 426793850ffSBarry Smith PetscInt m,n,M,N,i; 427793850ffSBarry Smith 428793850ffSBarry Smith PetscFunctionBegin; 429e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 430f3f84630SBarry Smith PetscValidPointer(mat,3); 431793850ffSBarry Smith 432793850ffSBarry Smith ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr); 433793850ffSBarry Smith ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr); 434793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 435793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 436793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 437793850ffSBarry Smith for (i=0; i<nmat; i++) { 438793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 439793850ffSBarry Smith } 440b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 441b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 442793850ffSBarry Smith PetscFunctionReturn(0); 443793850ffSBarry Smith } 444793850ffSBarry Smith 445793850ffSBarry Smith #undef __FUNCT__ 446793850ffSBarry Smith #define __FUNCT__ "MatCompositeAddMat" 447793850ffSBarry Smith /*@ 448793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 449793850ffSBarry Smith 450793850ffSBarry Smith Collective on Mat 451793850ffSBarry Smith 452793850ffSBarry Smith Input Parameters: 453793850ffSBarry Smith + mat - the composite matrix 454793850ffSBarry Smith - smat - the partial matrix 455793850ffSBarry Smith 456793850ffSBarry Smith Level: advanced 457793850ffSBarry Smith 458793850ffSBarry Smith .seealso: MatCreateComposite() 459793850ffSBarry Smith @*/ 4607087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 461793850ffSBarry Smith { 46238f2d2fdSLisandro Dalcin Mat_Composite *shell; 463793850ffSBarry Smith PetscErrorCode ierr; 464793850ffSBarry Smith Mat_CompositeLink ilink,next; 465793850ffSBarry Smith 466793850ffSBarry Smith PetscFunctionBegin; 4670700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4680700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 46938f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr); 470793850ffSBarry Smith ilink->next = 0; 471793850ffSBarry Smith ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 472c3122656SLisandro Dalcin ilink->mat = smat; 473793850ffSBarry Smith 47438f2d2fdSLisandro Dalcin shell = (Mat_Composite*)mat->data; 475793850ffSBarry Smith next = shell->head; 476*2205254eSKarl Rupp if (!next) shell->head = ilink; 477*2205254eSKarl Rupp else { 478793850ffSBarry Smith while (next->next) { 479793850ffSBarry Smith next = next->next; 480793850ffSBarry Smith } 481793850ffSBarry Smith next->next = ilink; 4826d7c1e57SBarry Smith ilink->prev = next; 4836d7c1e57SBarry Smith } 4846d7c1e57SBarry Smith shell->tail = ilink; 4856d7c1e57SBarry Smith PetscFunctionReturn(0); 4866d7c1e57SBarry Smith } 4876d7c1e57SBarry Smith 4886d7c1e57SBarry Smith #undef __FUNCT__ 4896d7c1e57SBarry Smith #define __FUNCT__ "MatCompositeSetType" 4906d7c1e57SBarry Smith /*@C 4916d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 4926d7c1e57SBarry Smith 4936d7c1e57SBarry Smith Collective on MPI_Comm 4946d7c1e57SBarry Smith 4956d7c1e57SBarry Smith Input Parameters: 4966d7c1e57SBarry Smith . mat - the composite matrix 4976d7c1e57SBarry Smith 4986d7c1e57SBarry Smith 4996d7c1e57SBarry Smith Level: advanced 5006d7c1e57SBarry Smith 5016d7c1e57SBarry Smith Notes: 5026d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 5036d7c1e57SBarry Smith matrix in the composite matrix. 5046d7c1e57SBarry Smith 5056d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 5066d7c1e57SBarry Smith 5076d7c1e57SBarry Smith @*/ 5087087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 5096d7c1e57SBarry Smith { 5106d7c1e57SBarry Smith Mat_Composite *b = (Mat_Composite*)mat->data; 511ace3abfcSBarry Smith PetscBool flg; 5126d7c1e57SBarry Smith PetscErrorCode ierr; 5136d7c1e57SBarry Smith 5146d7c1e57SBarry Smith PetscFunctionBegin; 515251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATCOMPOSITE,&flg);CHKERRQ(ierr); 516e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Can only use with composite matrix"); 5176d7c1e57SBarry Smith if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 5186d7c1e57SBarry Smith mat->ops->getdiagonal = 0; 5196d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite_Multiplicative; 5206d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 5216d7c1e57SBarry Smith b->type = MAT_COMPOSITE_MULTIPLICATIVE; 5226d7c1e57SBarry Smith } else { 5236d7c1e57SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Composite; 5246d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite; 5256d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite; 5266d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 527793850ffSBarry Smith } 528793850ffSBarry Smith PetscFunctionReturn(0); 529793850ffSBarry Smith } 530793850ffSBarry Smith 5316d7c1e57SBarry Smith 532b52f573bSBarry Smith #undef __FUNCT__ 533b52f573bSBarry Smith #define __FUNCT__ "MatCompositeMerge" 534b52f573bSBarry Smith /*@C 535b52f573bSBarry Smith MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 536b52f573bSBarry Smith by summing all the matrices inside the composite matrix. 537b52f573bSBarry Smith 538b52f573bSBarry Smith Collective on MPI_Comm 539b52f573bSBarry Smith 540b52f573bSBarry Smith Input Parameters: 541b52f573bSBarry Smith . mat - the composite matrix 542b52f573bSBarry Smith 543b52f573bSBarry Smith 544b52f573bSBarry Smith Options Database: 545b52f573bSBarry Smith . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 546b52f573bSBarry Smith 547b52f573bSBarry Smith Level: advanced 548b52f573bSBarry Smith 549b52f573bSBarry Smith Notes: 550b52f573bSBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 551b52f573bSBarry Smith matrix in the composite matrix. 552b52f573bSBarry Smith 553b52f573bSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 554b52f573bSBarry Smith 555b52f573bSBarry Smith @*/ 5567087cfbeSBarry Smith PetscErrorCode MatCompositeMerge(Mat mat) 557b52f573bSBarry Smith { 558b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 5596d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 560b52f573bSBarry Smith PetscErrorCode ierr; 5616d7c1e57SBarry Smith Mat tmat,newmat; 562b52f573bSBarry Smith 563b52f573bSBarry Smith PetscFunctionBegin; 564e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 565b52f573bSBarry Smith 566b52f573bSBarry Smith PetscFunctionBegin; 5676d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 568b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 569b52f573bSBarry Smith while ((next = next->next)) { 570b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 571b52f573bSBarry Smith } 5726d7c1e57SBarry Smith } else { 5736d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 5746d7c1e57SBarry Smith while ((prev = prev->prev)) { 5756d7c1e57SBarry Smith ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 5766bf464f9SBarry Smith ierr = MatDestroy(&tmat);CHKERRQ(ierr); 5776d7c1e57SBarry Smith tmat = newmat; 5786d7c1e57SBarry Smith } 5796d7c1e57SBarry Smith } 5806d7c1e57SBarry Smith ierr = MatHeaderReplace(mat,tmat);CHKERRQ(ierr); 581e4fc5df0SBarry Smith ierr = MatDiagonalScale(mat,shell->left,shell->right);CHKERRQ(ierr); 582e4fc5df0SBarry Smith ierr = MatScale(mat,shell->scale);CHKERRQ(ierr); 583b52f573bSBarry Smith PetscFunctionReturn(0); 584b52f573bSBarry Smith } 585