16f79c3a4SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 406be10caSBarry Smith /*@ 521c89e3eSBarry Smith MatAXPY - Computes Y = a*X + Y. 66f79c3a4SBarry Smith 73f9fe445SBarry Smith Logically Collective on Mat 8fee21e36SBarry Smith 998a79cdbSBarry Smith Input Parameters: 10607cd303SBarry Smith + a - the scalar multiplier 11607cd303SBarry Smith . X - the first matrix 12607cd303SBarry Smith . Y - the second matrix 13407f6b05SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 14407f6b05SHong Zhang or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 1598a79cdbSBarry Smith 162860a424SLois Curfman McInnes Level: intermediate 172860a424SLois Curfman McInnes 189cf4f1e8SLois Curfman McInnes .keywords: matrix, add 19d4bb536fSBarry Smith 202860a424SLois Curfman McInnes .seealso: MatAYPX() 2106be10caSBarry Smith @*/ 227087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 236f79c3a4SBarry Smith { 246849ba73SBarry Smith PetscErrorCode ierr; 25c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 260ff8bee4SStefano Zampini PetscBool sametype; 276f79c3a4SBarry Smith 283a40ed3dSBarry Smith PetscFunctionBegin; 290700a824SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 300700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 31c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 32273d9f13SBarry Smith ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 33273d9f13SBarry Smith ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 34e32f2f54SBarry 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); 351987afe7SBarry Smith 360ff8bee4SStefano Zampini ierr = PetscStrcmp(((PetscObject)X)->type_name,((PetscObject)Y)->type_name,&sametype);CHKERRQ(ierr); 37e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 380ff8bee4SStefano Zampini if (Y->ops->axpy && sametype) { 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 } 48c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VIENNACL) 49d67ff14aSKarl Rupp if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 50d67ff14aSKarl Rupp Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 51d67ff14aSKarl Rupp } 52c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VECCUDA) 53c41cb2e2SAlejandro Lamas Daviña if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 54c41cb2e2SAlejandro Lamas Daviña Y->valid_GPU_matrix = PETSC_CUDA_CPU; 55c41cb2e2SAlejandro Lamas Daviña } 56d67ff14aSKarl Rupp #endif 57607cd303SBarry Smith PetscFunctionReturn(0); 58607cd303SBarry Smith } 59607cd303SBarry Smith 60f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 61607cd303SBarry Smith { 6238baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 636849ba73SBarry Smith PetscErrorCode ierr; 6438baddfdSBarry Smith const PetscInt *row; 65b3cc6726SBarry Smith PetscScalar *val; 66b3cc6726SBarry Smith const PetscScalar *vals; 67607cd303SBarry Smith 68607cd303SBarry Smith PetscFunctionBegin; 698dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 7090f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 71f4df32b1SMatthew Knepley if (a == 1.0) { 72d4bb536fSBarry Smith for (i = start; i < end; i++) { 73d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 74d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 75d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 76d4bb536fSBarry Smith } 77d4bb536fSBarry Smith } else { 78854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 7906be10caSBarry Smith for (i=start; i<end; i++) { 80b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 8106be10caSBarry Smith for (j=0; j<ncols; j++) { 82f4df32b1SMatthew Knepley val[j] = a*vals[j]; 836f79c3a4SBarry Smith } 84b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 85b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 866f79c3a4SBarry Smith } 87b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 88d4bb536fSBarry Smith } 896d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 906d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 913a40ed3dSBarry Smith PetscFunctionReturn(0); 926f79c3a4SBarry Smith } 93052efed2SBarry Smith 94ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 95ec7775f6SShri Abhyankar { 96ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 97ec7775f6SShri Abhyankar PetscErrorCode ierr; 98ec7775f6SShri Abhyankar const PetscInt *row; 99ec7775f6SShri Abhyankar PetscScalar *val; 100ec7775f6SShri Abhyankar const PetscScalar *vals; 101ec7775f6SShri Abhyankar 102ec7775f6SShri Abhyankar PetscFunctionBegin; 103ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 104ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 105ec7775f6SShri Abhyankar if (a == 1.0) { 106ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 107ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 108ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 109ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 110ec7775f6SShri Abhyankar 111ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 112ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 113ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 114ec7775f6SShri Abhyankar } 115ec7775f6SShri Abhyankar } else { 116854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 117ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 118ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 119ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 120ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 121ec7775f6SShri Abhyankar 122ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 123ec7775f6SShri Abhyankar for (j=0; j<ncols; j++) { 124ec7775f6SShri Abhyankar val[j] = a*vals[j]; 125ec7775f6SShri Abhyankar } 126ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 127ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 128ec7775f6SShri Abhyankar } 129ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 130ec7775f6SShri Abhyankar } 131ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133ec7775f6SShri Abhyankar PetscFunctionReturn(0); 134ec7775f6SShri Abhyankar } 135ec7775f6SShri Abhyankar 136052efed2SBarry Smith /*@ 13787828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 138052efed2SBarry Smith 1393f9fe445SBarry Smith Neighbor-wise Collective on Mat 140fee21e36SBarry Smith 14198a79cdbSBarry Smith Input Parameters: 14298a79cdbSBarry Smith + Y - the matrices 14387828ca2SBarry Smith - a - the PetscScalar 14498a79cdbSBarry Smith 1452860a424SLois Curfman McInnes Level: intermediate 1462860a424SLois Curfman McInnes 1476f33a894SBarry Smith Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 1486f33a894SBarry 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 1496f33a894SBarry Smith entry). 1506f33a894SBarry Smith 1516f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 1526f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 1536f33a894SBarry Smith 154052efed2SBarry Smith .keywords: matrix, add, shift 1556b9ee512SLois Curfman McInnes 156f56f2b3fSBarry Smith .seealso: MatDiagonalSet() 157052efed2SBarry Smith @*/ 1587087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 159052efed2SBarry Smith { 1606849ba73SBarry Smith PetscErrorCode ierr; 161052efed2SBarry Smith 1623a40ed3dSBarry Smith PetscFunctionBegin; 1630700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 164ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 165ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1664994cf47SJed Brown MatCheckPreallocated(Y,1); 167b50b34bdSBarry Smith 1681c738b47SStefano Zampini if (Y->ops->shift) { 169f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 1701c738b47SStefano Zampini } else { 1711c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 1721c738b47SStefano Zampini } 1737d68702bSBarry Smith 1742eea6e16SBarry Smith #if defined(PETSC_HAVE_CUSP) 1752eea6e16SBarry Smith if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 1762eea6e16SBarry Smith Y->valid_GPU_matrix = PETSC_CUSP_CPU; 1772eea6e16SBarry Smith } 178c41cb2e2SAlejandro Lamas Daviña #elif 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 } 182c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VECCUDA) 183c41cb2e2SAlejandro Lamas Daviña if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 184c41cb2e2SAlejandro Lamas Daviña Y->valid_GPU_matrix = PETSC_CUDA_CPU; 185c41cb2e2SAlejandro Lamas Daviña } 186d67ff14aSKarl Rupp #endif 1873a40ed3dSBarry Smith PetscFunctionReturn(0); 188052efed2SBarry Smith } 1896d84be18SBarry Smith 1907087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 19109f38230SBarry Smith { 19209f38230SBarry Smith PetscErrorCode ierr; 19367576b8bSBarry Smith PetscInt i,start,end; 19409f38230SBarry Smith PetscScalar *v; 19509f38230SBarry Smith 19609f38230SBarry Smith PetscFunctionBegin; 19709f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 19809f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 19909f38230SBarry Smith for (i=start; i<end; i++) { 20009f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 20109f38230SBarry Smith } 20209f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 20309f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20409f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20509f38230SBarry Smith PetscFunctionReturn(0); 20609f38230SBarry Smith } 20709f38230SBarry Smith 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 249d4bb536fSBarry Smith /*@ 25004aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 251d4bb536fSBarry Smith 2523f9fe445SBarry Smith Logically on Mat 253fee21e36SBarry Smith 25498a79cdbSBarry Smith Input Parameters: 25504aac2b0SHong Zhang + a - the PetscScalar multiplier 25604aac2b0SHong Zhang . Y - the first matrix 25704aac2b0SHong Zhang . X - the second matrix 25804aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 259d4bb536fSBarry Smith 2602860a424SLois Curfman McInnes Level: intermediate 2612860a424SLois Curfman McInnes 262d4bb536fSBarry Smith .keywords: matrix, add 263d4bb536fSBarry Smith 2642860a424SLois Curfman McInnes .seealso: MatAXPY() 265d4bb536fSBarry Smith @*/ 2667087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 267d4bb536fSBarry Smith { 26887828ca2SBarry Smith PetscScalar one = 1.0; 2696849ba73SBarry Smith PetscErrorCode ierr; 27038baddfdSBarry Smith PetscInt mX,mY,nX,nY; 271d4bb536fSBarry Smith 2723a40ed3dSBarry Smith PetscFunctionBegin; 273c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 2740700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 275c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 276329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 277329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 278e32f2f54SBarry 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); 279d4bb536fSBarry Smith 280f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 281cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 2823a40ed3dSBarry Smith PetscFunctionReturn(0); 283d4bb536fSBarry Smith } 284b0a32e0cSBarry Smith 285b0a32e0cSBarry Smith /*@ 286b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 287b0a32e0cSBarry Smith 288b0a32e0cSBarry Smith Collective on Mat 289b0a32e0cSBarry Smith 290b0a32e0cSBarry Smith Input Parameter: 291b0a32e0cSBarry Smith . inmat - the matrix 292b0a32e0cSBarry Smith 293b0a32e0cSBarry Smith Output Parameter: 294*f3b1f45cSBarry Smith . mat - the explict operator 295b0a32e0cSBarry Smith 296b0a32e0cSBarry Smith Notes: 297b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 298b0a32e0cSBarry Smith identity matrix. 299b0a32e0cSBarry Smith 300b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 301b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 302b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 303b0a32e0cSBarry Smith 304b0a32e0cSBarry Smith Level: advanced 305b0a32e0cSBarry Smith 306b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 307b0a32e0cSBarry Smith @*/ 3087087cfbeSBarry Smith PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 309b0a32e0cSBarry Smith { 310b0a32e0cSBarry Smith Vec in,out; 311dfbe8321SBarry Smith PetscErrorCode ierr; 3120b12b109SJed Brown PetscInt i,m,n,M,N,*rows,start,end; 313b0a32e0cSBarry Smith MPI_Comm comm; 31487828ca2SBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 31538baddfdSBarry Smith PetscMPIInt size; 316b0a32e0cSBarry Smith 317b0a32e0cSBarry Smith PetscFunctionBegin; 3180700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 3194482741eSBarry Smith PetscValidPointer(mat,2); 320b0a32e0cSBarry Smith 321ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 322b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 323b0a32e0cSBarry Smith 3240b12b109SJed Brown ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 3250b12b109SJed Brown ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 3262a7a6963SBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 3270b12b109SJed Brown ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 3280b12b109SJed Brown ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 329785e854fSJed Brown ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 3308865f1eaSKarl Rupp for (i=0; i<m; i++) rows[i] = start + i; 331b0a32e0cSBarry Smith 332f69a0ea3SMatthew Knepley ierr = MatCreate(comm,mat);CHKERRQ(ierr); 3330b12b109SJed Brown ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 334b0a32e0cSBarry Smith if (size == 1) { 335be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 3360298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 337b0a32e0cSBarry Smith } else { 338be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 3390298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 340b0a32e0cSBarry Smith } 341b0a32e0cSBarry Smith 3420b12b109SJed Brown for (i=0; i<N; i++) { 343b0a32e0cSBarry Smith 3442dcb1b2aSMatthew Knepley ierr = VecSet(in,zero);CHKERRQ(ierr); 345b0a32e0cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 346b0a32e0cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 347b0a32e0cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 348b0a32e0cSBarry Smith 349b0a32e0cSBarry Smith ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 350b0a32e0cSBarry Smith 351b0a32e0cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 352b0a32e0cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 353b0a32e0cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 354b0a32e0cSBarry Smith 355b0a32e0cSBarry Smith } 356b0a32e0cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 3576bf464f9SBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 3586bf464f9SBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 359b0a32e0cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 360b0a32e0cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 36124f910e3SHong Zhang PetscFunctionReturn(0); 36224f910e3SHong Zhang } 3634325cce7SMatthew G Knepley 3644325cce7SMatthew G Knepley /*@ 365*f3b1f45cSBarry Smith MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of 366*f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 367*f3b1f45cSBarry Smith 368*f3b1f45cSBarry Smith Collective on Mat 369*f3b1f45cSBarry Smith 370*f3b1f45cSBarry Smith Input Parameter: 371*f3b1f45cSBarry Smith . inmat - the matrix 372*f3b1f45cSBarry Smith 373*f3b1f45cSBarry Smith Output Parameter: 374*f3b1f45cSBarry Smith . mat - the explict operator transposed 375*f3b1f45cSBarry Smith 376*f3b1f45cSBarry Smith Notes: 377*f3b1f45cSBarry Smith This computation is done by applying the transpose of the operator to columns of the 378*f3b1f45cSBarry Smith identity matrix. 379*f3b1f45cSBarry Smith 380*f3b1f45cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 381*f3b1f45cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 382*f3b1f45cSBarry Smith and is recommended for use only with relatively small systems. 383*f3b1f45cSBarry Smith 384*f3b1f45cSBarry Smith Level: advanced 385*f3b1f45cSBarry Smith 386*f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator 387*f3b1f45cSBarry Smith @*/ 388*f3b1f45cSBarry Smith PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 389*f3b1f45cSBarry Smith { 390*f3b1f45cSBarry Smith Vec in,out; 391*f3b1f45cSBarry Smith PetscErrorCode ierr; 392*f3b1f45cSBarry Smith PetscInt i,m,n,M,N,*rows,start,end; 393*f3b1f45cSBarry Smith MPI_Comm comm; 394*f3b1f45cSBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 395*f3b1f45cSBarry Smith PetscMPIInt size; 396*f3b1f45cSBarry Smith 397*f3b1f45cSBarry Smith PetscFunctionBegin; 398*f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 399*f3b1f45cSBarry Smith PetscValidPointer(mat,2); 400*f3b1f45cSBarry Smith 401*f3b1f45cSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 402*f3b1f45cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 403*f3b1f45cSBarry Smith 404*f3b1f45cSBarry Smith ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 405*f3b1f45cSBarry Smith ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 406*f3b1f45cSBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 407*f3b1f45cSBarry Smith ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 408*f3b1f45cSBarry Smith ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 409*f3b1f45cSBarry Smith ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 410*f3b1f45cSBarry Smith for (i=0; i<m; i++) rows[i] = start + i; 411*f3b1f45cSBarry Smith 412*f3b1f45cSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 413*f3b1f45cSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 414*f3b1f45cSBarry Smith if (size == 1) { 415*f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 416*f3b1f45cSBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 417*f3b1f45cSBarry Smith } else { 418*f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 419*f3b1f45cSBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 420*f3b1f45cSBarry Smith } 421*f3b1f45cSBarry Smith 422*f3b1f45cSBarry Smith for (i=0; i<N; i++) { 423*f3b1f45cSBarry Smith 424*f3b1f45cSBarry Smith ierr = VecSet(in,zero);CHKERRQ(ierr); 425*f3b1f45cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 426*f3b1f45cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 427*f3b1f45cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 428*f3b1f45cSBarry Smith 429*f3b1f45cSBarry Smith ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 430*f3b1f45cSBarry Smith 431*f3b1f45cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 432*f3b1f45cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 433*f3b1f45cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 434*f3b1f45cSBarry Smith 435*f3b1f45cSBarry Smith } 436*f3b1f45cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 437*f3b1f45cSBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 438*f3b1f45cSBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 439*f3b1f45cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 440*f3b1f45cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 441*f3b1f45cSBarry Smith PetscFunctionReturn(0); 442*f3b1f45cSBarry Smith } 443*f3b1f45cSBarry Smith 444*f3b1f45cSBarry Smith /*@ 4454325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4464325cce7SMatthew G Knepley 4474325cce7SMatthew G Knepley Input Parameters: 4484325cce7SMatthew G Knepley + A - The matrix 4494325cce7SMatthew G Knepley - tol - The zero tolerance 4504325cce7SMatthew G Knepley 4514325cce7SMatthew G Knepley Output Parameters: 4524325cce7SMatthew G Knepley . A - The chopped matrix 4534325cce7SMatthew G Knepley 4544325cce7SMatthew G Knepley Level: intermediate 4554325cce7SMatthew G Knepley 4564325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 4573fc99919SSatish Balay @*/ 4584325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 4594325cce7SMatthew G Knepley { 4604325cce7SMatthew G Knepley PetscScalar *newVals; 4614325cce7SMatthew G Knepley PetscInt *newCols; 4624325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 4634325cce7SMatthew G Knepley PetscErrorCode ierr; 4644325cce7SMatthew G Knepley 4654325cce7SMatthew G Knepley PetscFunctionBegin; 4664325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 4674325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 4684325cce7SMatthew G Knepley PetscInt ncols; 4694325cce7SMatthew G Knepley 4700298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 4714325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 4720298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 4734325cce7SMatthew G Knepley } 4744325cce7SMatthew G Knepley numRows = rEnd - rStart; 475b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 476dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 4774325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 4784325cce7SMatthew G Knepley const PetscScalar *vals; 4794325cce7SMatthew G Knepley const PetscInt *cols; 480fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 4814325cce7SMatthew G Knepley 4824325cce7SMatthew G Knepley if (r < rEnd) { 4834325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 4844325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 4854325cce7SMatthew G Knepley newCols[c] = cols[c]; 4864325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 4874325cce7SMatthew G Knepley } 488fad45679SMatthew G. Knepley newcols = ncols; 4894325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 490fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 4914325cce7SMatthew G Knepley } 4924325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4934325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4944325cce7SMatthew G Knepley } 4954325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 4964325cce7SMatthew G Knepley PetscFunctionReturn(0); 4974325cce7SMatthew G Knepley } 498