16f79c3a4SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 54a2ae208SSatish Balay #define __FUNCT__ "MatAXPY" 606be10caSBarry Smith /*@ 721c89e3eSBarry Smith MatAXPY - Computes Y = a*X + Y. 86f79c3a4SBarry Smith 93f9fe445SBarry Smith Logically Collective on Mat 10fee21e36SBarry Smith 1198a79cdbSBarry Smith Input Parameters: 12607cd303SBarry Smith + a - the scalar multiplier 13607cd303SBarry Smith . X - the first matrix 14607cd303SBarry Smith . Y - the second matrix 15407f6b05SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 16407f6b05SHong Zhang or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 1798a79cdbSBarry Smith 182860a424SLois Curfman McInnes Level: intermediate 192860a424SLois Curfman McInnes 209cf4f1e8SLois Curfman McInnes .keywords: matrix, add 21d4bb536fSBarry Smith 222860a424SLois Curfman McInnes .seealso: MatAYPX() 2306be10caSBarry Smith @*/ 247087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 256f79c3a4SBarry Smith { 266849ba73SBarry Smith PetscErrorCode ierr; 27c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 286f79c3a4SBarry Smith 293a40ed3dSBarry Smith PetscFunctionBegin; 300700a824SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 310700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 32c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 33273d9f13SBarry Smith ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 34273d9f13SBarry Smith ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 35e32f2f54SBarry Smith if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %D %D %D %D",m1,m2,n1,n2); 361987afe7SBarry Smith 37e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 38ab784542SHong Zhang if (Y->ops->axpy) { 39f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 40d4bb536fSBarry Smith } else { 41f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 42607cd303SBarry Smith } 43e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 442eea6e16SBarry Smith #if defined(PETSC_HAVE_CUSP) 452eea6e16SBarry Smith if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 462eea6e16SBarry Smith Y->valid_GPU_matrix = PETSC_CUSP_CPU; 472eea6e16SBarry Smith } 482eea6e16SBarry Smith #endif 49d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 50d67ff14aSKarl Rupp if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 51d67ff14aSKarl Rupp Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 52d67ff14aSKarl Rupp } 53d67ff14aSKarl Rupp #endif 54607cd303SBarry Smith PetscFunctionReturn(0); 55607cd303SBarry Smith } 56607cd303SBarry Smith 57607cd303SBarry Smith #undef __FUNCT__ 58607cd303SBarry Smith #define __FUNCT__ "MatAXPY_Basic" 59f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 60607cd303SBarry Smith { 6138baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 626849ba73SBarry Smith PetscErrorCode ierr; 6338baddfdSBarry Smith const PetscInt *row; 64b3cc6726SBarry Smith PetscScalar *val; 65b3cc6726SBarry Smith const PetscScalar *vals; 66607cd303SBarry Smith 67607cd303SBarry Smith PetscFunctionBegin; 688dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 6990f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 70f4df32b1SMatthew Knepley if (a == 1.0) { 71d4bb536fSBarry Smith for (i = start; i < end; i++) { 72d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 73d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 74d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 75d4bb536fSBarry Smith } 76d4bb536fSBarry Smith } else { 77854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 7806be10caSBarry Smith for (i=start; i<end; i++) { 79b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 8006be10caSBarry Smith for (j=0; j<ncols; j++) { 81f4df32b1SMatthew Knepley val[j] = a*vals[j]; 826f79c3a4SBarry Smith } 83b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 84b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 856f79c3a4SBarry Smith } 86b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 87d4bb536fSBarry Smith } 886d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 896d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 903a40ed3dSBarry Smith PetscFunctionReturn(0); 916f79c3a4SBarry Smith } 92052efed2SBarry Smith 934a2ae208SSatish Balay #undef __FUNCT__ 94ec7775f6SShri Abhyankar #define __FUNCT__ "MatAXPY_BasicWithPreallocation" 95ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 96ec7775f6SShri Abhyankar { 97ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 98ec7775f6SShri Abhyankar PetscErrorCode ierr; 99ec7775f6SShri Abhyankar const PetscInt *row; 100ec7775f6SShri Abhyankar PetscScalar *val; 101ec7775f6SShri Abhyankar const PetscScalar *vals; 102ec7775f6SShri Abhyankar 103ec7775f6SShri Abhyankar PetscFunctionBegin; 104ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 105ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 106ec7775f6SShri Abhyankar if (a == 1.0) { 107ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 108ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 109ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 110ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 111ec7775f6SShri Abhyankar 112ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 113ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 114ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 115ec7775f6SShri Abhyankar } 116ec7775f6SShri Abhyankar } else { 117854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 118ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 119ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 120ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 121ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 122ec7775f6SShri Abhyankar 123ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 124ec7775f6SShri Abhyankar for (j=0; j<ncols; j++) { 125ec7775f6SShri Abhyankar val[j] = a*vals[j]; 126ec7775f6SShri Abhyankar } 127ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 128ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 129ec7775f6SShri Abhyankar } 130ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 131ec7775f6SShri Abhyankar } 132ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 134ec7775f6SShri Abhyankar PetscFunctionReturn(0); 135ec7775f6SShri Abhyankar } 136ec7775f6SShri Abhyankar 137ec7775f6SShri Abhyankar #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatShift" 139052efed2SBarry Smith /*@ 14087828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 141052efed2SBarry Smith 1423f9fe445SBarry Smith Neighbor-wise Collective on Mat 143fee21e36SBarry Smith 14498a79cdbSBarry Smith Input Parameters: 14598a79cdbSBarry Smith + Y - the matrices 14687828ca2SBarry Smith - a - the PetscScalar 14798a79cdbSBarry Smith 1482860a424SLois Curfman McInnes Level: intermediate 1492860a424SLois Curfman McInnes 1506f33a894SBarry Smith Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 1516f33a894SBarry Smith fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 1526f33a894SBarry Smith entry). 1536f33a894SBarry Smith 1546f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 1556f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 1566f33a894SBarry Smith 157052efed2SBarry Smith .keywords: matrix, add, shift 1586b9ee512SLois Curfman McInnes 159f56f2b3fSBarry Smith .seealso: MatDiagonalSet() 160052efed2SBarry Smith @*/ 1617087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 162052efed2SBarry Smith { 1636849ba73SBarry Smith PetscErrorCode ierr; 164052efed2SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 1660700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 167ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 168ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1694994cf47SJed Brown MatCheckPreallocated(Y,1); 170b50b34bdSBarry Smith 171f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 1727d68702bSBarry Smith 1732eea6e16SBarry Smith #if defined(PETSC_HAVE_CUSP) 1742eea6e16SBarry Smith if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 1752eea6e16SBarry Smith Y->valid_GPU_matrix = PETSC_CUSP_CPU; 1762eea6e16SBarry Smith } 1772eea6e16SBarry Smith #endif 178d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 179d67ff14aSKarl Rupp if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 180d67ff14aSKarl Rupp Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 181d67ff14aSKarl Rupp } 182d67ff14aSKarl Rupp #endif 1833a40ed3dSBarry Smith PetscFunctionReturn(0); 184052efed2SBarry Smith } 1856d84be18SBarry Smith 1864a2ae208SSatish Balay #undef __FUNCT__ 18709f38230SBarry Smith #define __FUNCT__ "MatDiagonalSet_Default" 1887087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 18909f38230SBarry Smith { 19009f38230SBarry Smith PetscErrorCode ierr; 19167576b8bSBarry Smith PetscInt i,start,end; 19209f38230SBarry Smith PetscScalar *v; 19309f38230SBarry Smith 19409f38230SBarry Smith PetscFunctionBegin; 19509f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 19609f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 19709f38230SBarry Smith for (i=start; i<end; i++) { 19809f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 19909f38230SBarry Smith } 20009f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 20109f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20209f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20309f38230SBarry Smith PetscFunctionReturn(0); 20409f38230SBarry Smith } 20509f38230SBarry Smith 20609f38230SBarry Smith #undef __FUNCT__ 2074a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalSet" 2086d84be18SBarry Smith /*@ 209f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 210f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 211f56f2b3fSBarry Smith INSERT_VALUES. 2126d84be18SBarry Smith 2136d84be18SBarry Smith Input Parameters: 21498a79cdbSBarry Smith + Y - the input matrix 215f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 216f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 2176d84be18SBarry Smith 2183f9fe445SBarry Smith Neighbor-wise Collective on Mat and Vec 219fee21e36SBarry Smith 2206f33a894SBarry Smith Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2216f33a894SBarry Smith fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 2226f33a894SBarry Smith entry). 2236f33a894SBarry Smith 2242860a424SLois Curfman McInnes Level: intermediate 2252860a424SLois Curfman McInnes 2266b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal 2276b9ee512SLois Curfman McInnes 2286b9ee512SLois Curfman McInnes .seealso: MatShift() 2296d84be18SBarry Smith @*/ 2307087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 2316d84be18SBarry Smith { 2326849ba73SBarry Smith PetscErrorCode ierr; 23367576b8bSBarry Smith PetscInt matlocal,veclocal; 2346d84be18SBarry Smith 2353a40ed3dSBarry Smith PetscFunctionBegin; 2360700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 2370700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 23867576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 23967576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 24067576b8bSBarry Smith if (matlocal != veclocal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number local rows of matrix %D does not match that of vector for diagonal %D",matlocal,veclocal); 241f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 242f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 24394d884c6SBarry Smith } else { 24409f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 2456d84be18SBarry Smith } 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2476d84be18SBarry Smith } 248d4bb536fSBarry Smith 2494a2ae208SSatish Balay #undef __FUNCT__ 2504a2ae208SSatish Balay #define __FUNCT__ "MatAYPX" 251d4bb536fSBarry Smith /*@ 25204aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 253d4bb536fSBarry Smith 2543f9fe445SBarry Smith Logically on Mat 255fee21e36SBarry Smith 25698a79cdbSBarry Smith Input Parameters: 25704aac2b0SHong Zhang + a - the PetscScalar multiplier 25804aac2b0SHong Zhang . Y - the first matrix 25904aac2b0SHong Zhang . X - the second matrix 26004aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 261d4bb536fSBarry Smith 2622860a424SLois Curfman McInnes Level: intermediate 2632860a424SLois Curfman McInnes 264d4bb536fSBarry Smith .keywords: matrix, add 265d4bb536fSBarry Smith 2662860a424SLois Curfman McInnes .seealso: MatAXPY() 267d4bb536fSBarry Smith @*/ 2687087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 269d4bb536fSBarry Smith { 27087828ca2SBarry Smith PetscScalar one = 1.0; 2716849ba73SBarry Smith PetscErrorCode ierr; 27238baddfdSBarry Smith PetscInt mX,mY,nX,nY; 273d4bb536fSBarry Smith 2743a40ed3dSBarry Smith PetscFunctionBegin; 275c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 2760700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 277c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 278329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 279329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 280e32f2f54SBarry Smith if (mX != mY || nX != nY) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrices: %D %D first %D %D second",mX,mY,nX,nY); 281d4bb536fSBarry Smith 282f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 283cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 2843a40ed3dSBarry Smith PetscFunctionReturn(0); 285d4bb536fSBarry Smith } 286b0a32e0cSBarry Smith 2874a2ae208SSatish Balay #undef __FUNCT__ 2884a2ae208SSatish Balay #define __FUNCT__ "MatComputeExplicitOperator" 289b0a32e0cSBarry Smith /*@ 290b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 291b0a32e0cSBarry Smith 292b0a32e0cSBarry Smith Collective on Mat 293b0a32e0cSBarry Smith 294b0a32e0cSBarry Smith Input Parameter: 295b0a32e0cSBarry Smith . inmat - the matrix 296b0a32e0cSBarry Smith 297b0a32e0cSBarry Smith Output Parameter: 298b0a32e0cSBarry Smith . mat - the explict preconditioned operator 299b0a32e0cSBarry Smith 300b0a32e0cSBarry Smith Notes: 301b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 302b0a32e0cSBarry Smith identity matrix. 303b0a32e0cSBarry Smith 304b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 305b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 306b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 307b0a32e0cSBarry Smith 308b0a32e0cSBarry Smith Level: advanced 309b0a32e0cSBarry Smith 310b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 311b0a32e0cSBarry Smith @*/ 3127087cfbeSBarry Smith PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 313b0a32e0cSBarry Smith { 314b0a32e0cSBarry Smith Vec in,out; 315dfbe8321SBarry Smith PetscErrorCode ierr; 3160b12b109SJed Brown PetscInt i,m,n,M,N,*rows,start,end; 317b0a32e0cSBarry Smith MPI_Comm comm; 31887828ca2SBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 31938baddfdSBarry Smith PetscMPIInt size; 320b0a32e0cSBarry Smith 321b0a32e0cSBarry Smith PetscFunctionBegin; 3220700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 3234482741eSBarry Smith PetscValidPointer(mat,2); 324b0a32e0cSBarry Smith 325ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 326b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 327b0a32e0cSBarry Smith 3280b12b109SJed Brown ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 3290b12b109SJed Brown ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 3302a7a6963SBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 3310b12b109SJed Brown ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 3320b12b109SJed Brown ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 333785e854fSJed Brown ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 3348865f1eaSKarl Rupp for (i=0; i<m; i++) rows[i] = start + i; 335b0a32e0cSBarry Smith 336f69a0ea3SMatthew Knepley ierr = MatCreate(comm,mat);CHKERRQ(ierr); 3370b12b109SJed Brown ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 338b0a32e0cSBarry Smith if (size == 1) { 339be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 3400298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 341b0a32e0cSBarry Smith } else { 342be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 3430298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 344b0a32e0cSBarry Smith } 345b0a32e0cSBarry Smith 3460b12b109SJed Brown for (i=0; i<N; i++) { 347b0a32e0cSBarry Smith 3482dcb1b2aSMatthew Knepley ierr = VecSet(in,zero);CHKERRQ(ierr); 349b0a32e0cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 350b0a32e0cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 351b0a32e0cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 352b0a32e0cSBarry Smith 353b0a32e0cSBarry Smith ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 354b0a32e0cSBarry Smith 355b0a32e0cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 356b0a32e0cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 357b0a32e0cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 358b0a32e0cSBarry Smith 359b0a32e0cSBarry Smith } 360b0a32e0cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 3616bf464f9SBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 3626bf464f9SBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 363b0a32e0cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 364b0a32e0cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 36524f910e3SHong Zhang PetscFunctionReturn(0); 36624f910e3SHong Zhang } 3674325cce7SMatthew G Knepley 3684325cce7SMatthew G Knepley #undef __FUNCT__ 3694325cce7SMatthew G Knepley #define __FUNCT__ "MatChop" 3704325cce7SMatthew G Knepley /*@ 3714325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 3724325cce7SMatthew G Knepley 3734325cce7SMatthew G Knepley Input Parameters: 3744325cce7SMatthew G Knepley + A - The matrix 3754325cce7SMatthew G Knepley - tol - The zero tolerance 3764325cce7SMatthew G Knepley 3774325cce7SMatthew G Knepley Output Parameters: 3784325cce7SMatthew G Knepley . A - The chopped matrix 3794325cce7SMatthew G Knepley 3804325cce7SMatthew G Knepley Level: intermediate 3814325cce7SMatthew G Knepley 3824325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 3833fc99919SSatish Balay @*/ 3844325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 3854325cce7SMatthew G Knepley { 3864325cce7SMatthew G Knepley PetscScalar *newVals; 3874325cce7SMatthew G Knepley PetscInt *newCols; 3884325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 3894325cce7SMatthew G Knepley PetscErrorCode ierr; 3904325cce7SMatthew G Knepley 3914325cce7SMatthew G Knepley PetscFunctionBegin; 3924325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 3934325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 3944325cce7SMatthew G Knepley PetscInt ncols; 3954325cce7SMatthew G Knepley 3960298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 3974325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 3980298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 3994325cce7SMatthew G Knepley } 4004325cce7SMatthew G Knepley numRows = rEnd - rStart; 401*b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 402dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 4034325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 4044325cce7SMatthew G Knepley const PetscScalar *vals; 4054325cce7SMatthew G Knepley const PetscInt *cols; 406fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 4074325cce7SMatthew G Knepley 4084325cce7SMatthew G Knepley if (r < rEnd) { 4094325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 4104325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 4114325cce7SMatthew G Knepley newCols[c] = cols[c]; 4124325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 4134325cce7SMatthew G Knepley } 414fad45679SMatthew G. Knepley newcols = ncols; 4154325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 416fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 4174325cce7SMatthew G Knepley } 4184325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4194325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4204325cce7SMatthew G Knepley } 4214325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 4224325cce7SMatthew G Knepley PetscFunctionReturn(0); 4234325cce7SMatthew G Knepley } 424