1793850ffSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3793850ffSBarry Smith 48c8613bfSJakub Kruzik const char *const MatCompositeMergeTypes[] = {"left","right","MatCompositeMergeType","MAT_COMPOSITE_",0}; 58c8613bfSJakub Kruzik 6793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink; 7793850ffSBarry Smith struct _Mat_CompositeLink { 8793850ffSBarry Smith Mat mat; 96d7c1e57SBarry Smith Vec work; 106d7c1e57SBarry Smith Mat_CompositeLink next,prev; 11793850ffSBarry Smith }; 12793850ffSBarry Smith 13793850ffSBarry Smith typedef struct { 146d7c1e57SBarry Smith MatCompositeType type; 156d7c1e57SBarry Smith Mat_CompositeLink head,tail; 16793850ffSBarry Smith Vec work; 17e4fc5df0SBarry Smith PetscScalar scale; /* scale factor supplied with MatScale() */ 18e4fc5df0SBarry Smith Vec left,right; /* left and right diagonal scaling provided with MatDiagonalScale() */ 19e4fc5df0SBarry Smith Vec leftwork,rightwork; 206a7ede75SJakub Kruzik PetscInt nmat; 214b2558d6SJakub Kruzik PetscBool merge; 228c8613bfSJakub Kruzik MatCompositeMergeType mergetype; 23793850ffSBarry Smith } Mat_Composite; 24793850ffSBarry Smith 25793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat) 26793850ffSBarry Smith { 27793850ffSBarry Smith PetscErrorCode ierr; 282c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite*)mat->data; 296d7c1e57SBarry Smith Mat_CompositeLink next = shell->head,oldnext; 30793850ffSBarry Smith 31793850ffSBarry Smith PetscFunctionBegin; 32793850ffSBarry Smith while (next) { 336bf464f9SBarry Smith ierr = MatDestroy(&next->mat);CHKERRQ(ierr); 346d7c1e57SBarry Smith if (next->work && (!next->next || next->work != next->next->work)) { 356bf464f9SBarry Smith ierr = VecDestroy(&next->work);CHKERRQ(ierr); 366d7c1e57SBarry Smith } 376d7c1e57SBarry Smith oldnext = next; 38793850ffSBarry Smith next = next->next; 396d7c1e57SBarry Smith ierr = PetscFree(oldnext);CHKERRQ(ierr); 40793850ffSBarry Smith } 416bf464f9SBarry Smith ierr = VecDestroy(&shell->work);CHKERRQ(ierr); 426bf464f9SBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 436bf464f9SBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 446bf464f9SBarry Smith ierr = VecDestroy(&shell->leftwork);CHKERRQ(ierr); 456bf464f9SBarry Smith ierr = VecDestroy(&shell->rightwork);CHKERRQ(ierr); 466bf464f9SBarry Smith ierr = PetscFree(mat->data);CHKERRQ(ierr); 47793850ffSBarry Smith PetscFunctionReturn(0); 48793850ffSBarry Smith } 49793850ffSBarry Smith 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 */ 692a7a6963SBarry Smith ierr = MatCreateVecs(next->mat,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 PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y) 856d7c1e57SBarry Smith { 866d7c1e57SBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 876d7c1e57SBarry Smith Mat_CompositeLink tail = shell->tail; 886d7c1e57SBarry Smith PetscErrorCode ierr; 896d7c1e57SBarry Smith Vec in,out; 906d7c1e57SBarry Smith 916d7c1e57SBarry Smith PetscFunctionBegin; 92e32f2f54SBarry Smith if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 936d7c1e57SBarry Smith in = x; 94e4fc5df0SBarry Smith if (shell->left) { 95e4fc5df0SBarry Smith if (!shell->leftwork) { 96e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 97e4fc5df0SBarry Smith } 98e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 99e4fc5df0SBarry Smith in = shell->leftwork; 100e4fc5df0SBarry Smith } 1016d7c1e57SBarry Smith while (tail->prev) { 1026d7c1e57SBarry Smith if (!tail->prev->work) { /* should reuse previous work if the same size */ 1032a7a6963SBarry Smith ierr = MatCreateVecs(tail->mat,NULL,&tail->prev->work);CHKERRQ(ierr); 1046d7c1e57SBarry Smith } 1056d7c1e57SBarry Smith out = tail->prev->work; 1066d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr); 1076d7c1e57SBarry Smith in = out; 1086d7c1e57SBarry Smith tail = tail->prev; 1096d7c1e57SBarry Smith } 1106d7c1e57SBarry Smith ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr); 111e4fc5df0SBarry Smith if (shell->right) { 112e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 113e4fc5df0SBarry Smith } 114e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 1156d7c1e57SBarry Smith PetscFunctionReturn(0); 1166d7c1e57SBarry Smith } 1176d7c1e57SBarry Smith 118793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 119793850ffSBarry Smith { 120793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 121793850ffSBarry Smith Mat_CompositeLink next = shell->head; 122793850ffSBarry Smith PetscErrorCode ierr; 123e4fc5df0SBarry Smith Vec in; 124793850ffSBarry Smith 125793850ffSBarry Smith PetscFunctionBegin; 126e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 127e4fc5df0SBarry Smith in = x; 128e4fc5df0SBarry Smith if (shell->right) { 129e4fc5df0SBarry Smith if (!shell->rightwork) { 130e4fc5df0SBarry Smith ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr); 131793850ffSBarry Smith } 132e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr); 133e4fc5df0SBarry Smith in = shell->rightwork; 134e4fc5df0SBarry Smith } 135e4fc5df0SBarry Smith ierr = MatMult(next->mat,in,y);CHKERRQ(ierr); 136e4fc5df0SBarry Smith while ((next = next->next)) { 137e4fc5df0SBarry Smith ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr); 138e4fc5df0SBarry Smith } 139e4fc5df0SBarry Smith if (shell->left) { 140e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr); 141e4fc5df0SBarry Smith } 142e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 143793850ffSBarry Smith PetscFunctionReturn(0); 144793850ffSBarry Smith } 145793850ffSBarry Smith 146793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 147793850ffSBarry Smith { 148793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 149793850ffSBarry Smith Mat_CompositeLink next = shell->head; 150793850ffSBarry Smith PetscErrorCode ierr; 151e4fc5df0SBarry Smith Vec in; 152793850ffSBarry Smith 153793850ffSBarry Smith PetscFunctionBegin; 154e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 155e4fc5df0SBarry Smith in = x; 156e4fc5df0SBarry Smith if (shell->left) { 157e4fc5df0SBarry Smith if (!shell->leftwork) { 158e4fc5df0SBarry Smith ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr); 159793850ffSBarry Smith } 160e4fc5df0SBarry Smith ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr); 161e4fc5df0SBarry Smith in = shell->leftwork; 162e4fc5df0SBarry Smith } 163e4fc5df0SBarry Smith ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr); 164e4fc5df0SBarry Smith while ((next = next->next)) { 165e4fc5df0SBarry Smith ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr); 166e4fc5df0SBarry Smith } 167e4fc5df0SBarry Smith if (shell->right) { 168e4fc5df0SBarry Smith ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr); 169e4fc5df0SBarry Smith } 170e4fc5df0SBarry Smith ierr = VecScale(y,shell->scale);CHKERRQ(ierr); 171793850ffSBarry Smith PetscFunctionReturn(0); 172793850ffSBarry Smith } 173793850ffSBarry Smith 1747bf3a71bSJakub Kruzik PetscErrorCode MatMultAdd_Composite(Mat A,Vec x,Vec y,Vec z) 1757bf3a71bSJakub Kruzik { 1767bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 1777bf3a71bSJakub Kruzik PetscErrorCode ierr; 1787bf3a71bSJakub Kruzik 1797bf3a71bSJakub Kruzik PetscFunctionBegin; 1807bf3a71bSJakub Kruzik if (y != z) { 1817bf3a71bSJakub Kruzik ierr = MatMult(A,x,z);CHKERRQ(ierr); 1827bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 1837bf3a71bSJakub Kruzik } else { 1847bf3a71bSJakub Kruzik if (!shell->leftwork) { 1857bf3a71bSJakub Kruzik ierr = VecDuplicate(z,&shell->leftwork);CHKERRQ(ierr); 1867bf3a71bSJakub Kruzik } 1877bf3a71bSJakub Kruzik ierr = MatMult(A,x,shell->leftwork);CHKERRQ(ierr); 1887bf3a71bSJakub Kruzik ierr = VecCopy(y,z);CHKERRQ(ierr); 1897bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,shell->leftwork);CHKERRQ(ierr); 1907bf3a71bSJakub Kruzik } 1917bf3a71bSJakub Kruzik PetscFunctionReturn(0); 1927bf3a71bSJakub Kruzik } 1937bf3a71bSJakub Kruzik 1947bf3a71bSJakub Kruzik PetscErrorCode MatMultTransposeAdd_Composite(Mat A,Vec x,Vec y, Vec z) 1957bf3a71bSJakub Kruzik { 1967bf3a71bSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)A->data; 1977bf3a71bSJakub Kruzik PetscErrorCode ierr; 1987bf3a71bSJakub Kruzik 1997bf3a71bSJakub Kruzik PetscFunctionBegin; 2007bf3a71bSJakub Kruzik if (y != z) { 2017bf3a71bSJakub Kruzik ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 2027bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2037bf3a71bSJakub Kruzik } else { 2047bf3a71bSJakub Kruzik if (!shell->rightwork) { 2057bf3a71bSJakub Kruzik ierr = VecDuplicate(z,&shell->rightwork);CHKERRQ(ierr); 2067bf3a71bSJakub Kruzik } 2077bf3a71bSJakub Kruzik ierr = MatMultTranspose(A,x,shell->rightwork);CHKERRQ(ierr); 2087bf3a71bSJakub Kruzik ierr = VecCopy(y,z);CHKERRQ(ierr); 2097bf3a71bSJakub Kruzik ierr = VecAXPY(z,1.0,shell->rightwork);CHKERRQ(ierr); 2107bf3a71bSJakub Kruzik } 2117bf3a71bSJakub Kruzik PetscFunctionReturn(0); 2127bf3a71bSJakub Kruzik } 2137bf3a71bSJakub Kruzik 214793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 215793850ffSBarry Smith { 216793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 217793850ffSBarry Smith Mat_CompositeLink next = shell->head; 218793850ffSBarry Smith PetscErrorCode ierr; 219793850ffSBarry Smith 220793850ffSBarry Smith PetscFunctionBegin; 221e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 222e32f2f54SBarry Smith if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling"); 223e4fc5df0SBarry Smith 224793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 225793850ffSBarry Smith if (next->next && !shell->work) { 226793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 227793850ffSBarry Smith } 228793850ffSBarry Smith while ((next = next->next)) { 229793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 230793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 231793850ffSBarry Smith } 232e4fc5df0SBarry Smith ierr = VecScale(v,shell->scale);CHKERRQ(ierr); 233793850ffSBarry Smith PetscFunctionReturn(0); 234793850ffSBarry Smith } 235793850ffSBarry Smith 236793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 237793850ffSBarry Smith { 2384b2558d6SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)Y->data; 239b52f573bSBarry Smith PetscErrorCode ierr; 240b52f573bSBarry Smith 241793850ffSBarry Smith PetscFunctionBegin; 2424b2558d6SJakub Kruzik if (shell->merge) { 243b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 244b52f573bSBarry Smith } 245793850ffSBarry Smith PetscFunctionReturn(0); 246793850ffSBarry Smith } 247793850ffSBarry Smith 248e4fc5df0SBarry Smith PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha) 249e4fc5df0SBarry Smith { 250e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 251e4fc5df0SBarry Smith 252e4fc5df0SBarry Smith PetscFunctionBegin; 253321429dbSBarry Smith a->scale *= alpha; 254e4fc5df0SBarry Smith PetscFunctionReturn(0); 255e4fc5df0SBarry Smith } 256e4fc5df0SBarry Smith 257e4fc5df0SBarry Smith PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right) 258e4fc5df0SBarry Smith { 259e4fc5df0SBarry Smith Mat_Composite *a = (Mat_Composite*)inA->data; 260e4fc5df0SBarry Smith PetscErrorCode ierr; 261e4fc5df0SBarry Smith 262e4fc5df0SBarry Smith PetscFunctionBegin; 263e4fc5df0SBarry Smith if (left) { 264321429dbSBarry Smith if (!a->left) { 265e4fc5df0SBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 266e4fc5df0SBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 267321429dbSBarry Smith } else { 268321429dbSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 269321429dbSBarry Smith } 270e4fc5df0SBarry Smith } 271e4fc5df0SBarry Smith if (right) { 272321429dbSBarry Smith if (!a->right) { 273e4fc5df0SBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 274e4fc5df0SBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 275321429dbSBarry Smith } else { 276321429dbSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 277321429dbSBarry Smith } 278e4fc5df0SBarry Smith } 279e4fc5df0SBarry Smith PetscFunctionReturn(0); 280e4fc5df0SBarry Smith } 281793850ffSBarry Smith 2824b2558d6SJakub Kruzik PetscErrorCode MatSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,Mat A) 2834b2558d6SJakub Kruzik { 2844b2558d6SJakub Kruzik Mat_Composite *a = (Mat_Composite*)A->data; 2854b2558d6SJakub Kruzik PetscErrorCode ierr; 2864b2558d6SJakub Kruzik 2874b2558d6SJakub Kruzik PetscFunctionBegin; 2884b2558d6SJakub Kruzik ierr = PetscOptionsHead(PetscOptionsObject,"MATCOMPOSITE options");CHKERRQ(ierr); 2894b2558d6SJakub Kruzik ierr = PetscOptionsBool("-mat_composite_merge","Merge at MatAssemblyEnd","MatCompositeMerge",a->merge,&a->merge,NULL);CHKERRQ(ierr); 2908c8613bfSJakub Kruzik ierr = PetscOptionsEnum("-mat_composite_merge_type","Set composite merge direction","MatCompositeSetMergeType",MatCompositeMergeTypes,(PetscEnum)a->mergetype,(PetscEnum*)&a->mergetype,NULL);CHKERRQ(ierr); 2914b2558d6SJakub Kruzik ierr = PetscOptionsTail();CHKERRQ(ierr); 2924b2558d6SJakub Kruzik PetscFunctionReturn(0); 2934b2558d6SJakub Kruzik } 2944b2558d6SJakub Kruzik 2952c0821f3SBarry Smith /*@ 2968bd739bdSJakub Kruzik MatCreateComposite - Creates a matrix as the sum or product of one or more matrices 297793850ffSBarry Smith 298793850ffSBarry Smith Collective on MPI_Comm 299793850ffSBarry Smith 300793850ffSBarry Smith Input Parameters: 301793850ffSBarry Smith + comm - MPI communicator 302793850ffSBarry Smith . nmat - number of matrices to put in 303793850ffSBarry Smith - mats - the matrices 304793850ffSBarry Smith 305793850ffSBarry Smith Output Parameter: 306db36af27SMatthew G. Knepley . mat - the matrix 307793850ffSBarry Smith 3084b2558d6SJakub Kruzik Options Database Keys: 3094b2558d6SJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 310b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 3114b2558d6SJakub Kruzik 312793850ffSBarry Smith Level: advanced 313793850ffSBarry Smith 314793850ffSBarry Smith Notes: 315793850ffSBarry Smith Alternative construction 316793850ffSBarry Smith $ MatCreate(comm,&mat); 317793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 318793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 319793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 320793850ffSBarry Smith $ .... 321793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 322b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 323b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 324793850ffSBarry Smith 3256d7c1e57SBarry Smith For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0] 3266d7c1e57SBarry Smith 3276dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeGetMat(), MatCompositeMerge(), MatCompositeSetType(), MATCOMPOSITE 328793850ffSBarry Smith 329793850ffSBarry Smith @*/ 3307087cfbeSBarry Smith PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 331793850ffSBarry Smith { 332793850ffSBarry Smith PetscErrorCode ierr; 333793850ffSBarry Smith PetscInt m,n,M,N,i; 334793850ffSBarry Smith 335793850ffSBarry Smith PetscFunctionBegin; 336e32f2f54SBarry Smith if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 337f3f84630SBarry Smith PetscValidPointer(mat,3); 338793850ffSBarry Smith 339c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[0],PETSC_IGNORE,&n);CHKERRQ(ierr); 340c4f13f43SJakub Kruzik ierr = MatGetLocalSize(mats[nmat-1],&m,PETSC_IGNORE);CHKERRQ(ierr); 341c4f13f43SJakub Kruzik ierr = MatGetSize(mats[0],PETSC_IGNORE,&N);CHKERRQ(ierr); 342c4f13f43SJakub Kruzik ierr = MatGetSize(mats[nmat-1],&M,PETSC_IGNORE);CHKERRQ(ierr); 343793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 344793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 345793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 346793850ffSBarry Smith for (i=0; i<nmat; i++) { 347793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 348793850ffSBarry Smith } 349b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 350b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351793850ffSBarry Smith PetscFunctionReturn(0); 352793850ffSBarry Smith } 353793850ffSBarry Smith 354d7f81bf2SJakub Kruzik 355d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeAddMat_Composite(Mat mat,Mat smat) 356d7f81bf2SJakub Kruzik { 357d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 358d7f81bf2SJakub Kruzik Mat_CompositeLink ilink,next = shell->head; 359d7f81bf2SJakub Kruzik PetscErrorCode ierr; 360d7f81bf2SJakub Kruzik 361d7f81bf2SJakub Kruzik PetscFunctionBegin; 362d7f81bf2SJakub Kruzik ierr = PetscNewLog(mat,&ilink);CHKERRQ(ierr); 363d7f81bf2SJakub Kruzik ilink->next = 0; 364d7f81bf2SJakub Kruzik ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 365d7f81bf2SJakub Kruzik ilink->mat = smat; 366d7f81bf2SJakub Kruzik 367d7f81bf2SJakub Kruzik if (!next) shell->head = ilink; 368d7f81bf2SJakub Kruzik else { 369d7f81bf2SJakub Kruzik while (next->next) { 370d7f81bf2SJakub Kruzik next = next->next; 371d7f81bf2SJakub Kruzik } 372d7f81bf2SJakub Kruzik next->next = ilink; 373d7f81bf2SJakub Kruzik ilink->prev = next; 374d7f81bf2SJakub Kruzik } 375d7f81bf2SJakub Kruzik shell->tail = ilink; 376d7f81bf2SJakub Kruzik shell->nmat += 1; 377d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 378d7f81bf2SJakub Kruzik } 379d7f81bf2SJakub Kruzik 380793850ffSBarry Smith /*@ 3818bd739bdSJakub Kruzik MatCompositeAddMat - Add another matrix to a composite matrix. 382793850ffSBarry Smith 383793850ffSBarry Smith Collective on Mat 384793850ffSBarry Smith 385793850ffSBarry Smith Input Parameters: 386793850ffSBarry Smith + mat - the composite matrix 387793850ffSBarry Smith - smat - the partial matrix 388793850ffSBarry Smith 389793850ffSBarry Smith Level: advanced 390793850ffSBarry Smith 3918bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE 392793850ffSBarry Smith @*/ 3937087cfbeSBarry Smith PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat) 394793850ffSBarry Smith { 395793850ffSBarry Smith PetscErrorCode ierr; 396793850ffSBarry Smith 397793850ffSBarry Smith PetscFunctionBegin; 3980700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 3990700a824SBarry Smith PetscValidHeaderSpecific(smat,MAT_CLASSID,2); 400d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeAddMat_C",(Mat,Mat),(mat,smat));CHKERRQ(ierr); 401d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 402d7f81bf2SJakub Kruzik } 403793850ffSBarry Smith 404d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeSetType_Composite(Mat mat,MatCompositeType type) 405d7f81bf2SJakub Kruzik { 406d7f81bf2SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 407d7f81bf2SJakub Kruzik 408d7f81bf2SJakub Kruzik PetscFunctionBegin; 409ced1ac25SJakub Kruzik b->type = type; 410d7f81bf2SJakub Kruzik if (type == MAT_COMPOSITE_MULTIPLICATIVE) { 411d7f81bf2SJakub Kruzik mat->ops->getdiagonal = 0; 412d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite_Multiplicative; 413d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative; 414d7f81bf2SJakub Kruzik } else { 415d7f81bf2SJakub Kruzik mat->ops->getdiagonal = MatGetDiagonal_Composite; 416d7f81bf2SJakub Kruzik mat->ops->mult = MatMult_Composite; 417d7f81bf2SJakub Kruzik mat->ops->multtranspose = MatMultTranspose_Composite; 418793850ffSBarry Smith } 4196d7c1e57SBarry Smith PetscFunctionReturn(0); 4206d7c1e57SBarry Smith } 4216d7c1e57SBarry Smith 4222c0821f3SBarry Smith /*@ 4238bd739bdSJakub Kruzik MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product. 4246d7c1e57SBarry Smith 425*b2b245abSJakub Kruzik Logically Collective on Mat 4266d7c1e57SBarry Smith 4276d7c1e57SBarry Smith Input Parameters: 4286d7c1e57SBarry Smith . mat - the composite matrix 4296d7c1e57SBarry Smith 4306d7c1e57SBarry Smith Level: advanced 4316d7c1e57SBarry Smith 4326dbc55e5SJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeGetType(), MATCOMPOSITE 4336d7c1e57SBarry Smith 4346d7c1e57SBarry Smith @*/ 4357087cfbeSBarry Smith PetscErrorCode MatCompositeSetType(Mat mat,MatCompositeType type) 4366d7c1e57SBarry Smith { 4376d7c1e57SBarry Smith PetscErrorCode ierr; 4386d7c1e57SBarry Smith 4396d7c1e57SBarry Smith PetscFunctionBegin; 440d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 441*b2b245abSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 442d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetType_C",(Mat,MatCompositeType),(mat,type));CHKERRQ(ierr); 443793850ffSBarry Smith PetscFunctionReturn(0); 444793850ffSBarry Smith } 445793850ffSBarry Smith 4466dbc55e5SJakub Kruzik static PetscErrorCode MatCompositeGetType_Composite(Mat mat,MatCompositeType *type) 4476dbc55e5SJakub Kruzik { 4486dbc55e5SJakub Kruzik Mat_Composite *b = (Mat_Composite*)mat->data; 4496dbc55e5SJakub Kruzik 4506dbc55e5SJakub Kruzik PetscFunctionBegin; 4516dbc55e5SJakub Kruzik *type = b->type; 4526dbc55e5SJakub Kruzik PetscFunctionReturn(0); 4536dbc55e5SJakub Kruzik } 4546dbc55e5SJakub Kruzik 4556dbc55e5SJakub Kruzik /*@ 4566dbc55e5SJakub Kruzik MatCompositeGetType - Returns type of composite. 4576dbc55e5SJakub Kruzik 4586dbc55e5SJakub Kruzik Not Collective 4596dbc55e5SJakub Kruzik 4606dbc55e5SJakub Kruzik Input Parameter: 4616dbc55e5SJakub Kruzik . mat - the composite matrix 4626dbc55e5SJakub Kruzik 4636dbc55e5SJakub Kruzik Output Parameter: 4646dbc55e5SJakub Kruzik . type - type of composite 4656dbc55e5SJakub Kruzik 4666dbc55e5SJakub Kruzik Level: advanced 4676dbc55e5SJakub Kruzik 4686dbc55e5SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeSetType(), MATCOMPOSITE 4696dbc55e5SJakub Kruzik 4706dbc55e5SJakub Kruzik @*/ 4716dbc55e5SJakub Kruzik PetscErrorCode MatCompositeGetType(Mat mat,MatCompositeType *type) 4726dbc55e5SJakub Kruzik { 4736dbc55e5SJakub Kruzik PetscErrorCode ierr; 4746dbc55e5SJakub Kruzik 4756dbc55e5SJakub Kruzik PetscFunctionBegin; 4766dbc55e5SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4776dbc55e5SJakub Kruzik PetscValidPointer(type,2); 4786dbc55e5SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetType_C",(Mat,MatCompositeType*),(mat,type));CHKERRQ(ierr); 4796dbc55e5SJakub Kruzik PetscFunctionReturn(0); 4806dbc55e5SJakub Kruzik } 4816dbc55e5SJakub Kruzik 4828c8613bfSJakub Kruzik static PetscErrorCode MatCompositeSetMergeType_Composite(Mat mat,MatCompositeMergeType type) 4838c8613bfSJakub Kruzik { 4848c8613bfSJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 4858c8613bfSJakub Kruzik 4868c8613bfSJakub Kruzik PetscFunctionBegin; 4878c8613bfSJakub Kruzik shell->mergetype = type; 4888c8613bfSJakub Kruzik PetscFunctionReturn(0); 4898c8613bfSJakub Kruzik } 4908c8613bfSJakub Kruzik 4918c8613bfSJakub Kruzik /*@ 4928c8613bfSJakub Kruzik MatCompositeSetMergeType - Sets order of MatCompositeMerge(). 4938c8613bfSJakub Kruzik 4948c8613bfSJakub Kruzik Logically Collective on Mat 4958c8613bfSJakub Kruzik 4968c8613bfSJakub Kruzik Input Parameter: 4978c8613bfSJakub Kruzik + mat - the composite matrix 4988c8613bfSJakub Kruzik - type - MAT_COMPOSITE_MERGE RIGHT (default) to start merge from right with the first added matrix (mat[0]), 4998c8613bfSJakub Kruzik MAT_COMPOSITE_MERGE_LEFT to start merge from left with the last added matrix (mat[nmat-1]) 5008c8613bfSJakub Kruzik 5018c8613bfSJakub Kruzik Level: advanced 5028c8613bfSJakub Kruzik 5038c8613bfSJakub Kruzik Notes: 5048c8613bfSJakub Kruzik The resulting matrix is the same regardles of the MergeType. Only the order of operation is changed. 5058c8613bfSJakub Kruzik If set to MAT_COMPOSITE_MERGE_RIGHT the order of the merge is mat[nmat-1]*(mat[nmat-2]*(...*(mat[1]*mat[0]))) 5068c8613bfSJakub Kruzik otherwise the order is (((mat[nmat-1]*mat[nmat-2])*mat[nmat-3])*...)*mat[0]. 5078c8613bfSJakub Kruzik 5088c8613bfSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeMerge(), MATCOMPOSITE 5098c8613bfSJakub Kruzik 5108c8613bfSJakub Kruzik @*/ 5118c8613bfSJakub Kruzik PetscErrorCode MatCompositeSetMergeType(Mat mat,MatCompositeMergeType type) 5128c8613bfSJakub Kruzik { 5138c8613bfSJakub Kruzik PetscErrorCode ierr; 5148c8613bfSJakub Kruzik 5158c8613bfSJakub Kruzik PetscFunctionBegin; 5168c8613bfSJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5178c8613bfSJakub Kruzik PetscValidLogicalCollectiveEnum(mat,type,2); 5188c8613bfSJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeSetMergeType_C",(Mat,MatCompositeMergeType),(mat,type));CHKERRQ(ierr); 5198c8613bfSJakub Kruzik PetscFunctionReturn(0); 5208c8613bfSJakub Kruzik } 5218c8613bfSJakub Kruzik 522d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeMerge_Composite(Mat mat) 523b52f573bSBarry Smith { 524b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 5256d7c1e57SBarry Smith Mat_CompositeLink next = shell->head, prev = shell->tail; 526b52f573bSBarry Smith PetscErrorCode ierr; 5276d7c1e57SBarry Smith Mat tmat,newmat; 5281ba60692SJed Brown Vec left,right; 5291ba60692SJed Brown PetscScalar scale; 530b52f573bSBarry Smith 531b52f573bSBarry Smith PetscFunctionBegin; 532e32f2f54SBarry Smith if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 533b52f573bSBarry Smith 534b52f573bSBarry Smith PetscFunctionBegin; 5356d7c1e57SBarry Smith if (shell->type == MAT_COMPOSITE_ADDITIVE) { 5368c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 537b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 538b52f573bSBarry Smith while ((next = next->next)) { 539b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 540b52f573bSBarry Smith } 5416d7c1e57SBarry Smith } else { 5428c8613bfSJakub Kruzik ierr = MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 5438c8613bfSJakub Kruzik while ((prev = prev->prev)) { 5448c8613bfSJakub Kruzik ierr = MatAXPY(tmat,1.0,prev->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 5458c8613bfSJakub Kruzik } 5468c8613bfSJakub Kruzik } 5478c8613bfSJakub Kruzik } else { 5488c8613bfSJakub Kruzik if (shell->mergetype == MAT_COMPOSITE_MERGE_RIGHT) { 5496d7c1e57SBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 550b6cfcaf5SJakub Kruzik while ((next = next->next)) { 551b6cfcaf5SJakub Kruzik ierr = MatMatMult(next->mat,tmat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 5526bf464f9SBarry Smith ierr = MatDestroy(&tmat);CHKERRQ(ierr); 5536d7c1e57SBarry Smith tmat = newmat; 5546d7c1e57SBarry Smith } 55504d534ceSJakub Kruzik } else { 55604d534ceSJakub Kruzik ierr = MatDuplicate(prev->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 55704d534ceSJakub Kruzik while ((prev = prev->prev)) { 55804d534ceSJakub Kruzik ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr); 55904d534ceSJakub Kruzik ierr = MatDestroy(&tmat);CHKERRQ(ierr); 56004d534ceSJakub Kruzik tmat = newmat; 56104d534ceSJakub Kruzik } 56204d534ceSJakub Kruzik } 5636d7c1e57SBarry Smith } 5641ba60692SJed Brown 5651ba60692SJed Brown scale = shell->scale; 5661ba60692SJed Brown if ((left = shell->left)) {ierr = PetscObjectReference((PetscObject)left);CHKERRQ(ierr);} 5671ba60692SJed Brown if ((right = shell->right)) {ierr = PetscObjectReference((PetscObject)right);CHKERRQ(ierr);} 5681ba60692SJed Brown 56928be2f97SBarry Smith ierr = MatHeaderReplace(mat,&tmat);CHKERRQ(ierr); 5701ba60692SJed Brown 5711ba60692SJed Brown ierr = MatDiagonalScale(mat,left,right);CHKERRQ(ierr); 5721ba60692SJed Brown ierr = MatScale(mat,scale);CHKERRQ(ierr); 5731ba60692SJed Brown ierr = VecDestroy(&left);CHKERRQ(ierr); 5741ba60692SJed Brown ierr = VecDestroy(&right);CHKERRQ(ierr); 575b52f573bSBarry Smith PetscFunctionReturn(0); 576b52f573bSBarry Smith } 5776a7ede75SJakub Kruzik 5786a7ede75SJakub Kruzik /*@ 579d7f81bf2SJakub Kruzik MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 5808bd739bdSJakub Kruzik by summing or computing the product of all the matrices inside the composite matrix. 581d7f81bf2SJakub Kruzik 582*b2b245abSJakub Kruzik Collective 583d7f81bf2SJakub Kruzik 584d7f81bf2SJakub Kruzik Input Parameters: 585d7f81bf2SJakub Kruzik . mat - the composite matrix 586d7f81bf2SJakub Kruzik 587d7f81bf2SJakub Kruzik 5884b2558d6SJakub Kruzik Options Database Keys: 589b28d0daaSJakub Kruzik + -mat_composite_merge - merge in MatAssemblyEnd() 590b28d0daaSJakub Kruzik - -mat_composite_merge_type - set merge direction 591d7f81bf2SJakub Kruzik 592d7f81bf2SJakub Kruzik Level: advanced 593d7f81bf2SJakub Kruzik 594d7f81bf2SJakub Kruzik Notes: 595d7f81bf2SJakub Kruzik The MatType of the resulting matrix will be the same as the MatType of the FIRST 596d7f81bf2SJakub Kruzik matrix in the composite matrix. 597d7f81bf2SJakub Kruzik 5988c8613bfSJakub Kruzik .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MatCompositeSetMergeType(), MATCOMPOSITE 599d7f81bf2SJakub Kruzik 600d7f81bf2SJakub Kruzik @*/ 601d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeMerge(Mat mat) 602d7f81bf2SJakub Kruzik { 603d7f81bf2SJakub Kruzik PetscErrorCode ierr; 604d7f81bf2SJakub Kruzik 605d7f81bf2SJakub Kruzik PetscFunctionBegin; 606d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 607d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeMerge_C",(Mat),(mat));CHKERRQ(ierr); 608d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 609d7f81bf2SJakub Kruzik } 610d7f81bf2SJakub Kruzik 6116d0add67SJakub Kruzik static PetscErrorCode MatCompositeGetNumberMat_Composite(Mat mat,PetscInt *nmat) 612d7f81bf2SJakub Kruzik { 613d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 614d7f81bf2SJakub Kruzik 615d7f81bf2SJakub Kruzik PetscFunctionBegin; 616d7f81bf2SJakub Kruzik *nmat = shell->nmat; 617d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 618d7f81bf2SJakub Kruzik } 619d7f81bf2SJakub Kruzik 620d7f81bf2SJakub Kruzik /*@ 6216d0add67SJakub Kruzik MatCompositeGetNumberMat - Returns the number of matrices in the composite matrix. 6226a7ede75SJakub Kruzik 6236a7ede75SJakub Kruzik Not Collective 6246a7ede75SJakub Kruzik 6256a7ede75SJakub Kruzik Input Parameter: 626d7f81bf2SJakub Kruzik . mat - the composite matrix 6276a7ede75SJakub Kruzik 6286a7ede75SJakub Kruzik Output Parameter: 6296d0add67SJakub Kruzik . nmat - number of matrices in the composite matrix 6306a7ede75SJakub Kruzik 6318b5c3584SJakub Kruzik Level: advanced 6326a7ede75SJakub Kruzik 6338bd739bdSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetMat(), MATCOMPOSITE 6346a7ede75SJakub Kruzik 6356a7ede75SJakub Kruzik @*/ 6366d0add67SJakub Kruzik PetscErrorCode MatCompositeGetNumberMat(Mat mat,PetscInt *nmat) 6376a7ede75SJakub Kruzik { 638d7f81bf2SJakub Kruzik PetscErrorCode ierr; 6396a7ede75SJakub Kruzik 6406a7ede75SJakub Kruzik PetscFunctionBegin; 641d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6426a7ede75SJakub Kruzik PetscValidPointer(nmat,2); 6436d0add67SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetNumberMat_C",(Mat,PetscInt*),(mat,nmat));CHKERRQ(ierr); 644d7f81bf2SJakub Kruzik PetscFunctionReturn(0); 645d7f81bf2SJakub Kruzik } 646d7f81bf2SJakub Kruzik 647d7f81bf2SJakub Kruzik static PetscErrorCode MatCompositeGetMat_Composite(Mat mat,PetscInt i,Mat *Ai) 648d7f81bf2SJakub Kruzik { 649d7f81bf2SJakub Kruzik Mat_Composite *shell = (Mat_Composite*)mat->data; 650d7f81bf2SJakub Kruzik Mat_CompositeLink ilink; 651d7f81bf2SJakub Kruzik PetscInt k; 652d7f81bf2SJakub Kruzik 653d7f81bf2SJakub Kruzik PetscFunctionBegin; 654d7f81bf2SJakub Kruzik if (i >= shell->nmat) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"index out of range: %d >= %d",i,shell->nmat); 655d7f81bf2SJakub Kruzik ilink = shell->head; 656d7f81bf2SJakub Kruzik for (k=0; k<i; k++) { 657d7f81bf2SJakub Kruzik ilink = ilink->next; 658d7f81bf2SJakub Kruzik } 659d7f81bf2SJakub Kruzik *Ai = ilink->mat; 6606a7ede75SJakub Kruzik PetscFunctionReturn(0); 6616a7ede75SJakub Kruzik } 6626a7ede75SJakub Kruzik 6638b5c3584SJakub Kruzik /*@ 6648bd739bdSJakub Kruzik MatCompositeGetMat - Returns the ith matrix from the composite matrix. 6658b5c3584SJakub Kruzik 6668b5c3584SJakub Kruzik Logically Collective on Mat 6678b5c3584SJakub Kruzik 6688b5c3584SJakub Kruzik Input Parameter: 669d7f81bf2SJakub Kruzik + mat - the composite matrix 6708b5c3584SJakub Kruzik - i - the number of requested matrix 6718b5c3584SJakub Kruzik 6728b5c3584SJakub Kruzik Output Parameter: 6738b5c3584SJakub Kruzik . Ai - ith matrix in composite 6748b5c3584SJakub Kruzik 6758b5c3584SJakub Kruzik Level: advanced 6768b5c3584SJakub Kruzik 6776d0add67SJakub Kruzik .seealso: MatCreateComposite(), MatCompositeGetNumberMat(), MatCompositeAddMat(), MATCOMPOSITE 6788b5c3584SJakub Kruzik 6798b5c3584SJakub Kruzik @*/ 680d7f81bf2SJakub Kruzik PetscErrorCode MatCompositeGetMat(Mat mat,PetscInt i,Mat *Ai) 6818b5c3584SJakub Kruzik { 682d7f81bf2SJakub Kruzik PetscErrorCode ierr; 6838b5c3584SJakub Kruzik 6848b5c3584SJakub Kruzik PetscFunctionBegin; 685d7f81bf2SJakub Kruzik PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 686d7f81bf2SJakub Kruzik PetscValidLogicalCollectiveInt(mat,i,2); 6878b5c3584SJakub Kruzik PetscValidPointer(Ai,3); 688d7f81bf2SJakub Kruzik ierr = PetscUseMethod(mat,"MatCompositeGetMat_C",(Mat,PetscInt,Mat*),(mat,i,Ai));CHKERRQ(ierr); 6898b5c3584SJakub Kruzik PetscFunctionReturn(0); 6908b5c3584SJakub Kruzik } 6918b5c3584SJakub Kruzik 69241cd0125SJakub Kruzik static struct _MatOps MatOps_Values = {0, 69341cd0125SJakub Kruzik 0, 69441cd0125SJakub Kruzik 0, 69541cd0125SJakub Kruzik MatMult_Composite, 6967bf3a71bSJakub Kruzik MatMultAdd_Composite, 69741cd0125SJakub Kruzik /* 5*/ MatMultTranspose_Composite, 6987bf3a71bSJakub Kruzik MatMultTransposeAdd_Composite, 69941cd0125SJakub Kruzik 0, 70041cd0125SJakub Kruzik 0, 70141cd0125SJakub Kruzik 0, 70241cd0125SJakub Kruzik /* 10*/ 0, 70341cd0125SJakub Kruzik 0, 70441cd0125SJakub Kruzik 0, 70541cd0125SJakub Kruzik 0, 70641cd0125SJakub Kruzik 0, 70741cd0125SJakub Kruzik /* 15*/ 0, 70841cd0125SJakub Kruzik 0, 70941cd0125SJakub Kruzik MatGetDiagonal_Composite, 71041cd0125SJakub Kruzik MatDiagonalScale_Composite, 71141cd0125SJakub Kruzik 0, 71241cd0125SJakub Kruzik /* 20*/ 0, 71341cd0125SJakub Kruzik MatAssemblyEnd_Composite, 71441cd0125SJakub Kruzik 0, 71541cd0125SJakub Kruzik 0, 71641cd0125SJakub Kruzik /* 24*/ 0, 71741cd0125SJakub Kruzik 0, 71841cd0125SJakub Kruzik 0, 71941cd0125SJakub Kruzik 0, 72041cd0125SJakub Kruzik 0, 72141cd0125SJakub Kruzik /* 29*/ 0, 72241cd0125SJakub Kruzik 0, 72341cd0125SJakub Kruzik 0, 72441cd0125SJakub Kruzik 0, 72541cd0125SJakub Kruzik 0, 72641cd0125SJakub Kruzik /* 34*/ 0, 72741cd0125SJakub Kruzik 0, 72841cd0125SJakub Kruzik 0, 72941cd0125SJakub Kruzik 0, 73041cd0125SJakub Kruzik 0, 73141cd0125SJakub Kruzik /* 39*/ 0, 73241cd0125SJakub Kruzik 0, 73341cd0125SJakub Kruzik 0, 73441cd0125SJakub Kruzik 0, 73541cd0125SJakub Kruzik 0, 73641cd0125SJakub Kruzik /* 44*/ 0, 73741cd0125SJakub Kruzik MatScale_Composite, 73841cd0125SJakub Kruzik MatShift_Basic, 73941cd0125SJakub Kruzik 0, 74041cd0125SJakub Kruzik 0, 74141cd0125SJakub Kruzik /* 49*/ 0, 74241cd0125SJakub Kruzik 0, 74341cd0125SJakub Kruzik 0, 74441cd0125SJakub Kruzik 0, 74541cd0125SJakub Kruzik 0, 74641cd0125SJakub Kruzik /* 54*/ 0, 74741cd0125SJakub Kruzik 0, 74841cd0125SJakub Kruzik 0, 74941cd0125SJakub Kruzik 0, 75041cd0125SJakub Kruzik 0, 75141cd0125SJakub Kruzik /* 59*/ 0, 75241cd0125SJakub Kruzik MatDestroy_Composite, 75341cd0125SJakub Kruzik 0, 75441cd0125SJakub Kruzik 0, 75541cd0125SJakub Kruzik 0, 75641cd0125SJakub Kruzik /* 64*/ 0, 75741cd0125SJakub Kruzik 0, 75841cd0125SJakub Kruzik 0, 75941cd0125SJakub Kruzik 0, 76041cd0125SJakub Kruzik 0, 76141cd0125SJakub Kruzik /* 69*/ 0, 76241cd0125SJakub Kruzik 0, 76341cd0125SJakub Kruzik 0, 76441cd0125SJakub Kruzik 0, 76541cd0125SJakub Kruzik 0, 76641cd0125SJakub Kruzik /* 74*/ 0, 76741cd0125SJakub Kruzik 0, 7684b2558d6SJakub Kruzik MatSetFromOptions_Composite, 76941cd0125SJakub Kruzik 0, 77041cd0125SJakub Kruzik 0, 77141cd0125SJakub Kruzik /* 79*/ 0, 77241cd0125SJakub Kruzik 0, 77341cd0125SJakub Kruzik 0, 77441cd0125SJakub Kruzik 0, 77541cd0125SJakub Kruzik 0, 77641cd0125SJakub Kruzik /* 84*/ 0, 77741cd0125SJakub Kruzik 0, 77841cd0125SJakub Kruzik 0, 77941cd0125SJakub Kruzik 0, 78041cd0125SJakub Kruzik 0, 78141cd0125SJakub Kruzik /* 89*/ 0, 78241cd0125SJakub Kruzik 0, 78341cd0125SJakub Kruzik 0, 78441cd0125SJakub Kruzik 0, 78541cd0125SJakub Kruzik 0, 78641cd0125SJakub Kruzik /* 94*/ 0, 78741cd0125SJakub Kruzik 0, 78841cd0125SJakub Kruzik 0, 78941cd0125SJakub Kruzik 0, 79041cd0125SJakub Kruzik 0, 79141cd0125SJakub Kruzik /*99*/ 0, 79241cd0125SJakub Kruzik 0, 79341cd0125SJakub Kruzik 0, 79441cd0125SJakub Kruzik 0, 79541cd0125SJakub Kruzik 0, 79641cd0125SJakub Kruzik /*104*/ 0, 79741cd0125SJakub Kruzik 0, 79841cd0125SJakub Kruzik 0, 79941cd0125SJakub Kruzik 0, 80041cd0125SJakub Kruzik 0, 80141cd0125SJakub Kruzik /*109*/ 0, 80241cd0125SJakub Kruzik 0, 80341cd0125SJakub Kruzik 0, 80441cd0125SJakub Kruzik 0, 80541cd0125SJakub Kruzik 0, 80641cd0125SJakub Kruzik /*114*/ 0, 80741cd0125SJakub Kruzik 0, 80841cd0125SJakub Kruzik 0, 80941cd0125SJakub Kruzik 0, 81041cd0125SJakub Kruzik 0, 81141cd0125SJakub Kruzik /*119*/ 0, 81241cd0125SJakub Kruzik 0, 81341cd0125SJakub Kruzik 0, 81441cd0125SJakub Kruzik 0, 81541cd0125SJakub Kruzik 0, 81641cd0125SJakub Kruzik /*124*/ 0, 81741cd0125SJakub Kruzik 0, 81841cd0125SJakub Kruzik 0, 81941cd0125SJakub Kruzik 0, 82041cd0125SJakub Kruzik 0, 82141cd0125SJakub Kruzik /*129*/ 0, 82241cd0125SJakub Kruzik 0, 82341cd0125SJakub Kruzik 0, 82441cd0125SJakub Kruzik 0, 82541cd0125SJakub Kruzik 0, 82641cd0125SJakub Kruzik /*134*/ 0, 82741cd0125SJakub Kruzik 0, 82841cd0125SJakub Kruzik 0, 82941cd0125SJakub Kruzik 0, 83041cd0125SJakub Kruzik 0, 83141cd0125SJakub Kruzik /*139*/ 0, 83241cd0125SJakub Kruzik 0, 83341cd0125SJakub Kruzik 0 83441cd0125SJakub Kruzik }; 83541cd0125SJakub Kruzik 83641cd0125SJakub Kruzik /*MC 83741cd0125SJakub Kruzik MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices. 83841cd0125SJakub Kruzik The matrices need to have a correct size and parallel layout for the sum or product to be valid. 83941cd0125SJakub Kruzik 84041cd0125SJakub Kruzik Notes: 84141cd0125SJakub Kruzik to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE); 84241cd0125SJakub Kruzik 84341cd0125SJakub Kruzik Level: advanced 84441cd0125SJakub Kruzik 8458c8613bfSJakub Kruzik .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeSetType(), MatCompositeGetType(), MatCompositeMerge(), MatCompositeSetMergeType(), MatCompositeGetNumberMat(), MatCompositeGetMat() 84641cd0125SJakub Kruzik M*/ 84741cd0125SJakub Kruzik 84841cd0125SJakub Kruzik PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat A) 84941cd0125SJakub Kruzik { 85041cd0125SJakub Kruzik Mat_Composite *b; 85141cd0125SJakub Kruzik PetscErrorCode ierr; 85241cd0125SJakub Kruzik 85341cd0125SJakub Kruzik PetscFunctionBegin; 85441cd0125SJakub Kruzik ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 85541cd0125SJakub Kruzik A->data = (void*)b; 85641cd0125SJakub Kruzik ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 85741cd0125SJakub Kruzik 85841cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 85941cd0125SJakub Kruzik ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 86041cd0125SJakub Kruzik 86141cd0125SJakub Kruzik A->assembled = PETSC_TRUE; 86241cd0125SJakub Kruzik A->preallocated = PETSC_TRUE; 86341cd0125SJakub Kruzik b->type = MAT_COMPOSITE_ADDITIVE; 86441cd0125SJakub Kruzik b->scale = 1.0; 86541cd0125SJakub Kruzik b->nmat = 0; 8664b2558d6SJakub Kruzik b->merge = PETSC_FALSE; 8678c8613bfSJakub Kruzik b->mergetype = MAT_COMPOSITE_MERGE_RIGHT; 86841cd0125SJakub Kruzik 86941cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeAddMat_C",MatCompositeAddMat_Composite);CHKERRQ(ierr); 87041cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetType_C",MatCompositeSetType_Composite);CHKERRQ(ierr); 8716dbc55e5SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetType_C",MatCompositeGetType_Composite);CHKERRQ(ierr); 8728c8613bfSJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeSetMergeType_C",MatCompositeSetMergeType_Composite);CHKERRQ(ierr); 87341cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeMerge_C",MatCompositeMerge_Composite);CHKERRQ(ierr); 8746d0add67SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetNumberMat_C",MatCompositeGetNumberMat_Composite);CHKERRQ(ierr); 87541cd0125SJakub Kruzik ierr = PetscObjectComposeFunction((PetscObject)A,"MatCompositeGetMat_C",MatCompositeGetMat_Composite);CHKERRQ(ierr); 87641cd0125SJakub Kruzik 87741cd0125SJakub Kruzik ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 87841cd0125SJakub Kruzik PetscFunctionReturn(0); 87941cd0125SJakub Kruzik } 88041cd0125SJakub Kruzik 881