1793850ffSBarry Smith #define PETSCMAT_DLL 2793850ffSBarry Smith 3b9147fbbSdalcinl #include "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; 8793850ffSBarry Smith Mat_CompositeLink next; 9793850ffSBarry Smith }; 10793850ffSBarry Smith 11793850ffSBarry Smith typedef struct { 12793850ffSBarry Smith Mat_CompositeLink head; 13793850ffSBarry Smith Vec work; 14793850ffSBarry Smith } Mat_Composite; 15793850ffSBarry Smith 16793850ffSBarry Smith #undef __FUNCT__ 17793850ffSBarry Smith #define __FUNCT__ "MatDestroy_Composite" 18793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat) 19793850ffSBarry Smith { 20793850ffSBarry Smith PetscErrorCode ierr; 212c33bbaeSSatish Balay Mat_Composite *shell = (Mat_Composite*)mat->data; 22793850ffSBarry Smith Mat_CompositeLink next = shell->head; 23793850ffSBarry Smith 24793850ffSBarry Smith PetscFunctionBegin; 25793850ffSBarry Smith while (next) { 26793850ffSBarry Smith ierr = MatDestroy(next->mat);CHKERRQ(ierr); 27793850ffSBarry Smith next = next->next; 28793850ffSBarry Smith } 29793850ffSBarry Smith if (shell->work) {ierr = VecDestroy(shell->work);CHKERRQ(ierr);} 30793850ffSBarry Smith ierr = PetscFree(shell);CHKERRQ(ierr); 31793850ffSBarry Smith mat->data = 0; 32793850ffSBarry Smith PetscFunctionReturn(0); 33793850ffSBarry Smith } 34793850ffSBarry Smith 35793850ffSBarry Smith #undef __FUNCT__ 36793850ffSBarry Smith #define __FUNCT__ "MatMult_Composite" 37793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y) 38793850ffSBarry Smith { 39793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 40793850ffSBarry Smith Mat_CompositeLink next = shell->head; 41793850ffSBarry Smith PetscErrorCode ierr; 42793850ffSBarry Smith 43793850ffSBarry Smith PetscFunctionBegin; 44793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 45793850ffSBarry Smith ierr = MatMult(next->mat,x,y);CHKERRQ(ierr); 46793850ffSBarry Smith while ((next = next->next)) { 47793850ffSBarry Smith ierr = MatMultAdd(next->mat,x,y,y);CHKERRQ(ierr); 48793850ffSBarry Smith } 49793850ffSBarry Smith PetscFunctionReturn(0); 50793850ffSBarry Smith } 51793850ffSBarry Smith 52793850ffSBarry Smith #undef __FUNCT__ 53793850ffSBarry Smith #define __FUNCT__ "MatMultTranspose_Composite" 54793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y) 55793850ffSBarry Smith { 56793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 57793850ffSBarry Smith Mat_CompositeLink next = shell->head; 58793850ffSBarry Smith PetscErrorCode ierr; 59793850ffSBarry Smith 60793850ffSBarry Smith PetscFunctionBegin; 61793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 62793850ffSBarry Smith ierr = MatMultTranspose(next->mat,x,y);CHKERRQ(ierr); 63793850ffSBarry Smith while ((next = next->next)) { 64793850ffSBarry Smith ierr = MatMultTransposeAdd(next->mat,x,y,y);CHKERRQ(ierr); 65793850ffSBarry Smith } 66793850ffSBarry Smith PetscFunctionReturn(0); 67793850ffSBarry Smith } 68793850ffSBarry Smith 69793850ffSBarry Smith #undef __FUNCT__ 70793850ffSBarry Smith #define __FUNCT__ "MatGetDiagonal_Composite" 71793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v) 72793850ffSBarry Smith { 73793850ffSBarry Smith Mat_Composite *shell = (Mat_Composite*)A->data; 74793850ffSBarry Smith Mat_CompositeLink next = shell->head; 75793850ffSBarry Smith PetscErrorCode ierr; 76793850ffSBarry Smith 77793850ffSBarry Smith PetscFunctionBegin; 78793850ffSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 79793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr); 80793850ffSBarry Smith if (next->next && !shell->work) { 81793850ffSBarry Smith ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr); 82793850ffSBarry Smith } 83793850ffSBarry Smith while ((next = next->next)) { 84793850ffSBarry Smith ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr); 85793850ffSBarry Smith ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr); 86793850ffSBarry Smith } 87793850ffSBarry Smith PetscFunctionReturn(0); 88793850ffSBarry Smith } 89793850ffSBarry Smith 90793850ffSBarry Smith #undef __FUNCT__ 91793850ffSBarry Smith #define __FUNCT__ "MatAssemblyEnd_Composite" 92793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t) 93793850ffSBarry Smith { 94b52f573bSBarry Smith PetscErrorCode ierr; 95b52f573bSBarry Smith PetscTruth flg; 96b52f573bSBarry Smith 97793850ffSBarry Smith PetscFunctionBegin; 98b52f573bSBarry Smith ierr = PetscOptionsHasName(Y->prefix,"-mat_composite_merge",&flg);CHKERRQ(ierr); 99b52f573bSBarry Smith if (flg) { 100b52f573bSBarry Smith ierr = MatCompositeMerge(Y);CHKERRQ(ierr); 101b52f573bSBarry Smith } 102793850ffSBarry Smith PetscFunctionReturn(0); 103793850ffSBarry Smith } 104793850ffSBarry Smith 105793850ffSBarry Smith 106793850ffSBarry Smith static struct _MatOps MatOps_Values = {0, 107793850ffSBarry Smith 0, 108793850ffSBarry Smith 0, 109793850ffSBarry Smith MatMult_Composite, 110793850ffSBarry Smith 0, 111793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite, 112793850ffSBarry Smith 0, 113793850ffSBarry Smith 0, 114793850ffSBarry Smith 0, 115793850ffSBarry Smith 0, 116793850ffSBarry Smith /*10*/ 0, 117793850ffSBarry Smith 0, 118793850ffSBarry Smith 0, 119793850ffSBarry Smith 0, 120793850ffSBarry Smith 0, 121793850ffSBarry Smith /*15*/ 0, 122793850ffSBarry Smith 0, 123793850ffSBarry Smith MatGetDiagonal_Composite, 124793850ffSBarry Smith 0, 125793850ffSBarry Smith 0, 126793850ffSBarry Smith /*20*/ 0, 127793850ffSBarry Smith MatAssemblyEnd_Composite, 128793850ffSBarry Smith 0, 129793850ffSBarry Smith 0, 130793850ffSBarry Smith 0, 131793850ffSBarry Smith /*25*/ 0, 132793850ffSBarry Smith 0, 133793850ffSBarry Smith 0, 134793850ffSBarry Smith 0, 135793850ffSBarry Smith 0, 136793850ffSBarry Smith /*30*/ 0, 137793850ffSBarry Smith 0, 138793850ffSBarry Smith 0, 139793850ffSBarry Smith 0, 140793850ffSBarry Smith 0, 141793850ffSBarry Smith /*35*/ 0, 142793850ffSBarry Smith 0, 143793850ffSBarry Smith 0, 144793850ffSBarry Smith 0, 145793850ffSBarry Smith 0, 146793850ffSBarry Smith /*40*/ 0, 147793850ffSBarry Smith 0, 148793850ffSBarry Smith 0, 149793850ffSBarry Smith 0, 150793850ffSBarry Smith 0, 151793850ffSBarry Smith /*45*/ 0, 152793850ffSBarry Smith 0, 153793850ffSBarry Smith 0, 154793850ffSBarry Smith 0, 155793850ffSBarry Smith 0, 156793850ffSBarry Smith /*50*/ 0, 157793850ffSBarry Smith 0, 158793850ffSBarry Smith 0, 159793850ffSBarry Smith 0, 160793850ffSBarry Smith 0, 161793850ffSBarry Smith /*55*/ 0, 162793850ffSBarry Smith 0, 163793850ffSBarry Smith 0, 164793850ffSBarry Smith 0, 165793850ffSBarry Smith 0, 166793850ffSBarry Smith /*60*/ 0, 167793850ffSBarry Smith MatDestroy_Composite, 168793850ffSBarry Smith 0, 169793850ffSBarry Smith 0, 170793850ffSBarry Smith 0, 171793850ffSBarry Smith /*65*/ 0, 172793850ffSBarry Smith 0, 173793850ffSBarry Smith 0, 174793850ffSBarry Smith 0, 175793850ffSBarry Smith 0, 176793850ffSBarry Smith /*70*/ 0, 177793850ffSBarry Smith 0, 178793850ffSBarry Smith 0, 179793850ffSBarry Smith 0, 180793850ffSBarry Smith 0, 181793850ffSBarry Smith /*75*/ 0, 182793850ffSBarry Smith 0, 183793850ffSBarry Smith 0, 184793850ffSBarry Smith 0, 185793850ffSBarry Smith 0, 186793850ffSBarry Smith /*80*/ 0, 187793850ffSBarry Smith 0, 188793850ffSBarry Smith 0, 189793850ffSBarry Smith 0, 190793850ffSBarry Smith 0, 191793850ffSBarry Smith /*85*/ 0, 192793850ffSBarry Smith 0, 193793850ffSBarry Smith 0, 194793850ffSBarry Smith 0, 195793850ffSBarry Smith 0, 196793850ffSBarry Smith /*90*/ 0, 197793850ffSBarry Smith 0, 198793850ffSBarry Smith 0, 199793850ffSBarry Smith 0, 200793850ffSBarry Smith 0, 201793850ffSBarry Smith /*95*/ 0, 202793850ffSBarry Smith 0, 203793850ffSBarry Smith 0, 204793850ffSBarry Smith 0}; 205793850ffSBarry Smith 206793850ffSBarry Smith /*MC 207793850ffSBarry Smith MATCOMPOSITE - A matrix defined by the sum of one or more matrices (all matrices are of same size and parallel layout) 208793850ffSBarry Smith 209793850ffSBarry Smith Level: advanced 210793850ffSBarry Smith 211b52f573bSBarry Smith .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge() 212793850ffSBarry Smith M*/ 213793850ffSBarry Smith 214793850ffSBarry Smith EXTERN_C_BEGIN 215793850ffSBarry Smith #undef __FUNCT__ 216793850ffSBarry Smith #define __FUNCT__ "MatCreate_Composite" 217793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Composite(Mat A) 218793850ffSBarry Smith { 219793850ffSBarry Smith Mat_Composite *b; 220793850ffSBarry Smith PetscErrorCode ierr; 221793850ffSBarry Smith 222793850ffSBarry Smith PetscFunctionBegin; 223*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Composite,&b);CHKERRQ(ierr); 224793850ffSBarry Smith A->data = (void*)b; 225*38f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 226793850ffSBarry Smith 2272813ca57SBarry Smith ierr = PetscMapSetBlockSize(&A->rmap,1);CHKERRQ(ierr); 2282813ca57SBarry Smith ierr = PetscMapSetBlockSize(&A->cmap,1);CHKERRQ(ierr); 2296148ca0dSBarry Smith ierr = PetscMapSetUp(&A->rmap);CHKERRQ(ierr); 2306148ca0dSBarry Smith ierr = PetscMapSetUp(&A->cmap);CHKERRQ(ierr); 231793850ffSBarry Smith 232793850ffSBarry Smith A->assembled = PETSC_TRUE; 233793850ffSBarry Smith A->preallocated = PETSC_FALSE; 234793850ffSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr); 235793850ffSBarry Smith PetscFunctionReturn(0); 236793850ffSBarry Smith } 237793850ffSBarry Smith EXTERN_C_END 238793850ffSBarry Smith 239793850ffSBarry Smith #undef __FUNCT__ 240793850ffSBarry Smith #define __FUNCT__ "MatCreateComposite" 241793850ffSBarry Smith /*@C 242793850ffSBarry Smith MatCreateComposite - Creates a matrix as the sum of zero or more matrices 243793850ffSBarry Smith 244793850ffSBarry Smith Collective on MPI_Comm 245793850ffSBarry Smith 246793850ffSBarry Smith Input Parameters: 247793850ffSBarry Smith + comm - MPI communicator 248793850ffSBarry Smith . nmat - number of matrices to put in 249793850ffSBarry Smith - mats - the matrices 250793850ffSBarry Smith 251793850ffSBarry Smith Output Parameter: 252793850ffSBarry Smith . A - the matrix 253793850ffSBarry Smith 254793850ffSBarry Smith Level: advanced 255793850ffSBarry Smith 256793850ffSBarry Smith Notes: 257793850ffSBarry Smith Alternative construction 258793850ffSBarry Smith $ MatCreate(comm,&mat); 259793850ffSBarry Smith $ MatSetSizes(mat,m,n,M,N); 260793850ffSBarry Smith $ MatSetType(mat,MATCOMPOSITE); 261793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[0]); 262793850ffSBarry Smith $ .... 263793850ffSBarry Smith $ MatCompositeAddMat(mat,mats[nmat-1]); 264b52f573bSBarry Smith $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY); 265b52f573bSBarry Smith $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY); 266793850ffSBarry Smith 267b52f573bSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge() 268793850ffSBarry Smith 269793850ffSBarry Smith @*/ 270793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat) 271793850ffSBarry Smith { 272793850ffSBarry Smith PetscErrorCode ierr; 273793850ffSBarry Smith PetscInt m,n,M,N,i; 274793850ffSBarry Smith 275793850ffSBarry Smith PetscFunctionBegin; 276793850ffSBarry Smith if (nmat < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix"); 277f3f84630SBarry Smith PetscValidPointer(mat,3); 278793850ffSBarry Smith 279793850ffSBarry Smith ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr); 280793850ffSBarry Smith ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr); 281793850ffSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 282793850ffSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 283793850ffSBarry Smith ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr); 284793850ffSBarry Smith for (i=0; i<nmat; i++) { 285793850ffSBarry Smith ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr); 286793850ffSBarry Smith } 287b52f573bSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 288b52f573bSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 289793850ffSBarry Smith PetscFunctionReturn(0); 290793850ffSBarry Smith } 291793850ffSBarry Smith 292793850ffSBarry Smith #undef __FUNCT__ 293793850ffSBarry Smith #define __FUNCT__ "MatCompositeAddMat" 294793850ffSBarry Smith /*@ 295793850ffSBarry Smith MatCompositeAddMat - add another matrix to a composite matrix 296793850ffSBarry Smith 297793850ffSBarry Smith Collective on Mat 298793850ffSBarry Smith 299793850ffSBarry Smith Input Parameters: 300793850ffSBarry Smith + mat - the composite matrix 301793850ffSBarry Smith - smat - the partial matrix 302793850ffSBarry Smith 303793850ffSBarry Smith Level: advanced 304793850ffSBarry Smith 305793850ffSBarry Smith .seealso: MatCreateComposite() 306793850ffSBarry Smith @*/ 307793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat mat,Mat smat) 308793850ffSBarry Smith { 309*38f2d2fdSLisandro Dalcin Mat_Composite *shell; 310793850ffSBarry Smith PetscErrorCode ierr; 311793850ffSBarry Smith Mat_CompositeLink ilink,next; 312793850ffSBarry Smith 313793850ffSBarry Smith PetscFunctionBegin; 314793850ffSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 315793850ffSBarry Smith PetscValidHeaderSpecific(smat,MAT_COOKIE,2); 316*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(mat,struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr); 317793850ffSBarry Smith ilink->next = 0; 318793850ffSBarry Smith ierr = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr); 319c3122656SLisandro Dalcin ilink->mat = smat; 320793850ffSBarry Smith 321*38f2d2fdSLisandro Dalcin shell = (Mat_Composite*)mat->data; 322793850ffSBarry Smith next = shell->head; 323793850ffSBarry Smith if (!next) { 324793850ffSBarry Smith shell->head = ilink; 325793850ffSBarry Smith } else { 326793850ffSBarry Smith while (next->next) { 327793850ffSBarry Smith next = next->next; 328793850ffSBarry Smith } 329793850ffSBarry Smith next->next = ilink; 330793850ffSBarry Smith } 331793850ffSBarry Smith PetscFunctionReturn(0); 332793850ffSBarry Smith } 333793850ffSBarry Smith 334b52f573bSBarry Smith #undef __FUNCT__ 335b52f573bSBarry Smith #define __FUNCT__ "MatCompositeMerge" 336b52f573bSBarry Smith /*@C 337b52f573bSBarry Smith MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix 338b52f573bSBarry Smith by summing all the matrices inside the composite matrix. 339b52f573bSBarry Smith 340b52f573bSBarry Smith Collective on MPI_Comm 341b52f573bSBarry Smith 342b52f573bSBarry Smith Input Parameters: 343b52f573bSBarry Smith . mat - the composite matrix 344b52f573bSBarry Smith 345b52f573bSBarry Smith 346b52f573bSBarry Smith Options Database: 347b52f573bSBarry Smith . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked) 348b52f573bSBarry Smith 349b52f573bSBarry Smith Level: advanced 350b52f573bSBarry Smith 351b52f573bSBarry Smith Notes: 352b52f573bSBarry Smith The MatType of the resulting matrix will be the same as the MatType of the FIRST 353b52f573bSBarry Smith matrix in the composite matrix. 354b52f573bSBarry Smith 355b52f573bSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE 356b52f573bSBarry Smith 357b52f573bSBarry Smith @*/ 358b52f573bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat mat) 359b52f573bSBarry Smith { 360b52f573bSBarry Smith Mat_Composite *shell = (Mat_Composite*)mat->data; 361b52f573bSBarry Smith Mat_CompositeLink next = shell->head; 362b52f573bSBarry Smith PetscErrorCode ierr; 363b52f573bSBarry Smith Mat tmat; 364b52f573bSBarry Smith 365b52f573bSBarry Smith PetscFunctionBegin; 366b52f573bSBarry Smith if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()"); 367b52f573bSBarry Smith 368b52f573bSBarry Smith PetscFunctionBegin; 369b52f573bSBarry Smith ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr); 370b52f573bSBarry Smith while ((next = next->next)) { 371b52f573bSBarry Smith ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 372b52f573bSBarry Smith } 373b52f573bSBarry Smith ierr = MatDestroy_Composite(mat);CHKERRQ(ierr); 374b52f573bSBarry Smith 375b52f573bSBarry Smith PetscFunctionReturn(0); 376b52f573bSBarry Smith } 377