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; 25646531bbSStefano Zampini PetscInt M1,M2,N1,N2; 26c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 270ff8bee4SStefano Zampini PetscBool sametype; 286f79c3a4SBarry Smith 293a40ed3dSBarry Smith PetscFunctionBegin; 300700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 31c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 32646531bbSStefano Zampini PetscValidHeaderSpecific(X,MAT_CLASSID,3); 33646531bbSStefano Zampini PetscCheckSameComm(Y,1,X,3); 34646531bbSStefano Zampini ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr); 35646531bbSStefano Zampini ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr); 36646531bbSStefano Zampini ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr); 37646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr); 38646531bbSStefano Zampini if (M1 != M2 || N1 != N2) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_SIZ,"Non conforming matrix add: global sizes %D x %D, %D x %D",M1,M2,N1,N2); 39646531bbSStefano Zampini if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: local sizes %D x %D, %D x %D",m1,m2,n1,n2); 401987afe7SBarry Smith 410ff8bee4SStefano Zampini ierr = PetscStrcmp(((PetscObject)X)->type_name,((PetscObject)Y)->type_name,&sametype);CHKERRQ(ierr); 42e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 430ff8bee4SStefano Zampini if (Y->ops->axpy && sametype) { 44f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 45d4bb536fSBarry Smith } else { 46646531bbSStefano Zampini if (str != DIFFERENT_NONZERO_PATTERN) { 47f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 48646531bbSStefano Zampini } else { 49646531bbSStefano Zampini Mat B; 50646531bbSStefano Zampini 51646531bbSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 52646531bbSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 53646531bbSStefano Zampini ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 54646531bbSStefano Zampini } 55607cd303SBarry Smith } 56e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 57fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 58b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 59b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 60c41cb2e2SAlejandro Lamas Daviña } 61d67ff14aSKarl Rupp #endif 62607cd303SBarry Smith PetscFunctionReturn(0); 63607cd303SBarry Smith } 64607cd303SBarry Smith 65646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 66646531bbSStefano Zampini { 67646531bbSStefano Zampini PetscErrorCode ierr; 6827afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 69646531bbSStefano Zampini 70646531bbSStefano Zampini PetscFunctionBegin; 71646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 72646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 73646531bbSStefano Zampini if (preall) { 74646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 75646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 76646531bbSStefano Zampini Mat preallocator; 77646531bbSStefano Zampini PetscInt r,rstart,rend; 78646531bbSStefano Zampini PetscInt m,n,M,N; 79646531bbSStefano Zampini 80*6eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 81*6eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 82646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 83646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 84646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 85646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 86646531bbSStefano Zampini ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 87646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 88646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 89646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 90646531bbSStefano Zampini PetscInt ncols; 91646531bbSStefano Zampini const PetscInt *row; 92646531bbSStefano Zampini const PetscScalar *vals; 93646531bbSStefano Zampini 94646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 95646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 96646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 97646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 98646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 99646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 100646531bbSStefano Zampini } 101646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 103*6eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 104*6eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 105646531bbSStefano Zampini 106646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 107646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 108646531bbSStefano Zampini ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr); 109646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 110646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 111646531bbSStefano Zampini } 112646531bbSStefano Zampini PetscFunctionReturn(0); 113646531bbSStefano Zampini } 114646531bbSStefano Zampini 115f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 116607cd303SBarry Smith { 11738baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 1186849ba73SBarry Smith PetscErrorCode ierr; 11938baddfdSBarry Smith const PetscInt *row; 120b3cc6726SBarry Smith PetscScalar *val; 121b3cc6726SBarry Smith const PetscScalar *vals; 122607cd303SBarry Smith 123607cd303SBarry Smith PetscFunctionBegin; 1248dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 12590f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 126*6eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 127f4df32b1SMatthew Knepley if (a == 1.0) { 128d4bb536fSBarry Smith for (i = start; i < end; i++) { 129d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 130d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 131d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 132d4bb536fSBarry Smith } 133d4bb536fSBarry Smith } else { 13421a3365eSStefano Zampini PetscInt vs = 100; 13521a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 13621a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 13706be10caSBarry Smith for (i=start; i<end; i++) { 138b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 13921a3365eSStefano Zampini if (vs < ncols) { 14021a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 14121a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 1426f79c3a4SBarry Smith } 14321a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 144b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 145b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 1466f79c3a4SBarry Smith } 147b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 148d4bb536fSBarry Smith } 149*6eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 1506d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1516d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1523a40ed3dSBarry Smith PetscFunctionReturn(0); 1536f79c3a4SBarry Smith } 154052efed2SBarry Smith 155ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 156ec7775f6SShri Abhyankar { 157ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 158ec7775f6SShri Abhyankar PetscErrorCode ierr; 159ec7775f6SShri Abhyankar const PetscInt *row; 160ec7775f6SShri Abhyankar PetscScalar *val; 161ec7775f6SShri Abhyankar const PetscScalar *vals; 162ec7775f6SShri Abhyankar 163ec7775f6SShri Abhyankar PetscFunctionBegin; 164ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 165ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 166*6eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 167*6eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 168ec7775f6SShri Abhyankar if (a == 1.0) { 169ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 170ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 171ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 172ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 173ec7775f6SShri Abhyankar 174ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 175ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 176ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 177ec7775f6SShri Abhyankar } 178ec7775f6SShri Abhyankar } else { 17921a3365eSStefano Zampini PetscInt vs = 100; 18021a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 18121a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 182ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 183ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 184ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 185ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 186ec7775f6SShri Abhyankar 187ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 18821a3365eSStefano Zampini if (vs < ncols) { 18921a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 19021a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 191ec7775f6SShri Abhyankar } 19221a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 193ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 194ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 195ec7775f6SShri Abhyankar } 196ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 197ec7775f6SShri Abhyankar } 198*6eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 199*6eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 200ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 201ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 202ec7775f6SShri Abhyankar PetscFunctionReturn(0); 203ec7775f6SShri Abhyankar } 204ec7775f6SShri Abhyankar 205052efed2SBarry Smith /*@ 20687828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 207052efed2SBarry Smith 2083f9fe445SBarry Smith Neighbor-wise Collective on Mat 209fee21e36SBarry Smith 21098a79cdbSBarry Smith Input Parameters: 21198a79cdbSBarry Smith + Y - the matrices 21287828ca2SBarry Smith - a - the PetscScalar 21398a79cdbSBarry Smith 2142860a424SLois Curfman McInnes Level: intermediate 2152860a424SLois Curfman McInnes 21695452b02SPatrick Sanan Notes: 21795452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2186f33a894SBarry 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 2196f33a894SBarry Smith entry). 2206f33a894SBarry Smith 2210c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2220c0fd78eSBarry Smith 2236f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 2246f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 2256f33a894SBarry Smith 226052efed2SBarry Smith .keywords: matrix, add, shift 2276b9ee512SLois Curfman McInnes 2280c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 229052efed2SBarry Smith @*/ 2307087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 231052efed2SBarry Smith { 2326849ba73SBarry Smith PetscErrorCode ierr; 233052efed2SBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 2350700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 236ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 237ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2384994cf47SJed Brown MatCheckPreallocated(Y,1); 239b50b34bdSBarry Smith 2401c738b47SStefano Zampini if (Y->ops->shift) { 241f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 2421c738b47SStefano Zampini } else { 2431c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 2441c738b47SStefano Zampini } 2457d68702bSBarry Smith 2465528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 247fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 248b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 249b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 250c41cb2e2SAlejandro Lamas Daviña } 251d67ff14aSKarl Rupp #endif 2523a40ed3dSBarry Smith PetscFunctionReturn(0); 253052efed2SBarry Smith } 2546d84be18SBarry Smith 2557087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 25609f38230SBarry Smith { 25709f38230SBarry Smith PetscErrorCode ierr; 25867576b8bSBarry Smith PetscInt i,start,end; 25909f38230SBarry Smith PetscScalar *v; 26009f38230SBarry Smith 26109f38230SBarry Smith PetscFunctionBegin; 26209f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 26309f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 26409f38230SBarry Smith for (i=start; i<end; i++) { 26509f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 26609f38230SBarry Smith } 26709f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 26809f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 26909f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 27009f38230SBarry Smith PetscFunctionReturn(0); 27109f38230SBarry Smith } 27209f38230SBarry Smith 2736d84be18SBarry Smith /*@ 274f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 275f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 276f56f2b3fSBarry Smith INSERT_VALUES. 2776d84be18SBarry Smith 2786d84be18SBarry Smith Input Parameters: 27998a79cdbSBarry Smith + Y - the input matrix 280f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 281f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 2826d84be18SBarry Smith 2833f9fe445SBarry Smith Neighbor-wise Collective on Mat and Vec 284fee21e36SBarry Smith 28595452b02SPatrick Sanan Notes: 28695452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2876f33a894SBarry 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 2886f33a894SBarry Smith entry). 2896f33a894SBarry Smith 2902860a424SLois Curfman McInnes Level: intermediate 2912860a424SLois Curfman McInnes 2926b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal 2936b9ee512SLois Curfman McInnes 2940c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 2956d84be18SBarry Smith @*/ 2967087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 2976d84be18SBarry Smith { 2986849ba73SBarry Smith PetscErrorCode ierr; 29967576b8bSBarry Smith PetscInt matlocal,veclocal; 3006d84be18SBarry Smith 3013a40ed3dSBarry Smith PetscFunctionBegin; 3020700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3030700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 30467576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 30567576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 30667576b8bSBarry 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); 307f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 308f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 30994d884c6SBarry Smith } else { 31009f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3116d84be18SBarry Smith } 3125528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 3146d84be18SBarry Smith } 315d4bb536fSBarry Smith 316d4bb536fSBarry Smith /*@ 31704aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 318d4bb536fSBarry Smith 3193f9fe445SBarry Smith Logically on Mat 320fee21e36SBarry Smith 32198a79cdbSBarry Smith Input Parameters: 32204aac2b0SHong Zhang + a - the PetscScalar multiplier 32304aac2b0SHong Zhang . Y - the first matrix 32404aac2b0SHong Zhang . X - the second matrix 32504aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 326d4bb536fSBarry Smith 3272860a424SLois Curfman McInnes Level: intermediate 3282860a424SLois Curfman McInnes 329d4bb536fSBarry Smith .keywords: matrix, add 330d4bb536fSBarry Smith 3312860a424SLois Curfman McInnes .seealso: MatAXPY() 332d4bb536fSBarry Smith @*/ 3337087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 334d4bb536fSBarry Smith { 33587828ca2SBarry Smith PetscScalar one = 1.0; 3366849ba73SBarry Smith PetscErrorCode ierr; 33738baddfdSBarry Smith PetscInt mX,mY,nX,nY; 338d4bb536fSBarry Smith 3393a40ed3dSBarry Smith PetscFunctionBegin; 340c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 3410700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 342c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 343329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 344329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 345e32f2f54SBarry 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); 346d4bb536fSBarry Smith 347f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 348cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 3493a40ed3dSBarry Smith PetscFunctionReturn(0); 350d4bb536fSBarry Smith } 351b0a32e0cSBarry Smith 352b0a32e0cSBarry Smith /*@ 353b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 354b0a32e0cSBarry Smith 355b0a32e0cSBarry Smith Collective on Mat 356b0a32e0cSBarry Smith 357b0a32e0cSBarry Smith Input Parameter: 358b0a32e0cSBarry Smith . inmat - the matrix 359b0a32e0cSBarry Smith 360b0a32e0cSBarry Smith Output Parameter: 361f3b1f45cSBarry Smith . mat - the explict operator 362b0a32e0cSBarry Smith 363b0a32e0cSBarry Smith Notes: 364b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 365b0a32e0cSBarry Smith identity matrix. 366b0a32e0cSBarry Smith 367b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 368b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 369b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 370b0a32e0cSBarry Smith 371b0a32e0cSBarry Smith Level: advanced 372b0a32e0cSBarry Smith 373b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 374b0a32e0cSBarry Smith @*/ 3757087cfbeSBarry Smith PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 376b0a32e0cSBarry Smith { 377dfbe8321SBarry Smith PetscErrorCode ierr; 378b0a32e0cSBarry Smith MPI_Comm comm; 37938baddfdSBarry Smith PetscMPIInt size; 380b0a32e0cSBarry Smith 381b0a32e0cSBarry Smith PetscFunctionBegin; 3820700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 3834482741eSBarry Smith PetscValidPointer(mat,2); 384b0a32e0cSBarry Smith 385ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 386b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 387b3d09e86SJed Brown ierr = MatConvert_Shell(inmat,size == 1 ? MATSEQDENSE : MATAIJ,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 38824f910e3SHong Zhang PetscFunctionReturn(0); 38924f910e3SHong Zhang } 3904325cce7SMatthew G Knepley 3914325cce7SMatthew G Knepley /*@ 392f3b1f45cSBarry Smith MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of 393f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 394f3b1f45cSBarry Smith 395f3b1f45cSBarry Smith Collective on Mat 396f3b1f45cSBarry Smith 397f3b1f45cSBarry Smith Input Parameter: 398f3b1f45cSBarry Smith . inmat - the matrix 399f3b1f45cSBarry Smith 400f3b1f45cSBarry Smith Output Parameter: 401f3b1f45cSBarry Smith . mat - the explict operator transposed 402f3b1f45cSBarry Smith 403f3b1f45cSBarry Smith Notes: 404f3b1f45cSBarry Smith This computation is done by applying the transpose of the operator to columns of the 405f3b1f45cSBarry Smith identity matrix. 406f3b1f45cSBarry Smith 407f3b1f45cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 408f3b1f45cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 409f3b1f45cSBarry Smith and is recommended for use only with relatively small systems. 410f3b1f45cSBarry Smith 411f3b1f45cSBarry Smith Level: advanced 412f3b1f45cSBarry Smith 413f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator 414f3b1f45cSBarry Smith @*/ 415f3b1f45cSBarry Smith PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 416f3b1f45cSBarry Smith { 417f3b1f45cSBarry Smith Vec in,out; 418f3b1f45cSBarry Smith PetscErrorCode ierr; 419f3b1f45cSBarry Smith PetscInt i,m,n,M,N,*rows,start,end; 420f3b1f45cSBarry Smith MPI_Comm comm; 421f3b1f45cSBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 422f3b1f45cSBarry Smith PetscMPIInt size; 423f3b1f45cSBarry Smith 424f3b1f45cSBarry Smith PetscFunctionBegin; 425f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 426f3b1f45cSBarry Smith PetscValidPointer(mat,2); 427f3b1f45cSBarry Smith 428f3b1f45cSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 429f3b1f45cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 430f3b1f45cSBarry Smith 431f3b1f45cSBarry Smith ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 432f3b1f45cSBarry Smith ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 433f3b1f45cSBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 434f3b1f45cSBarry Smith ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 435f3b1f45cSBarry Smith ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 436f3b1f45cSBarry Smith ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 437f3b1f45cSBarry Smith for (i=0; i<m; i++) rows[i] = start + i; 438f3b1f45cSBarry Smith 439f3b1f45cSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 440f3b1f45cSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 441f3b1f45cSBarry Smith if (size == 1) { 442f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 443f3b1f45cSBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 444f3b1f45cSBarry Smith } else { 445f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 446f3b1f45cSBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 447f3b1f45cSBarry Smith } 448f3b1f45cSBarry Smith 449f3b1f45cSBarry Smith for (i=0; i<N; i++) { 450f3b1f45cSBarry Smith 451f3b1f45cSBarry Smith ierr = VecSet(in,zero);CHKERRQ(ierr); 452f3b1f45cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 453f3b1f45cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 454f3b1f45cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 455f3b1f45cSBarry Smith 456f3b1f45cSBarry Smith ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 457f3b1f45cSBarry Smith 458f3b1f45cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 459f3b1f45cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 460f3b1f45cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 461f3b1f45cSBarry Smith 462f3b1f45cSBarry Smith } 463f3b1f45cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 464f3b1f45cSBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 465f3b1f45cSBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 466f3b1f45cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 467f3b1f45cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 468f3b1f45cSBarry Smith PetscFunctionReturn(0); 469f3b1f45cSBarry Smith } 470f3b1f45cSBarry Smith 471f3b1f45cSBarry Smith /*@ 4724325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4734325cce7SMatthew G Knepley 4744325cce7SMatthew G Knepley Input Parameters: 4754325cce7SMatthew G Knepley + A - The matrix 4764325cce7SMatthew G Knepley - tol - The zero tolerance 4774325cce7SMatthew G Knepley 4784325cce7SMatthew G Knepley Output Parameters: 4794325cce7SMatthew G Knepley . A - The chopped matrix 4804325cce7SMatthew G Knepley 4814325cce7SMatthew G Knepley Level: intermediate 4824325cce7SMatthew G Knepley 4834325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 4843fc99919SSatish Balay @*/ 4854325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 4864325cce7SMatthew G Knepley { 4874325cce7SMatthew G Knepley PetscScalar *newVals; 4884325cce7SMatthew G Knepley PetscInt *newCols; 4894325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 4904325cce7SMatthew G Knepley PetscErrorCode ierr; 4914325cce7SMatthew G Knepley 4924325cce7SMatthew G Knepley PetscFunctionBegin; 4934325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 4944325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 4954325cce7SMatthew G Knepley PetscInt ncols; 4964325cce7SMatthew G Knepley 4970298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 4984325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 4990298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5004325cce7SMatthew G Knepley } 5014325cce7SMatthew G Knepley numRows = rEnd - rStart; 502b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 503dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 5044325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 5054325cce7SMatthew G Knepley const PetscScalar *vals; 5064325cce7SMatthew G Knepley const PetscInt *cols; 507fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 5084325cce7SMatthew G Knepley 5094325cce7SMatthew G Knepley if (r < rEnd) { 5104325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 5114325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 5124325cce7SMatthew G Knepley newCols[c] = cols[c]; 5134325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 5144325cce7SMatthew G Knepley } 515fad45679SMatthew G. Knepley newcols = ncols; 5164325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 517fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 5184325cce7SMatthew G Knepley } 5194325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5204325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5214325cce7SMatthew G Knepley } 5224325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 5234325cce7SMatthew G Knepley PetscFunctionReturn(0); 5244325cce7SMatthew G Knepley } 525