1793850ffSBarry Smith 2*c6db04a5SJed Brown #include <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) { 30793850ffSBarry Smith ierr = MatDestroy(next->mat);CHKERRQ(ierr); 316d7c1e57SBarry Smith if (next->work && (!next->next || next->work != next->next->work)) { 326d7c1e57SBarry 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 } 38793850ffSBarry Smith if (shell->work) {ierr = VecDestroy(shell->work);CHKERRQ(ierr);} 39e4fc5df0SBarry Smith if (shell->left) {ierr = VecDestroy(shell->left);CHKERRQ(ierr);} 40e4fc5df0SBarry Smith if (shell->right) {ierr = VecDestroy(shell->right);CHKERRQ(ierr);} 41e4fc5df0SBarry Smith if (shell->leftwork) {ierr = VecDestroy(shell->leftwork);CHKERRQ(ierr);} 42e4fc5df0SBarry Smith if (shell->rightwork) {ierr = VecDestroy(shell->rightwork);CHKERRQ(ierr);} 43793850ffSBarry Smith ierr = PetscFree(shell);CHKERRQ(ierr); 44793850ffSBarry Smith mat->data = 0; 45793850ffSBarry Smith PetscFunctionReturn(0); 46793850ffSBarry Smith } 47793850ffSBarry Smith 48793850ffSBarry Smith #undef __FUNCT__ 496d7c1e57SBarry Smith #define __FUNCT__ "MatMult_Composite_Multiplicative" 506d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 516d7c1e57SBarry Smith { 526d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 536d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 546d7c1e57SBarry Smith PetscErrorCode ierr; 556d7c1e57SBarry Smith Vec in,out; 566d7c1e57SBarry Smith 576d7c1e57SBarry Smith PetscFunctionBegin; 58e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 596d7c1e57SBarry Smith in = x; 60e4fc5df0SBarry Smith if (shell->right) { 61e4fc5df0SBarry Smith if (!shell->rightwork) { 62e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 63e4fc5df0SBarry Smith } 64e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 65e4fc5df0SBarry Smith in = shell->rightwork; 66e4fc5df0SBarry Smith } 676d7c1e57SBarry Smith while (next->next) { 686d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 696d7c1e57SBarry Smith ierr = MatGetVecs(next->mat,PETSC_NULL,&next->work);CHKERRQ(ierr); 706d7c1e57SBarry Smith } 716d7c1e57SBarry Smith out = next->work; 726d7c1e57SBarry Smith ierr = MatMult(next->mat,in,out);CHKERRQ(ierr); 736d7c1e57SBarry Smith in = out; 746d7c1e57SBarry Smith next = next->next; 756d7c1e57SBarry Smith } 766d7c1e57SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 77e4fc5df0SBarry Smith if (shell->left) { 78e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 79e4fc5df0SBarry Smith } 80e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 816d7c1e57SBarry Smith PetscFunctionReturn(0); 826d7c1e57SBarry Smith } 836d7c1e57SBarry Smith 846d7c1e57SBarry Smith #undef __FUNCT__ 856d7c1e57SBarry Smith #define __FUNCT__ "MatMultTranspose_Composite_Multiplicative" 866d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 876d7c1e57SBarry Smith { 886d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 896d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 906d7c1e57SBarry Smith PetscErrorCode ierr; 916d7c1e57SBarry Smith Vec in,out; 926d7c1e57SBarry Smith 936d7c1e57SBarry Smith PetscFunctionBegin; 94e32f2f54SBarry Smith if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 956d7c1e57SBarry Smith in = x; 96e4fc5df0SBarry Smith if (shell->left) { 97e4fc5df0SBarry Smith if (!shell->leftwork) { 98e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 99e4fc5df0SBarry Smith } 100e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 101e4fc5df0SBarry Smith in = shell->leftwork; 102e4fc5df0SBarry Smith } 1036d7c1e57SBarry Smith while (tail->prev) { 1046d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1056d7c1e57SBarry Smith ierr = MatGetVecs(tail->mat,PETSC_NULL,&tail->prev->work);CHKERRQ(ierr); 1066d7c1e57SBarry Smith } 1076d7c1e57SBarry Smith out = tail->prev->work; 1086d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1096d7c1e57SBarry Smith in = out; 1106d7c1e57SBarry Smith tail = tail->prev; 1116d7c1e57SBarry Smith } 1126d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 113e4fc5df0SBarry Smith if (shell->right) { 114e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 115e4fc5df0SBarry Smith } 116e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1176d7c1e57SBarry Smith PetscFunctionReturn(0); 1186d7c1e57SBarry Smith } 1196d7c1e57SBarry Smith 1206d7c1e57SBarry Smith #undef __FUNCT__ 121793850ffSBarry Smith #define __FUNCT__ "MatMult_Composite" 122793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 123793850ffSBarry Smith { 124793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 125793850ffSBarry Smith Mat_CompositeLink next = shell->head; 126793850ffSBarry Smith PetscErrorCode ierr; 127e4fc5df0SBarry Smith Vec in; 128793850ffSBarry Smith 129793850ffSBarry Smith PetscFunctionBegin; 130e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 131e4fc5df0SBarry Smith in = x; 132e4fc5df0SBarry Smith if (shell->right) { 133e4fc5df0SBarry Smith if (!shell->rightwork) { 134e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 135793850ffSBarry Smith } 136e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 137e4fc5df0SBarry Smith in = shell->rightwork; 138e4fc5df0SBarry Smith } 139e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 140e4fc5df0SBarry Smith while ((next = next->next)) { 141e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 142e4fc5df0SBarry Smith } 143e4fc5df0SBarry Smith if (shell->left) { 144e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 145e4fc5df0SBarry Smith } 146e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 147793850ffSBarry Smith PetscFunctionReturn(0); 148793850ffSBarry Smith } 149793850ffSBarry Smith 150793850ffSBarry Smith #undef __FUNCT__ 151793850ffSBarry Smith #define __FUNCT__ "MatMultTranspose_Composite" 152793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 153793850ffSBarry Smith { 154793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 155793850ffSBarry Smith Mat_CompositeLink next = shell->head; 156793850ffSBarry Smith PetscErrorCode ierr; 157e4fc5df0SBarry Smith Vec in; 158793850ffSBarry Smith 159793850ffSBarry Smith PetscFunctionBegin; 160e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 161e4fc5df0SBarry Smith in = x; 162e4fc5df0SBarry Smith if (shell->left) { 163e4fc5df0SBarry Smith if (!shell->leftwork) { 164e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 165793850ffSBarry Smith } 166e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 167e4fc5df0SBarry Smith in = shell->leftwork; 168e4fc5df0SBarry Smith } 169e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 170e4fc5df0SBarry Smith while ((next = next->next)) { 171e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 172e4fc5df0SBarry Smith } 173e4fc5df0SBarry Smith if (shell->right) { 174e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 175e4fc5df0SBarry Smith } 176e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 177793850ffSBarry Smith PetscFunctionReturn(0); 178793850ffSBarry Smith } 179793850ffSBarry Smith 180793850ffSBarry Smith #undef __FUNCT__ 181793850ffSBarry Smith #define __FUNCT__ "MatGetDiagonal_Composite" 182793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 183793850ffSBarry Smith { 184793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 185793850ffSBarry Smith Mat_CompositeLink next = shell->head; 186793850ffSBarry Smith PetscErrorCode ierr; 187793850ffSBarry Smith 188793850ffSBarry Smith PetscFunctionBegin; 189e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 190e32f2f54SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 191e4fc5df0SBarry Smith 192793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 193793850ffSBarry Smith if (next->next && !shell->work) { 194793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 195793850ffSBarry Smith } 196793850ffSBarry Smith while ((next = next->next)) { 197793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 198793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 199793850ffSBarry Smith } 200e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 201793850ffSBarry Smith PetscFunctionReturn(0); 202793850ffSBarry Smith } 203793850ffSBarry Smith 204793850ffSBarry Smith #undef __FUNCT__ 205793850ffSBarry Smith #define __FUNCT__ "MatAssemblyEnd_Composite" 206793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 207793850ffSBarry Smith { 208b52f573bSBarry Smith PetscErrorCode ierr; 209ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 210b52f573bSBarry Smith 211793850ffSBarry Smith PetscFunctionBegin; 212acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)Y)->prefix,"-mat_composite_merge",&flg,PETSC_NULL);CHKERRQ(ierr); 213b52f573bSBarry Smith if (flg) { 214b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 215b52f573bSBarry Smith } 216793850ffSBarry Smith PetscFunctionReturn(0); 217793850ffSBarry Smith } 218793850ffSBarry Smith 219e4fc5df0SBarry Smith #undef __FUNCT__ 220e4fc5df0SBarry Smith #define __FUNCT__ "MatScale_Composite" 221e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 222e4fc5df0SBarry Smith { 223e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 224e4fc5df0SBarry Smith 225e4fc5df0SBarry Smith PetscFunctionBegin; 226321429dbSBarry Smith a->scale *= alpha; 227e4fc5df0SBarry Smith PetscFunctionReturn(0); 228e4fc5df0SBarry Smith } 229e4fc5df0SBarry Smith 230e4fc5df0SBarry Smith #undef __FUNCT__ 231e4fc5df0SBarry Smith #define __FUNCT__ "MatDiagonalScale_Composite" 232e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 233e4fc5df0SBarry Smith { 234e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 235e4fc5df0SBarry Smith PetscErrorCode ierr; 236e4fc5df0SBarry Smith 237e4fc5df0SBarry Smith PetscFunctionBegin; 238e4fc5df0SBarry Smith if (left) { 239321429dbSBarry Smith if (!a->left) { 240e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 241e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 242321429dbSBarry Smith } else { 243321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 244321429dbSBarry Smith } 245e4fc5df0SBarry Smith } 246e4fc5df0SBarry Smith if (right) { 247321429dbSBarry Smith if (!a->right) { 248e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 249e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 250321429dbSBarry Smith } else { 251321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 252321429dbSBarry Smith } 253e4fc5df0SBarry Smith } 254e4fc5df0SBarry Smith PetscFunctionReturn(0); 255e4fc5df0SBarry Smith } 256793850ffSBarry Smith 257793850ffSBarry Smith static struct _MatOps MatOps_Values = {0, 258793850ffSBarry Smith 0, 259793850ffSBarry Smith 0, 260793850ffSBarry Smith MatMult_Composite, 261793850ffSBarry Smith 0, 262793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite, 263793850ffSBarry Smith 0, 264793850ffSBarry Smith 0, 265793850ffSBarry Smith 0, 266793850ffSBarry Smith 0, 267793850ffSBarry Smith /*10*/ 0, 268793850ffSBarry Smith 0, 269793850ffSBarry Smith 0, 270793850ffSBarry Smith 0, 271793850ffSBarry Smith 0, 272793850ffSBarry Smith /*15*/ 0, 273793850ffSBarry Smith 0, 274793850ffSBarry Smith MatGetDiagonal_Composite, 275e4fc5df0SBarry Smith MatDiagonalScale_Composite, 276793850ffSBarry Smith 0, 277793850ffSBarry Smith /*20*/ 0, 278793850ffSBarry Smith MatAssemblyEnd_Composite, 279793850ffSBarry Smith 0, 280793850ffSBarry Smith 0, 281d519adbfSMatthew Knepley /*24*/ 0, 282793850ffSBarry Smith 0, 283793850ffSBarry Smith 0, 284793850ffSBarry Smith 0, 285793850ffSBarry Smith 0, 286d519adbfSMatthew Knepley /*29*/ 0, 287793850ffSBarry Smith 0, 288793850ffSBarry Smith 0, 289793850ffSBarry Smith 0, 290793850ffSBarry Smith 0, 291d519adbfSMatthew Knepley /*34*/ 0, 292793850ffSBarry Smith 0, 293793850ffSBarry Smith 0, 294793850ffSBarry Smith 0, 295793850ffSBarry Smith 0, 296d519adbfSMatthew Knepley /*39*/ 0, 297793850ffSBarry Smith 0, 298793850ffSBarry Smith 0, 299793850ffSBarry Smith 0, 300793850ffSBarry Smith 0, 301d519adbfSMatthew Knepley /*44*/ 0, 302e4fc5df0SBarry Smith MatScale_Composite, 303793850ffSBarry Smith 0, 304793850ffSBarry Smith 0, 305793850ffSBarry Smith 0, 306d519adbfSMatthew Knepley /*49*/ 0, 307793850ffSBarry Smith 0, 308793850ffSBarry Smith 0, 309793850ffSBarry Smith 0, 310793850ffSBarry Smith 0, 311d519adbfSMatthew Knepley /*54*/ 0, 312793850ffSBarry Smith 0, 313793850ffSBarry Smith 0, 314793850ffSBarry Smith 0, 315793850ffSBarry Smith 0, 316d519adbfSMatthew Knepley /*59*/ 0, 317793850ffSBarry Smith MatDestroy_Composite, 318793850ffSBarry Smith 0, 319793850ffSBarry Smith 0, 320793850ffSBarry Smith 0, 321d519adbfSMatthew Knepley /*64*/ 0, 322793850ffSBarry Smith 0, 323793850ffSBarry Smith 0, 324793850ffSBarry Smith 0, 325793850ffSBarry Smith 0, 326d519adbfSMatthew Knepley /*69*/ 0, 327793850ffSBarry Smith 0, 328793850ffSBarry Smith 0, 329793850ffSBarry Smith 0, 330793850ffSBarry Smith 0, 331d519adbfSMatthew Knepley /*74*/ 0, 332793850ffSBarry Smith 0, 333793850ffSBarry Smith 0, 334793850ffSBarry Smith 0, 335793850ffSBarry Smith 0, 336d519adbfSMatthew Knepley /*79*/ 0, 337793850ffSBarry Smith 0, 338793850ffSBarry Smith 0, 339793850ffSBarry Smith 0, 340793850ffSBarry Smith 0, 341d519adbfSMatthew Knepley /*84*/ 0, 342793850ffSBarry Smith 0, 343793850ffSBarry Smith 0, 344793850ffSBarry Smith 0, 345793850ffSBarry Smith 0, 346d519adbfSMatthew Knepley /*89*/ 0, 347793850ffSBarry Smith 0, 348793850ffSBarry Smith 0, 349793850ffSBarry Smith 0, 350793850ffSBarry Smith 0, 351d519adbfSMatthew Knepley /*94*/ 0, 352793850ffSBarry Smith 0, 353793850ffSBarry Smith 0, 354793850ffSBarry Smith 0}; 355793850ffSBarry Smith 356793850ffSBarry Smith /*MC 3576d7c1e57SBarry Smith MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices (all matrices are of same size and parallel layout). 3586d7c1e57SBarry Smith 3596d7c1e57SBarry Smith Notes: to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 360793850ffSBarry Smith 361793850ffSBarry Smith Level: advanced 362793850ffSBarry Smith 3636d7c1e57SBarry Smith .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 364793850ffSBarry Smith M*/ 365793850ffSBarry Smith 366793850ffSBarry Smith EXTERN_C_BEGIN 367793850ffSBarry Smith #undef __FUNCT__ 368793850ffSBarry Smith #define __FUNCT__ "MatCreate_Composite" 3697087cfbeSBarry Smith PetscErrorCode MatCreate_Composite(Mat A) 370793850ffSBarry Smith { 371793850ffSBarry Smith Mat_Composite *b; 372793850ffSBarry Smith PetscErrorCode ierr; 373793850ffSBarry Smith 374793850ffSBarry Smith PetscFunctionBegin; 37538f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Composite,&b);CHKERRQ(ierr); 376793850ffSBarry Smith A->data = (void*)b; 37738f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 378793850ffSBarry Smith 37926283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 38026283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 38126283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 38226283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 383793850ffSBarry Smith 384793850ffSBarry Smith A->assembled = PETSC_TRUE; 385793850ffSBarry Smith A->preallocated = PETSC_FALSE; 3866d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 387e4fc5df0SBarry Smith b->scale = 1.0; 388793850ffSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 389793850ffSBarry Smith PetscFunctionReturn(0); 390793850ffSBarry Smith } 391793850ffSBarry Smith EXTERN_C_END 392793850ffSBarry Smith 393793850ffSBarry Smith #undef __FUNCT__ 394793850ffSBarry Smith #define __FUNCT__ "MatCreateComposite" 395793850ffSBarry Smith /*@C 396793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 397793850ffSBarry Smith 398793850ffSBarry Smith Collective on MPI_Comm 399793850ffSBarry Smith 400793850ffSBarry Smith Input Parameters: 401793850ffSBarry Smith + comm - MPI communicator 402793850ffSBarry Smith . nmat - number of matrices to put in 403793850ffSBarry Smith - mats - the matrices 404793850ffSBarry Smith 405793850ffSBarry Smith Output Parameter: 406793850ffSBarry Smith . A - the matrix 407793850ffSBarry Smith 408793850ffSBarry Smith Level: advanced 409793850ffSBarry Smith 410793850ffSBarry Smith Notes: 411793850ffSBarry Smith Alternative construction 412793850ffSBarry Smith $ MatCreate(comm,&mat); 413793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 414793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 415793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 416793850ffSBarry Smith $ .... 417793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 418b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 419b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 420793850ffSBarry Smith 4216d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4226d7c1e57SBarry Smith 4236d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 424793850ffSBarry Smith 425793850ffSBarry Smith @*/ 4267087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 427793850ffSBarry Smith { 428793850ffSBarry Smith PetscErrorCode ierr; 429793850ffSBarry Smith PetscInt m,n,M,N,i; 430793850ffSBarry Smith 431793850ffSBarry Smith PetscFunctionBegin; 432e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 433f3f84630SBarry Smith PetscValidPointer(mat,3); 434793850ffSBarry Smith 435793850ffSBarry Smith ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr); 436793850ffSBarry Smith ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr); 437793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 438793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 439793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 440793850ffSBarry Smith for (i=0; i<nmat; i++) { 441793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 442793850ffSBarry Smith } 443b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 444b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 445793850ffSBarry Smith PetscFunctionReturn(0); 446793850ffSBarry Smith } 447793850ffSBarry Smith 448793850ffSBarry Smith #undef __FUNCT__ 449793850ffSBarry Smith #define __FUNCT__ "MatCompositeAddMat" 450793850ffSBarry Smith /*@ 451793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 452793850ffSBarry Smith 453793850ffSBarry Smith Collective on Mat 454793850ffSBarry Smith 455793850ffSBarry Smith Input Parameters: 456793850ffSBarry Smith + mat - the composite matrix 457793850ffSBarry Smith - smat - the partial matrix 458793850ffSBarry Smith 459793850ffSBarry Smith Level: advanced 460793850ffSBarry Smith 461793850ffSBarry Smith .seealso: MatCreateComposite() 462793850ffSBarry Smith @*/ 4637087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 464793850ffSBarry Smith { 46538f2d2fdSLisandro Dalcin Mat_Composite *shell; 466793850ffSBarry Smith PetscErrorCode ierr; 467793850ffSBarry Smith Mat_CompositeLink ilink,next; 468793850ffSBarry Smith 469793850ffSBarry Smith PetscFunctionBegin; 4700700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4710700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 47238f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr); 473793850ffSBarry Smith ilink->next = 0; 474793850ffSBarry Smith ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 475c3122656SLisandro Dalcin ilink->mat = smat; 476793850ffSBarry Smith 47738f2d2fdSLisandro Dalcin shell = (Mat_Composite*)mat->data; 478793850ffSBarry Smith next = shell->head; 479793850ffSBarry Smith if (!next) { 480793850ffSBarry Smith shell->head = ilink; 481793850ffSBarry Smith } else { 482793850ffSBarry Smith while (next->next) { 483793850ffSBarry Smith next = next->next; 484793850ffSBarry Smith } 485793850ffSBarry Smith next->next = ilink; 4866d7c1e57SBarry Smith ilink->prev = next; 4876d7c1e57SBarry Smith } 4886d7c1e57SBarry Smith shell->tail = ilink; 4896d7c1e57SBarry Smith PetscFunctionReturn(0); 4906d7c1e57SBarry Smith } 4916d7c1e57SBarry Smith 4926d7c1e57SBarry Smith #undef __FUNCT__ 4936d7c1e57SBarry Smith #define __FUNCT__ "MatCompositeSetType" 4946d7c1e57SBarry Smith /*@C 4956d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 4966d7c1e57SBarry Smith 4976d7c1e57SBarry Smith Collective on MPI_Comm 4986d7c1e57SBarry Smith 4996d7c1e57SBarry Smith Input Parameters: 5006d7c1e57SBarry Smith . mat - the composite matrix 5016d7c1e57SBarry Smith 5026d7c1e57SBarry Smith 5036d7c1e57SBarry Smith Level: advanced 5046d7c1e57SBarry Smith 5056d7c1e57SBarry Smith Notes: 5066d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 5076d7c1e57SBarry Smith matrix in the composite matrix. 5086d7c1e57SBarry Smith 5096d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 5106d7c1e57SBarry Smith 5116d7c1e57SBarry Smith @*/ 5127087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 5136d7c1e57SBarry Smith { 5146d7c1e57SBarry Smith Mat_Composite *b = (Mat_Composite*)mat->data; 515ace3abfcSBarry Smith PetscBool flg; 5166d7c1e57SBarry Smith PetscErrorCode ierr; 5176d7c1e57SBarry Smith 5186d7c1e57SBarry Smith PetscFunctionBegin; 5196d7c1e57SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATCOMPOSITE,&flg);CHKERRQ(ierr); 520e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Can only use with composite matrix"); 5216d7c1e57SBarry Smith if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 5226d7c1e57SBarry Smith mat->ops->getdiagonal = 0; 5236d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite_Multiplicative; 5246d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 5256d7c1e57SBarry Smith b->type = MAT_COMPOSITE_MULTIPLICATIVE; 5266d7c1e57SBarry Smith } else { 5276d7c1e57SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Composite; 5286d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite; 5296d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite; 5306d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 531793850ffSBarry Smith } 532793850ffSBarry Smith PetscFunctionReturn(0); 533793850ffSBarry Smith } 534793850ffSBarry Smith 5356d7c1e57SBarry Smith 536b52f573bSBarry Smith #undef __FUNCT__ 537b52f573bSBarry Smith #define __FUNCT__ "MatCompositeMerge" 538b52f573bSBarry Smith /*@C 539b52f573bSBarry Smith MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 540b52f573bSBarry Smith by summing all the matrices inside the composite matrix. 541b52f573bSBarry Smith 542b52f573bSBarry Smith Collective on MPI_Comm 543b52f573bSBarry Smith 544b52f573bSBarry Smith Input Parameters: 545b52f573bSBarry Smith . mat - the composite matrix 546b52f573bSBarry Smith 547b52f573bSBarry Smith 548b52f573bSBarry Smith Options Database: 549b52f573bSBarry Smith . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 550b52f573bSBarry Smith 551b52f573bSBarry Smith Level: advanced 552b52f573bSBarry Smith 553b52f573bSBarry Smith Notes: 554b52f573bSBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 555b52f573bSBarry Smith matrix in the composite matrix. 556b52f573bSBarry Smith 557b52f573bSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 558b52f573bSBarry Smith 559b52f573bSBarry Smith @*/ 5607087cfbeSBarry Smith PetscErrorCode MatCompositeMerge(Mat mat) 561b52f573bSBarry Smith { 562b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 5636d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 564b52f573bSBarry Smith PetscErrorCode ierr; 5656d7c1e57SBarry Smith Mat tmat,newmat; 566b52f573bSBarry Smith 567b52f573bSBarry Smith PetscFunctionBegin; 568e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 569b52f573bSBarry Smith 570b52f573bSBarry Smith PetscFunctionBegin; 5716d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 572b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 573b52f573bSBarry Smith while ((next = next->next)) { 574b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 575b52f573bSBarry Smith } 5766d7c1e57SBarry Smith } else { 5776d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 5786d7c1e57SBarry Smith while ((prev = prev->prev)) { 5796d7c1e57SBarry Smith ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 5806d7c1e57SBarry Smith ierr = MatDestroy(tmat);CHKERRQ(ierr); 5816d7c1e57SBarry Smith tmat = newmat; 5826d7c1e57SBarry Smith } 5836d7c1e57SBarry Smith } 5846d7c1e57SBarry Smith ierr = MatHeaderReplace(mat,tmat);CHKERRQ(ierr); 585e4fc5df0SBarry Smith ierr = MatDiagonalScale(mat,shell->left,shell->right);CHKERRQ(ierr); 586e4fc5df0SBarry Smith ierr = MatScale(mat,shell->scale);CHKERRQ(ierr); 587b52f573bSBarry Smith PetscFunctionReturn(0); 588b52f573bSBarry Smith } 589