16f79c3a4SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 4d54c938cSPierre Jolivet static PetscErrorCode MatTransposeAXPY_Private(Mat Y,PetscScalar a,Mat X,MatStructure str,Mat T) 5d54c938cSPierre Jolivet { 6d54c938cSPierre Jolivet PetscErrorCode ierr,(*f)(Mat,Mat*); 7d54c938cSPierre Jolivet Mat A,F; 8d54c938cSPierre Jolivet 9d54c938cSPierre Jolivet PetscFunctionBegin; 10d54c938cSPierre Jolivet ierr = PetscObjectQueryFunction((PetscObject)T,"MatTransposeGetMat_C",&f);CHKERRQ(ierr); 11d54c938cSPierre Jolivet if (f) { 12d54c938cSPierre Jolivet ierr = MatTransposeGetMat(T,&A);CHKERRQ(ierr); 13d54c938cSPierre Jolivet if (T == X) { 14d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Explicitly transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 15d54c938cSPierre Jolivet ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 16d54c938cSPierre Jolivet A = Y; 17d54c938cSPierre Jolivet } else { 18d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 19d54c938cSPierre Jolivet ierr = MatTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 20d54c938cSPierre Jolivet } 21d54c938cSPierre Jolivet } else { 22d54c938cSPierre Jolivet ierr = MatHermitianTransposeGetMat(T,&A);CHKERRQ(ierr); 23d54c938cSPierre Jolivet if (T == X) { 24d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Explicitly Hermitian transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 25d54c938cSPierre Jolivet ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 26d54c938cSPierre Jolivet A = Y; 27d54c938cSPierre Jolivet } else { 28d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Hermitian transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 29d54c938cSPierre Jolivet ierr = MatHermitianTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 30d54c938cSPierre Jolivet } 31d54c938cSPierre Jolivet } 32d54c938cSPierre Jolivet ierr = MatAXPY(A,a,F,str);CHKERRQ(ierr); 33d54c938cSPierre Jolivet ierr = MatDestroy(&F);CHKERRQ(ierr); 34d54c938cSPierre Jolivet PetscFunctionReturn(0); 35d54c938cSPierre Jolivet } 36d54c938cSPierre Jolivet 3706be10caSBarry Smith /*@ 3821c89e3eSBarry Smith MatAXPY - Computes Y = a*X + Y. 396f79c3a4SBarry Smith 403f9fe445SBarry Smith Logically Collective on Mat 41fee21e36SBarry Smith 4298a79cdbSBarry Smith Input Parameters: 43607cd303SBarry Smith + a - the scalar multiplier 44607cd303SBarry Smith . X - the first matrix 45607cd303SBarry Smith . Y - the second matrix 46407f6b05SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 47407f6b05SHong Zhang or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 4898a79cdbSBarry Smith 492860a424SLois Curfman McInnes Level: intermediate 502860a424SLois Curfman McInnes 512860a424SLois Curfman McInnes .seealso: MatAYPX() 5206be10caSBarry Smith @*/ 537087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 546f79c3a4SBarry Smith { 55d54c938cSPierre Jolivet PetscErrorCode ierr; 56646531bbSStefano Zampini PetscInt M1,M2,N1,N2; 57c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 58a683ea4aSPierre Jolivet MatType t1,t2; 59a683ea4aSPierre Jolivet PetscBool sametype,transpose; 606f79c3a4SBarry Smith 613a40ed3dSBarry Smith PetscFunctionBegin; 620700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 63c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 64646531bbSStefano Zampini PetscValidHeaderSpecific(X,MAT_CLASSID,3); 65646531bbSStefano Zampini PetscCheckSameComm(Y,1,X,3); 66646531bbSStefano Zampini ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr); 67646531bbSStefano Zampini ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr); 68646531bbSStefano Zampini ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr); 69646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr); 70646531bbSStefano 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); 71646531bbSStefano 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); 72*1ccb7edbSStefano Zampini if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (Y)"); 73*1ccb7edbSStefano Zampini if (!X->assembled) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (X)"); 741987afe7SBarry Smith 75a683ea4aSPierre Jolivet ierr = MatGetType(X,&t1);CHKERRQ(ierr); 76a683ea4aSPierre Jolivet ierr = MatGetType(Y,&t2);CHKERRQ(ierr); 77a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,t2,&sametype);CHKERRQ(ierr); 78e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 790ff8bee4SStefano Zampini if (Y->ops->axpy && sametype) { 80f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 81d4bb536fSBarry Smith } else { 82a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 83a683ea4aSPierre Jolivet if (transpose) { 84d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,X);CHKERRQ(ierr); 85a683ea4aSPierre Jolivet } else { 86a683ea4aSPierre Jolivet ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 87a683ea4aSPierre Jolivet if (transpose) { 88d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,Y);CHKERRQ(ierr); 89a683ea4aSPierre Jolivet } else { 90646531bbSStefano Zampini if (str != DIFFERENT_NONZERO_PATTERN) { 91f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 92646531bbSStefano Zampini } else { 93646531bbSStefano Zampini Mat B; 94646531bbSStefano Zampini 95646531bbSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 96646531bbSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 97646531bbSStefano Zampini ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 98646531bbSStefano Zampini } 99607cd303SBarry Smith } 100a683ea4aSPierre Jolivet } 101a683ea4aSPierre Jolivet } 102e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 103607cd303SBarry Smith PetscFunctionReturn(0); 104607cd303SBarry Smith } 105607cd303SBarry Smith 106646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 107646531bbSStefano Zampini { 108646531bbSStefano Zampini PetscErrorCode ierr; 10927afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 110646531bbSStefano Zampini 111646531bbSStefano Zampini PetscFunctionBegin; 112646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 113646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 114646531bbSStefano Zampini if (preall) { 115646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 116646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 117646531bbSStefano Zampini Mat preallocator; 118646531bbSStefano Zampini PetscInt r,rstart,rend; 119646531bbSStefano Zampini PetscInt m,n,M,N; 120646531bbSStefano Zampini 1216eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 1226eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 123646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 124646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 125646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 126646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 127646531bbSStefano Zampini ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 128646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 129646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 130646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 131646531bbSStefano Zampini PetscInt ncols; 132646531bbSStefano Zampini const PetscInt *row; 133646531bbSStefano Zampini const PetscScalar *vals; 134646531bbSStefano Zampini 135646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 136646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 137646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 138646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 139646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 140646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 141646531bbSStefano Zampini } 142646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1446eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 1456eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 146646531bbSStefano Zampini 147646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 148646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 149646531bbSStefano Zampini ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr); 150646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 151646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 152646531bbSStefano Zampini } 153646531bbSStefano Zampini PetscFunctionReturn(0); 154646531bbSStefano Zampini } 155646531bbSStefano Zampini 156f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 157607cd303SBarry Smith { 15838baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 1596849ba73SBarry Smith PetscErrorCode ierr; 16038baddfdSBarry Smith const PetscInt *row; 161b3cc6726SBarry Smith PetscScalar *val; 162b3cc6726SBarry Smith const PetscScalar *vals; 163607cd303SBarry Smith 164607cd303SBarry Smith PetscFunctionBegin; 1658dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 16690f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 1676eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 168f4df32b1SMatthew Knepley if (a == 1.0) { 169d4bb536fSBarry Smith for (i = start; i < end; i++) { 170d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 171d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 172d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 173d4bb536fSBarry Smith } 174d4bb536fSBarry Smith } else { 17521a3365eSStefano Zampini PetscInt vs = 100; 17621a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 17721a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 17806be10caSBarry Smith for (i=start; i<end; i++) { 179b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 18021a3365eSStefano Zampini if (vs < ncols) { 18121a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 18221a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 1836f79c3a4SBarry Smith } 18421a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 185b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 186b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 1876f79c3a4SBarry Smith } 188b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 189d4bb536fSBarry Smith } 1906eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 1916d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1926d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 1946f79c3a4SBarry Smith } 195052efed2SBarry Smith 196ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 197ec7775f6SShri Abhyankar { 198ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 199ec7775f6SShri Abhyankar PetscErrorCode ierr; 200ec7775f6SShri Abhyankar const PetscInt *row; 201ec7775f6SShri Abhyankar PetscScalar *val; 202ec7775f6SShri Abhyankar const PetscScalar *vals; 203ec7775f6SShri Abhyankar 204ec7775f6SShri Abhyankar PetscFunctionBegin; 205ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 206ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 2076eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 2086eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 209ec7775f6SShri Abhyankar if (a == 1.0) { 210ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 211ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 212ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 213ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 214ec7775f6SShri Abhyankar 215ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 216ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 217ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 218ec7775f6SShri Abhyankar } 219ec7775f6SShri Abhyankar } else { 22021a3365eSStefano Zampini PetscInt vs = 100; 22121a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 22221a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 223ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 224ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 225ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 226ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 227ec7775f6SShri Abhyankar 228ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 22921a3365eSStefano Zampini if (vs < ncols) { 23021a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 23121a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 232ec7775f6SShri Abhyankar } 23321a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 234ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 235ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 236ec7775f6SShri Abhyankar } 237ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 238ec7775f6SShri Abhyankar } 2396eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 2406eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 241ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 242ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243ec7775f6SShri Abhyankar PetscFunctionReturn(0); 244ec7775f6SShri Abhyankar } 245ec7775f6SShri Abhyankar 246052efed2SBarry Smith /*@ 24787828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 248052efed2SBarry Smith 2493f9fe445SBarry Smith Neighbor-wise Collective on Mat 250fee21e36SBarry Smith 25198a79cdbSBarry Smith Input Parameters: 25298a79cdbSBarry Smith + Y - the matrices 25387828ca2SBarry Smith - a - the PetscScalar 25498a79cdbSBarry Smith 2552860a424SLois Curfman McInnes Level: intermediate 2562860a424SLois Curfman McInnes 25795452b02SPatrick Sanan Notes: 25895452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2596f33a894SBarry 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 2606f33a894SBarry Smith entry). 2616f33a894SBarry Smith 2620c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2630c0fd78eSBarry Smith 2646f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 2656f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 2666f33a894SBarry Smith 2670c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 268052efed2SBarry Smith @*/ 2697087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 270052efed2SBarry Smith { 2716849ba73SBarry Smith PetscErrorCode ierr; 272052efed2SBarry Smith 2733a40ed3dSBarry Smith PetscFunctionBegin; 2740700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 275ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 276ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2774994cf47SJed Brown MatCheckPreallocated(Y,1); 2781b71c1a7SBarry Smith if (a == 0.0) PetscFunctionReturn(0); 279b50b34bdSBarry Smith 2801c738b47SStefano Zampini if (Y->ops->shift) { 281f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 2821c738b47SStefano Zampini } else { 2831c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 2841c738b47SStefano Zampini } 2857d68702bSBarry Smith 2865528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 2873a40ed3dSBarry Smith PetscFunctionReturn(0); 288052efed2SBarry Smith } 2896d84be18SBarry Smith 2907087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 29109f38230SBarry Smith { 29209f38230SBarry Smith PetscErrorCode ierr; 29367576b8bSBarry Smith PetscInt i,start,end; 29409f38230SBarry Smith PetscScalar *v; 29509f38230SBarry Smith 29609f38230SBarry Smith PetscFunctionBegin; 29709f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 29809f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 29909f38230SBarry Smith for (i=start; i<end; i++) { 30009f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 30109f38230SBarry Smith } 30209f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 30309f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30409f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30509f38230SBarry Smith PetscFunctionReturn(0); 30609f38230SBarry Smith } 30709f38230SBarry Smith 3086d84be18SBarry Smith /*@ 309f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 310f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 311f56f2b3fSBarry Smith INSERT_VALUES. 3126d84be18SBarry Smith 3136d84be18SBarry Smith Input Parameters: 31498a79cdbSBarry Smith + Y - the input matrix 315f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 316f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3176d84be18SBarry Smith 318d083f849SBarry Smith Neighbor-wise Collective on Mat 319fee21e36SBarry Smith 32095452b02SPatrick Sanan Notes: 32195452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3226f33a894SBarry 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 3236f33a894SBarry Smith entry). 3246f33a894SBarry Smith 3252860a424SLois Curfman McInnes Level: intermediate 3262860a424SLois Curfman McInnes 3270c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 3286d84be18SBarry Smith @*/ 3297087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 3306d84be18SBarry Smith { 3316849ba73SBarry Smith PetscErrorCode ierr; 33267576b8bSBarry Smith PetscInt matlocal,veclocal; 3336d84be18SBarry Smith 3343a40ed3dSBarry Smith PetscFunctionBegin; 3350700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3360700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 33767576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 33867576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 33967576b8bSBarry 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); 340f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 341f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 34294d884c6SBarry Smith } else { 34309f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3446d84be18SBarry Smith } 3455528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3463a40ed3dSBarry Smith PetscFunctionReturn(0); 3476d84be18SBarry Smith } 348d4bb536fSBarry Smith 349d4bb536fSBarry Smith /*@ 35004aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 351d4bb536fSBarry Smith 3523f9fe445SBarry Smith Logically on Mat 353fee21e36SBarry Smith 35498a79cdbSBarry Smith Input Parameters: 35504aac2b0SHong Zhang + a - the PetscScalar multiplier 35604aac2b0SHong Zhang . Y - the first matrix 35704aac2b0SHong Zhang . X - the second matrix 35804aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 359d4bb536fSBarry Smith 3602860a424SLois Curfman McInnes Level: intermediate 3612860a424SLois Curfman McInnes 3622860a424SLois Curfman McInnes .seealso: MatAXPY() 363d4bb536fSBarry Smith @*/ 3647087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 365d4bb536fSBarry Smith { 36687828ca2SBarry Smith PetscScalar one = 1.0; 3676849ba73SBarry Smith PetscErrorCode ierr; 36838baddfdSBarry Smith PetscInt mX,mY,nX,nY; 369d4bb536fSBarry Smith 3703a40ed3dSBarry Smith PetscFunctionBegin; 371c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 3720700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 373c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 374329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 375329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 376e32f2f54SBarry 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); 377f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 378cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 3793a40ed3dSBarry Smith PetscFunctionReturn(0); 380d4bb536fSBarry Smith } 381b0a32e0cSBarry Smith 382b0a32e0cSBarry Smith /*@ 3830bacdadaSStefano Zampini MatComputeOperator - Computes the explicit matrix 384b0a32e0cSBarry Smith 385b0a32e0cSBarry Smith Collective on Mat 386b0a32e0cSBarry Smith 387b0a32e0cSBarry Smith Input Parameter: 388186a3e20SStefano Zampini + inmat - the matrix 389186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator 390b0a32e0cSBarry Smith 391b0a32e0cSBarry Smith Output Parameter: 392f3b1f45cSBarry Smith . mat - the explict operator 393b0a32e0cSBarry Smith 394b0a32e0cSBarry Smith Notes: 395186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 396186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 397186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 398b0a32e0cSBarry Smith 399b0a32e0cSBarry Smith Level: advanced 400b0a32e0cSBarry Smith 401b0a32e0cSBarry Smith @*/ 4020bacdadaSStefano Zampini PetscErrorCode MatComputeOperator(Mat inmat,MatType mattype,Mat *mat) 403b0a32e0cSBarry Smith { 404dfbe8321SBarry Smith PetscErrorCode ierr; 405b0a32e0cSBarry Smith 406b0a32e0cSBarry Smith PetscFunctionBegin; 4070700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 408186a3e20SStefano Zampini PetscValidPointer(mat,3); 409186a3e20SStefano Zampini ierr = MatConvert_Shell(inmat,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 41024f910e3SHong Zhang PetscFunctionReturn(0); 41124f910e3SHong Zhang } 4124325cce7SMatthew G Knepley 4134325cce7SMatthew G Knepley /*@ 4140bacdadaSStefano Zampini MatComputeOperatorTranspose - Computes the explicit matrix representation of 415f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 416f3b1f45cSBarry Smith 417f3b1f45cSBarry Smith Collective on Mat 418f3b1f45cSBarry Smith 419f3b1f45cSBarry Smith Input Parameter: 420f3b1f45cSBarry Smith . inmat - the matrix 421f3b1f45cSBarry Smith 422f3b1f45cSBarry Smith Output Parameter: 423f3b1f45cSBarry Smith . mat - the explict operator transposed 424f3b1f45cSBarry Smith 425f3b1f45cSBarry Smith Notes: 426186a3e20SStefano Zampini This computation is done by applying the transpose of the operator to columns of the identity matrix. 427186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 428186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 429f3b1f45cSBarry Smith 430f3b1f45cSBarry Smith Level: advanced 431f3b1f45cSBarry Smith 432f3b1f45cSBarry Smith @*/ 4330bacdadaSStefano Zampini PetscErrorCode MatComputeOperatorTranspose(Mat inmat,MatType mattype,Mat *mat) 434f3b1f45cSBarry Smith { 435186a3e20SStefano Zampini Mat A; 436f3b1f45cSBarry Smith PetscErrorCode ierr; 437f3b1f45cSBarry Smith 438f3b1f45cSBarry Smith PetscFunctionBegin; 439f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 440186a3e20SStefano Zampini PetscValidPointer(mat,3); 441186a3e20SStefano Zampini ierr = MatCreateTranspose(inmat,&A);CHKERRQ(ierr); 442186a3e20SStefano Zampini ierr = MatConvert_Shell(A,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 443186a3e20SStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 444f3b1f45cSBarry Smith PetscFunctionReturn(0); 445f3b1f45cSBarry Smith } 446f3b1f45cSBarry Smith 447f3b1f45cSBarry Smith /*@ 4484325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4494325cce7SMatthew G Knepley 4504325cce7SMatthew G Knepley Input Parameters: 4514325cce7SMatthew G Knepley + A - The matrix 4524325cce7SMatthew G Knepley - tol - The zero tolerance 4534325cce7SMatthew G Knepley 4544325cce7SMatthew G Knepley Output Parameters: 4554325cce7SMatthew G Knepley . A - The chopped matrix 4564325cce7SMatthew G Knepley 4574325cce7SMatthew G Knepley Level: intermediate 4584325cce7SMatthew G Knepley 4594325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 4603fc99919SSatish Balay @*/ 4614325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 4624325cce7SMatthew G Knepley { 4634325cce7SMatthew G Knepley PetscScalar *newVals; 4644325cce7SMatthew G Knepley PetscInt *newCols; 4654325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 4664325cce7SMatthew G Knepley PetscErrorCode ierr; 4674325cce7SMatthew G Knepley 4684325cce7SMatthew G Knepley PetscFunctionBegin; 4694325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 4704325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 4714325cce7SMatthew G Knepley PetscInt ncols; 4724325cce7SMatthew G Knepley 4730298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 4744325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 4750298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 4764325cce7SMatthew G Knepley } 4774325cce7SMatthew G Knepley numRows = rEnd - rStart; 478b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 479dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 4804325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 4814325cce7SMatthew G Knepley const PetscScalar *vals; 4824325cce7SMatthew G Knepley const PetscInt *cols; 483fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 4844325cce7SMatthew G Knepley 4854325cce7SMatthew G Knepley if (r < rEnd) { 4864325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 4874325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 4884325cce7SMatthew G Knepley newCols[c] = cols[c]; 4894325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 4904325cce7SMatthew G Knepley } 491fad45679SMatthew G. Knepley newcols = ncols; 4924325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 493fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 4944325cce7SMatthew G Knepley } 4954325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4964325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4974325cce7SMatthew G Knepley } 4984325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 4994325cce7SMatthew G Knepley PetscFunctionReturn(0); 5004325cce7SMatthew G Knepley } 501