1793850ffSBarry Smith #define PETSCMAT_DLL 2793850ffSBarry Smith 37c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 4793850ffSBarry Smith 5793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink; 6793850ffSBarry Smith struct _Mat_CompositeLink { 7793850ffSBarry Smith Mat mat; 86d7c1e57SBarry Smith Vec work; 96d7c1e57SBarry Smith Mat_CompositeLink next,prev; 10793850ffSBarry Smith }; 11793850ffSBarry Smith 12793850ffSBarry Smith typedef struct { 136d7c1e57SBarry Smith MatCompositeType type; 146d7c1e57SBarry Smith Mat_CompositeLink head,tail; 15793850ffSBarry Smith Vec work; 16e4fc5df0SBarry Smith PetscScalar scale; /* scale factor supplied with MatScale() */ 17e4fc5df0SBarry Smith Vec left,right; /* left and right diagonal scaling provided with MatDiagonalScale() */ 18e4fc5df0SBarry Smith Vec leftwork,rightwork; 19793850ffSBarry Smith } Mat_Composite; 20793850ffSBarry Smith 21793850ffSBarry Smith #undef __FUNCT__ 22793850ffSBarry Smith #define __FUNCT__ "MatDestroy_Composite" 23793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat) 24793850ffSBarry Smith { 25793850ffSBarry Smith PetscErrorCode ierr; 262c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite*)mat->data; 276d7c1e57SBarry Smith Mat_CompositeLink next = shell->head,oldnext; 28793850ffSBarry Smith 29793850ffSBarry Smith PetscFunctionBegin; 30793850ffSBarry Smith while (next) { 31793850ffSBarry Smith ierr = MatDestroy(next->mat);CHKERRQ(ierr); 326d7c1e57SBarry Smith if (next->work && (!next->next || next->work != next->next->work)) { 336d7c1e57SBarry Smith ierr = VecDestroy(next->work);CHKERRQ(ierr); 346d7c1e57SBarry Smith } 356d7c1e57SBarry Smith oldnext = next; 36793850ffSBarry Smith next = next->next; 376d7c1e57SBarry Smith ierr = PetscFree(oldnext);CHKERRQ(ierr); 38793850ffSBarry Smith } 39793850ffSBarry Smith if (shell->work) {ierr = VecDestroy(shell->work);CHKERRQ(ierr);} 40e4fc5df0SBarry Smith if (shell->left) {ierr = VecDestroy(shell->left);CHKERRQ(ierr);} 41e4fc5df0SBarry Smith if (shell->right) {ierr = VecDestroy(shell->right);CHKERRQ(ierr);} 42e4fc5df0SBarry Smith if (shell->leftwork) {ierr = VecDestroy(shell->leftwork);CHKERRQ(ierr);} 43e4fc5df0SBarry Smith if (shell->rightwork) {ierr = VecDestroy(shell->rightwork);CHKERRQ(ierr);} 44793850ffSBarry Smith ierr = PetscFree(shell);CHKERRQ(ierr); 45793850ffSBarry Smith mat->data = 0; 46793850ffSBarry Smith PetscFunctionReturn(0); 47793850ffSBarry Smith } 48793850ffSBarry Smith 49793850ffSBarry Smith #undef __FUNCT__ 506d7c1e57SBarry Smith #define __FUNCT__ "MatMult_Composite_Multiplicative" 516d7c1e57SBarry Smith PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y) 526d7c1e57SBarry Smith { 536d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 546d7c1e57SBarry Smith Mat_CompositeLink next = shell->head; 556d7c1e57SBarry Smith PetscErrorCode ierr; 566d7c1e57SBarry Smith Vec in,out; 576d7c1e57SBarry Smith 586d7c1e57SBarry Smith PetscFunctionBegin; 596d7c1e57SBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 606d7c1e57SBarry Smith in = x; 61e4fc5df0SBarry Smith if (shell->right) { 62e4fc5df0SBarry Smith if (!shell->rightwork) { 63e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 64e4fc5df0SBarry Smith } 65e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 66e4fc5df0SBarry Smith in = shell->rightwork; 67e4fc5df0SBarry Smith } 686d7c1e57SBarry Smith while (next->next) { 696d7c1e57SBarry Smith if (!next->work) { /* should reuse previous work if the same size */ 706d7c1e57SBarry Smith ierr = MatGetVecs(next->mat,PETSC_NULL,&next->work);CHKERRQ(ierr); 716d7c1e57SBarry Smith } 726d7c1e57SBarry Smith out = next->work; 736d7c1e57SBarry Smith ierr = MatMult(next->mat,in,out);CHKERRQ(ierr); 746d7c1e57SBarry Smith in = out; 756d7c1e57SBarry Smith next = next->next; 766d7c1e57SBarry Smith } 776d7c1e57SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 78e4fc5df0SBarry Smith if (shell->left) { 79e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 80e4fc5df0SBarry Smith } 81e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 826d7c1e57SBarry Smith PetscFunctionReturn(0); 836d7c1e57SBarry Smith } 846d7c1e57SBarry Smith 856d7c1e57SBarry Smith #undef __FUNCT__ 866d7c1e57SBarry Smith #define __FUNCT__ "MatMultTranspose_Composite_Multiplicative" 876d7c1e57SBarry Smith PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 886d7c1e57SBarry Smith { 896d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 906d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 916d7c1e57SBarry Smith PetscErrorCode ierr; 926d7c1e57SBarry Smith Vec in,out; 936d7c1e57SBarry Smith 946d7c1e57SBarry Smith PetscFunctionBegin; 956d7c1e57SBarry Smith if (!tail) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 966d7c1e57SBarry Smith in = x; 97e4fc5df0SBarry Smith if (shell->left) { 98e4fc5df0SBarry Smith if (!shell->leftwork) { 99e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 100e4fc5df0SBarry Smith } 101e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 102e4fc5df0SBarry Smith in = shell->leftwork; 103e4fc5df0SBarry Smith } 1046d7c1e57SBarry Smith while (tail->prev) { 1056d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1066d7c1e57SBarry Smith ierr = MatGetVecs(tail->mat,PETSC_NULL,&tail->prev->work);CHKERRQ(ierr); 1076d7c1e57SBarry Smith } 1086d7c1e57SBarry Smith out = tail->prev->work; 1096d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1106d7c1e57SBarry Smith in = out; 1116d7c1e57SBarry Smith tail = tail->prev; 1126d7c1e57SBarry Smith } 1136d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 114e4fc5df0SBarry Smith if (shell->right) { 115e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 116e4fc5df0SBarry Smith } 117e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1186d7c1e57SBarry Smith PetscFunctionReturn(0); 1196d7c1e57SBarry Smith } 1206d7c1e57SBarry Smith 1216d7c1e57SBarry Smith #undef __FUNCT__ 122793850ffSBarry Smith #define __FUNCT__ "MatMult_Composite" 123793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 124793850ffSBarry Smith { 125793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 126793850ffSBarry Smith Mat_CompositeLink next = shell->head; 127793850ffSBarry Smith PetscErrorCode ierr; 128e4fc5df0SBarry Smith Vec in; 129793850ffSBarry Smith 130793850ffSBarry Smith PetscFunctionBegin; 131793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 132e4fc5df0SBarry Smith in = x; 133e4fc5df0SBarry Smith if (shell->right) { 134e4fc5df0SBarry Smith if (!shell->rightwork) { 135e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 136793850ffSBarry Smith } 137e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 138e4fc5df0SBarry Smith in = shell->rightwork; 139e4fc5df0SBarry Smith } 140e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 141e4fc5df0SBarry Smith while ((next = next->next)) { 142e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 143e4fc5df0SBarry Smith } 144e4fc5df0SBarry Smith if (shell->left) { 145e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 146e4fc5df0SBarry Smith } 147e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 148793850ffSBarry Smith PetscFunctionReturn(0); 149793850ffSBarry Smith } 150793850ffSBarry Smith 151793850ffSBarry Smith #undef __FUNCT__ 152793850ffSBarry Smith #define __FUNCT__ "MatMultTranspose_Composite" 153793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 154793850ffSBarry Smith { 155793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 156793850ffSBarry Smith Mat_CompositeLink next = shell->head; 157793850ffSBarry Smith PetscErrorCode ierr; 158e4fc5df0SBarry Smith Vec in; 159793850ffSBarry Smith 160793850ffSBarry Smith PetscFunctionBegin; 161793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 162e4fc5df0SBarry Smith in = x; 163e4fc5df0SBarry Smith if (shell->left) { 164e4fc5df0SBarry Smith if (!shell->leftwork) { 165e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 166793850ffSBarry Smith } 167e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 168e4fc5df0SBarry Smith in = shell->leftwork; 169e4fc5df0SBarry Smith } 170e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 171e4fc5df0SBarry Smith while ((next = next->next)) { 172e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 173e4fc5df0SBarry Smith } 174e4fc5df0SBarry Smith if (shell->right) { 175e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 176e4fc5df0SBarry Smith } 177e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 178793850ffSBarry Smith PetscFunctionReturn(0); 179793850ffSBarry Smith } 180793850ffSBarry Smith 181793850ffSBarry Smith #undef __FUNCT__ 182793850ffSBarry Smith #define __FUNCT__ "MatGetDiagonal_Composite" 183793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 184793850ffSBarry Smith { 185793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 186793850ffSBarry Smith Mat_CompositeLink next = shell->head; 187793850ffSBarry Smith PetscErrorCode ierr; 188793850ffSBarry Smith 189793850ffSBarry Smith PetscFunctionBegin; 190793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 191e4fc5df0SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 192e4fc5df0SBarry Smith 193793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 194793850ffSBarry Smith if (next->next && !shell->work) { 195793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 196793850ffSBarry Smith } 197793850ffSBarry Smith while ((next = next->next)) { 198793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 199793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 200793850ffSBarry Smith } 201e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 202793850ffSBarry Smith PetscFunctionReturn(0); 203793850ffSBarry Smith } 204793850ffSBarry Smith 205793850ffSBarry Smith #undef __FUNCT__ 206793850ffSBarry Smith #define __FUNCT__ "MatAssemblyEnd_Composite" 207793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 208793850ffSBarry Smith { 209b52f573bSBarry Smith PetscErrorCode ierr; 210b52f573bSBarry Smith PetscTruth flg; 211b52f573bSBarry Smith 212793850ffSBarry Smith PetscFunctionBegin; 2137adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)Y)->prefix,"-mat_composite_merge",&flg);CHKERRQ(ierr); 214b52f573bSBarry Smith if (flg) { 215b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 216b52f573bSBarry Smith } 217793850ffSBarry Smith PetscFunctionReturn(0); 218793850ffSBarry Smith } 219793850ffSBarry Smith 220e4fc5df0SBarry Smith #undef __FUNCT__ 221e4fc5df0SBarry Smith #define __FUNCT__ "MatScale_Composite" 222e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 223e4fc5df0SBarry Smith { 224e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 225e4fc5df0SBarry Smith 226e4fc5df0SBarry Smith PetscFunctionBegin; 227*321429dbSBarry Smith a->scale *= alpha; 228e4fc5df0SBarry Smith PetscFunctionReturn(0); 229e4fc5df0SBarry Smith } 230e4fc5df0SBarry Smith 231e4fc5df0SBarry Smith #undef __FUNCT__ 232e4fc5df0SBarry Smith #define __FUNCT__ "MatDiagonalScale_Composite" 233e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 234e4fc5df0SBarry Smith { 235e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 236e4fc5df0SBarry Smith PetscErrorCode ierr; 237e4fc5df0SBarry Smith 238e4fc5df0SBarry Smith PetscFunctionBegin; 239e4fc5df0SBarry Smith if (left) { 240*321429dbSBarry Smith if (!a->left) { 241e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 242e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 243*321429dbSBarry Smith } else { 244*321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 245*321429dbSBarry Smith } 246e4fc5df0SBarry Smith } 247e4fc5df0SBarry Smith if (right) { 248*321429dbSBarry Smith if (!a->right) { 249e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 250e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 251*321429dbSBarry Smith } else { 252*321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 253*321429dbSBarry Smith } 254e4fc5df0SBarry Smith } 255e4fc5df0SBarry Smith PetscFunctionReturn(0); 256e4fc5df0SBarry Smith } 257793850ffSBarry Smith 258793850ffSBarry Smith static struct _MatOps MatOps_Values = {0, 259793850ffSBarry Smith 0, 260793850ffSBarry Smith 0, 261793850ffSBarry Smith MatMult_Composite, 262793850ffSBarry Smith 0, 263793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite, 264793850ffSBarry Smith 0, 265793850ffSBarry Smith 0, 266793850ffSBarry Smith 0, 267793850ffSBarry Smith 0, 268793850ffSBarry Smith /*10*/ 0, 269793850ffSBarry Smith 0, 270793850ffSBarry Smith 0, 271793850ffSBarry Smith 0, 272793850ffSBarry Smith 0, 273793850ffSBarry Smith /*15*/ 0, 274793850ffSBarry Smith 0, 275793850ffSBarry Smith MatGetDiagonal_Composite, 276e4fc5df0SBarry Smith MatDiagonalScale_Composite, 277793850ffSBarry Smith 0, 278793850ffSBarry Smith /*20*/ 0, 279793850ffSBarry Smith MatAssemblyEnd_Composite, 280793850ffSBarry Smith 0, 281793850ffSBarry Smith 0, 282793850ffSBarry Smith 0, 283793850ffSBarry Smith /*25*/ 0, 284793850ffSBarry Smith 0, 285793850ffSBarry Smith 0, 286793850ffSBarry Smith 0, 287793850ffSBarry Smith 0, 288793850ffSBarry Smith /*30*/ 0, 289793850ffSBarry Smith 0, 290793850ffSBarry Smith 0, 291793850ffSBarry Smith 0, 292793850ffSBarry Smith 0, 293793850ffSBarry Smith /*35*/ 0, 294793850ffSBarry Smith 0, 295793850ffSBarry Smith 0, 296793850ffSBarry Smith 0, 297793850ffSBarry Smith 0, 298793850ffSBarry Smith /*40*/ 0, 299793850ffSBarry Smith 0, 300793850ffSBarry Smith 0, 301793850ffSBarry Smith 0, 302793850ffSBarry Smith 0, 303793850ffSBarry Smith /*45*/ 0, 304e4fc5df0SBarry Smith MatScale_Composite, 305793850ffSBarry Smith 0, 306793850ffSBarry Smith 0, 307793850ffSBarry Smith 0, 308793850ffSBarry Smith /*50*/ 0, 309793850ffSBarry Smith 0, 310793850ffSBarry Smith 0, 311793850ffSBarry Smith 0, 312793850ffSBarry Smith 0, 313793850ffSBarry Smith /*55*/ 0, 314793850ffSBarry Smith 0, 315793850ffSBarry Smith 0, 316793850ffSBarry Smith 0, 317793850ffSBarry Smith 0, 318793850ffSBarry Smith /*60*/ 0, 319793850ffSBarry Smith MatDestroy_Composite, 320793850ffSBarry Smith 0, 321793850ffSBarry Smith 0, 322793850ffSBarry Smith 0, 323793850ffSBarry Smith /*65*/ 0, 324793850ffSBarry Smith 0, 325793850ffSBarry Smith 0, 326793850ffSBarry Smith 0, 327793850ffSBarry Smith 0, 328793850ffSBarry Smith /*70*/ 0, 329793850ffSBarry Smith 0, 330793850ffSBarry Smith 0, 331793850ffSBarry Smith 0, 332793850ffSBarry Smith 0, 333793850ffSBarry Smith /*75*/ 0, 334793850ffSBarry Smith 0, 335793850ffSBarry Smith 0, 336793850ffSBarry Smith 0, 337793850ffSBarry Smith 0, 338793850ffSBarry Smith /*80*/ 0, 339793850ffSBarry Smith 0, 340793850ffSBarry Smith 0, 341793850ffSBarry Smith 0, 342793850ffSBarry Smith 0, 343793850ffSBarry Smith /*85*/ 0, 344793850ffSBarry Smith 0, 345793850ffSBarry Smith 0, 346793850ffSBarry Smith 0, 347793850ffSBarry Smith 0, 348793850ffSBarry Smith /*90*/ 0, 349793850ffSBarry Smith 0, 350793850ffSBarry Smith 0, 351793850ffSBarry Smith 0, 352793850ffSBarry Smith 0, 353793850ffSBarry Smith /*95*/ 0, 354793850ffSBarry Smith 0, 355793850ffSBarry Smith 0, 356793850ffSBarry Smith 0}; 357793850ffSBarry Smith 358793850ffSBarry Smith /*MC 3596d7c1e57SBarry Smith MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices (all matrices are of same size and parallel layout). 3606d7c1e57SBarry Smith 3616d7c1e57SBarry Smith Notes: to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 362793850ffSBarry Smith 363793850ffSBarry Smith Level: advanced 364793850ffSBarry Smith 3656d7c1e57SBarry Smith .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 366793850ffSBarry Smith M*/ 367793850ffSBarry Smith 368793850ffSBarry Smith EXTERN_C_BEGIN 369793850ffSBarry Smith #undef __FUNCT__ 370793850ffSBarry Smith #define __FUNCT__ "MatCreate_Composite" 371793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Composite(Mat A) 372793850ffSBarry Smith { 373793850ffSBarry Smith Mat_Composite *b; 374793850ffSBarry Smith PetscErrorCode ierr; 375793850ffSBarry Smith 376793850ffSBarry Smith PetscFunctionBegin; 37738f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Composite,&b);CHKERRQ(ierr); 378793850ffSBarry Smith A->data = (void*)b; 37938f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 380793850ffSBarry Smith 381d0f46423SBarry Smith ierr = PetscMapSetBlockSize(A->rmap,1);CHKERRQ(ierr); 382d0f46423SBarry Smith ierr = PetscMapSetBlockSize(A->cmap,1);CHKERRQ(ierr); 383d0f46423SBarry Smith ierr = PetscMapSetUp(A->rmap);CHKERRQ(ierr); 384d0f46423SBarry Smith ierr = PetscMapSetUp(A->cmap);CHKERRQ(ierr); 385793850ffSBarry Smith 386793850ffSBarry Smith A->assembled = PETSC_TRUE; 387793850ffSBarry Smith A->preallocated = PETSC_FALSE; 3886d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 389e4fc5df0SBarry Smith b->scale = 1.0; 390793850ffSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 391793850ffSBarry Smith PetscFunctionReturn(0); 392793850ffSBarry Smith } 393793850ffSBarry Smith EXTERN_C_END 394793850ffSBarry Smith 395793850ffSBarry Smith #undef __FUNCT__ 396793850ffSBarry Smith #define __FUNCT__ "MatCreateComposite" 397793850ffSBarry Smith /*@C 398793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 399793850ffSBarry Smith 400793850ffSBarry Smith Collective on MPI_Comm 401793850ffSBarry Smith 402793850ffSBarry Smith Input Parameters: 403793850ffSBarry Smith + comm - MPI communicator 404793850ffSBarry Smith . nmat - number of matrices to put in 405793850ffSBarry Smith - mats - the matrices 406793850ffSBarry Smith 407793850ffSBarry Smith Output Parameter: 408793850ffSBarry Smith . A - the matrix 409793850ffSBarry Smith 410793850ffSBarry Smith Level: advanced 411793850ffSBarry Smith 412793850ffSBarry Smith Notes: 413793850ffSBarry Smith Alternative construction 414793850ffSBarry Smith $ MatCreate(comm,&mat); 415793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 416793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 417793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 418793850ffSBarry Smith $ .... 419793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 420b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 421b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 422793850ffSBarry Smith 4236d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 4246d7c1e57SBarry Smith 4256d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType 426793850ffSBarry Smith 427793850ffSBarry Smith @*/ 428793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 429793850ffSBarry Smith { 430793850ffSBarry Smith PetscErrorCode ierr; 431793850ffSBarry Smith PetscInt m,n,M,N,i; 432793850ffSBarry Smith 433793850ffSBarry Smith PetscFunctionBegin; 434793850ffSBarry Smith if (nmat < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 435f3f84630SBarry Smith PetscValidPointer(mat,3); 436793850ffSBarry Smith 437793850ffSBarry Smith ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr); 438793850ffSBarry Smith ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr); 439793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 440793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 441793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 442793850ffSBarry Smith for (i=0; i<nmat; i++) { 443793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 444793850ffSBarry Smith } 445b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 446b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 447793850ffSBarry Smith PetscFunctionReturn(0); 448793850ffSBarry Smith } 449793850ffSBarry Smith 450793850ffSBarry Smith #undef __FUNCT__ 451793850ffSBarry Smith #define __FUNCT__ "MatCompositeAddMat" 452793850ffSBarry Smith /*@ 453793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 454793850ffSBarry Smith 455793850ffSBarry Smith Collective on Mat 456793850ffSBarry Smith 457793850ffSBarry Smith Input Parameters: 458793850ffSBarry Smith + mat - the composite matrix 459793850ffSBarry Smith - smat - the partial matrix 460793850ffSBarry Smith 461793850ffSBarry Smith Level: advanced 462793850ffSBarry Smith 463793850ffSBarry Smith .seealso: MatCreateComposite() 464793850ffSBarry Smith @*/ 465793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat mat,Mat smat) 466793850ffSBarry Smith { 46738f2d2fdSLisandro Dalcin Mat_Composite *shell; 468793850ffSBarry Smith PetscErrorCode ierr; 469793850ffSBarry Smith Mat_CompositeLink ilink,next; 470793850ffSBarry Smith 471793850ffSBarry Smith PetscFunctionBegin; 472793850ffSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 473793850ffSBarry Smith PetscValidHeaderSpecific(smat,MAT_COOKIE,2); 47438f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr); 475793850ffSBarry Smith ilink->next = 0; 476793850ffSBarry Smith ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 477c3122656SLisandro Dalcin ilink->mat = smat; 478793850ffSBarry Smith 47938f2d2fdSLisandro Dalcin shell = (Mat_Composite*)mat->data; 480793850ffSBarry Smith next = shell->head; 481793850ffSBarry Smith if (!next) { 482793850ffSBarry Smith shell->head = ilink; 483793850ffSBarry Smith } else { 484793850ffSBarry Smith while (next->next) { 485793850ffSBarry Smith next = next->next; 486793850ffSBarry Smith } 487793850ffSBarry Smith next->next = ilink; 4886d7c1e57SBarry Smith ilink->prev = next; 4896d7c1e57SBarry Smith } 4906d7c1e57SBarry Smith shell->tail = ilink; 4916d7c1e57SBarry Smith PetscFunctionReturn(0); 4926d7c1e57SBarry Smith } 4936d7c1e57SBarry Smith 4946d7c1e57SBarry Smith #undef __FUNCT__ 4956d7c1e57SBarry Smith #define __FUNCT__ "MatCompositeSetType" 4966d7c1e57SBarry Smith /*@C 4976d7c1e57SBarry Smith MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product 4986d7c1e57SBarry Smith 4996d7c1e57SBarry Smith Collective on MPI_Comm 5006d7c1e57SBarry Smith 5016d7c1e57SBarry Smith Input Parameters: 5026d7c1e57SBarry Smith . mat - the composite matrix 5036d7c1e57SBarry Smith 5046d7c1e57SBarry Smith 5056d7c1e57SBarry Smith Level: advanced 5066d7c1e57SBarry Smith 5076d7c1e57SBarry Smith Notes: 5086d7c1e57SBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 5096d7c1e57SBarry Smith matrix in the composite matrix. 5106d7c1e57SBarry Smith 5116d7c1e57SBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 5126d7c1e57SBarry Smith 5136d7c1e57SBarry Smith @*/ 5146d7c1e57SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat mat,MatCompositeType type) 5156d7c1e57SBarry Smith { 5166d7c1e57SBarry Smith Mat_Composite *b = (Mat_Composite*)mat->data; 5176d7c1e57SBarry Smith PetscTruth flg; 5186d7c1e57SBarry Smith PetscErrorCode ierr; 5196d7c1e57SBarry Smith 5206d7c1e57SBarry Smith PetscFunctionBegin; 5216d7c1e57SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATCOMPOSITE,&flg);CHKERRQ(ierr); 5226d7c1e57SBarry Smith if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Can only use with composite matrix"); 5236d7c1e57SBarry Smith if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 5246d7c1e57SBarry Smith mat->ops->getdiagonal = 0; 5256d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite_Multiplicative; 5266d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 5276d7c1e57SBarry Smith b->type = MAT_COMPOSITE_MULTIPLICATIVE; 5286d7c1e57SBarry Smith } else { 5296d7c1e57SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Composite; 5306d7c1e57SBarry Smith mat->ops->mult = MatMult_Composite; 5316d7c1e57SBarry Smith mat->ops->multtranspose = MatMultTranspose_Composite; 5326d7c1e57SBarry Smith b->type = MAT_COMPOSITE_ADDITIVE; 533793850ffSBarry Smith } 534793850ffSBarry Smith PetscFunctionReturn(0); 535793850ffSBarry Smith } 536793850ffSBarry Smith 5376d7c1e57SBarry Smith 538b52f573bSBarry Smith #undef __FUNCT__ 539b52f573bSBarry Smith #define __FUNCT__ "MatCompositeMerge" 540b52f573bSBarry Smith /*@C 541b52f573bSBarry Smith MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 542b52f573bSBarry Smith by summing all the matrices inside the composite matrix. 543b52f573bSBarry Smith 544b52f573bSBarry Smith Collective on MPI_Comm 545b52f573bSBarry Smith 546b52f573bSBarry Smith Input Parameters: 547b52f573bSBarry Smith . mat - the composite matrix 548b52f573bSBarry Smith 549b52f573bSBarry Smith 550b52f573bSBarry Smith Options Database: 551b52f573bSBarry Smith . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 552b52f573bSBarry Smith 553b52f573bSBarry Smith Level: advanced 554b52f573bSBarry Smith 555b52f573bSBarry Smith Notes: 556b52f573bSBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 557b52f573bSBarry Smith matrix in the composite matrix. 558b52f573bSBarry Smith 559b52f573bSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 560b52f573bSBarry Smith 561b52f573bSBarry Smith @*/ 562b52f573bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat mat) 563b52f573bSBarry Smith { 564b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 5656d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 566b52f573bSBarry Smith PetscErrorCode ierr; 5676d7c1e57SBarry Smith Mat tmat,newmat; 568b52f573bSBarry Smith 569b52f573bSBarry Smith PetscFunctionBegin; 570b52f573bSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 571b52f573bSBarry Smith 572b52f573bSBarry Smith PetscFunctionBegin; 5736d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 574b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 575b52f573bSBarry Smith while ((next = next->next)) { 576b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 577b52f573bSBarry Smith } 5786d7c1e57SBarry Smith } else { 5796d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 5806d7c1e57SBarry Smith while ((prev = prev->prev)) { 5816d7c1e57SBarry Smith ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 5826d7c1e57SBarry Smith ierr = MatDestroy(tmat);CHKERRQ(ierr); 5836d7c1e57SBarry Smith tmat = newmat; 5846d7c1e57SBarry Smith } 5856d7c1e57SBarry Smith } 5866d7c1e57SBarry Smith ierr = MatHeaderReplace(mat,tmat);CHKERRQ(ierr); 587e4fc5df0SBarry Smith ierr = MatDiagonalScale(mat,shell->left,shell->right);CHKERRQ(ierr); 588e4fc5df0SBarry Smith ierr = MatScale(mat,shell->scale);CHKERRQ(ierr); 589b52f573bSBarry Smith PetscFunctionReturn(0); 590b52f573bSBarry Smith } 591