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 { 24*a683ea4aSPierre Jolivet PetscErrorCode ierr,(*f)(Mat,Mat*); 25646531bbSStefano Zampini PetscInt M1,M2,N1,N2; 26c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 27*a683ea4aSPierre Jolivet MatType t1,t2; 28*a683ea4aSPierre Jolivet Mat T,F; 29*a683ea4aSPierre Jolivet PetscBool sametype,transpose; 306f79c3a4SBarry Smith 313a40ed3dSBarry Smith PetscFunctionBegin; 320700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 33c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 34646531bbSStefano Zampini PetscValidHeaderSpecific(X,MAT_CLASSID,3); 35646531bbSStefano Zampini PetscCheckSameComm(Y,1,X,3); 36646531bbSStefano Zampini ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr); 37646531bbSStefano Zampini ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr); 38646531bbSStefano Zampini ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr); 39646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr); 40646531bbSStefano 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); 41646531bbSStefano 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); 421987afe7SBarry Smith 43*a683ea4aSPierre Jolivet ierr = MatGetType(X,&t1);CHKERRQ(ierr); 44*a683ea4aSPierre Jolivet ierr = MatGetType(Y,&t2);CHKERRQ(ierr); 45*a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,t2,&sametype);CHKERRQ(ierr); 46e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 470ff8bee4SStefano Zampini if (Y->ops->axpy && sametype) { 48f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 49d4bb536fSBarry Smith } else { 50*a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 51*a683ea4aSPierre Jolivet if (transpose) { 52*a683ea4aSPierre Jolivet PetscObjectQueryFunction((PetscObject)X,"MatTransposeGetMat_C",&f); 53*a683ea4aSPierre Jolivet if (f) { 54*a683ea4aSPierre Jolivet ierr = f(X,&T);CHKERRQ(ierr); 55*a683ea4aSPierre Jolivet ierr = MatTranspose(T,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 56*a683ea4aSPierre Jolivet } else { 57*a683ea4aSPierre Jolivet ierr = MatHermitianTransposeGetMat(X,&T);CHKERRQ(ierr); 58*a683ea4aSPierre Jolivet ierr = MatHermitianTranspose(T,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 59*a683ea4aSPierre Jolivet } 60*a683ea4aSPierre Jolivet ierr = MatAXPY(Y,a,F,str);CHKERRQ(ierr); 61*a683ea4aSPierre Jolivet ierr = MatDestroy(&F);CHKERRQ(ierr); 62*a683ea4aSPierre Jolivet } else { 63*a683ea4aSPierre Jolivet ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 64*a683ea4aSPierre Jolivet if (transpose) { 65*a683ea4aSPierre Jolivet PetscObjectQueryFunction((PetscObject)Y,"MatTransposeGetMat_C",&f); 66*a683ea4aSPierre Jolivet if (f) { 67*a683ea4aSPierre Jolivet ierr = f(Y,&T);CHKERRQ(ierr); 68*a683ea4aSPierre Jolivet ierr = MatTranspose(T,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 69*a683ea4aSPierre Jolivet } else { 70*a683ea4aSPierre Jolivet ierr = MatHermitianTransposeGetMat(Y,&T);CHKERRQ(ierr); 71*a683ea4aSPierre Jolivet ierr = MatHermitianTranspose(T,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 72*a683ea4aSPierre Jolivet } 73*a683ea4aSPierre Jolivet ierr = MatAXPY(F,a,X,str);CHKERRQ(ierr); 74*a683ea4aSPierre Jolivet ierr = MatDestroy(&F);CHKERRQ(ierr); 75*a683ea4aSPierre Jolivet } else { 76646531bbSStefano Zampini if (str != DIFFERENT_NONZERO_PATTERN) { 77f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 78646531bbSStefano Zampini } else { 79646531bbSStefano Zampini Mat B; 80646531bbSStefano Zampini 81646531bbSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 82646531bbSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 83646531bbSStefano Zampini ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 84646531bbSStefano Zampini } 85607cd303SBarry Smith } 86*a683ea4aSPierre Jolivet } 87*a683ea4aSPierre Jolivet } 88e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 89fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 90b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 91b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 92c41cb2e2SAlejandro Lamas Daviña } 93d67ff14aSKarl Rupp #endif 94607cd303SBarry Smith PetscFunctionReturn(0); 95607cd303SBarry Smith } 96607cd303SBarry Smith 97646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 98646531bbSStefano Zampini { 99646531bbSStefano Zampini PetscErrorCode ierr; 10027afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 101646531bbSStefano Zampini 102646531bbSStefano Zampini PetscFunctionBegin; 103646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 104646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 105646531bbSStefano Zampini if (preall) { 106646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 107646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 108646531bbSStefano Zampini Mat preallocator; 109646531bbSStefano Zampini PetscInt r,rstart,rend; 110646531bbSStefano Zampini PetscInt m,n,M,N; 111646531bbSStefano Zampini 112646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 113646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 114646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 115646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 116646531bbSStefano Zampini ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 117646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 118646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 119646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 120646531bbSStefano Zampini PetscInt ncols; 121646531bbSStefano Zampini const PetscInt *row; 122646531bbSStefano Zampini const PetscScalar *vals; 123646531bbSStefano Zampini 124646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 125646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 126646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 127646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 128646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 129646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 130646531bbSStefano Zampini } 131646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133646531bbSStefano Zampini 134646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 135646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 136646531bbSStefano Zampini ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr); 137646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 138646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 139646531bbSStefano Zampini } 140646531bbSStefano Zampini PetscFunctionReturn(0); 141646531bbSStefano Zampini } 142646531bbSStefano Zampini 143f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 144607cd303SBarry Smith { 14538baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 1466849ba73SBarry Smith PetscErrorCode ierr; 14738baddfdSBarry Smith const PetscInt *row; 148b3cc6726SBarry Smith PetscScalar *val; 149b3cc6726SBarry Smith const PetscScalar *vals; 150607cd303SBarry Smith 151607cd303SBarry Smith PetscFunctionBegin; 1528dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 15390f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 154f4df32b1SMatthew Knepley if (a == 1.0) { 155d4bb536fSBarry Smith for (i = start; i < end; i++) { 156d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 157d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 158d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 159d4bb536fSBarry Smith } 160d4bb536fSBarry Smith } else { 16121a3365eSStefano Zampini PetscInt vs = 100; 16221a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 16321a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 16406be10caSBarry Smith for (i=start; i<end; i++) { 165b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 16621a3365eSStefano Zampini if (vs < ncols) { 16721a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 16821a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 1696f79c3a4SBarry Smith } 17021a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 171b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 172b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 1736f79c3a4SBarry Smith } 174b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 175d4bb536fSBarry Smith } 1766d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1776d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1783a40ed3dSBarry Smith PetscFunctionReturn(0); 1796f79c3a4SBarry Smith } 180052efed2SBarry Smith 181ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 182ec7775f6SShri Abhyankar { 183ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 184ec7775f6SShri Abhyankar PetscErrorCode ierr; 185ec7775f6SShri Abhyankar const PetscInt *row; 186ec7775f6SShri Abhyankar PetscScalar *val; 187ec7775f6SShri Abhyankar const PetscScalar *vals; 188ec7775f6SShri Abhyankar 189ec7775f6SShri Abhyankar PetscFunctionBegin; 190ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 191ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 192ec7775f6SShri Abhyankar if (a == 1.0) { 193ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 194ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 195ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 196ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 197ec7775f6SShri Abhyankar 198ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 199ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 200ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 201ec7775f6SShri Abhyankar } 202ec7775f6SShri Abhyankar } else { 20321a3365eSStefano Zampini PetscInt vs = 100; 20421a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 20521a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 206ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 207ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 208ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 209ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 210ec7775f6SShri Abhyankar 211ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 21221a3365eSStefano Zampini if (vs < ncols) { 21321a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 21421a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 215ec7775f6SShri Abhyankar } 21621a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 217ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 218ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 219ec7775f6SShri Abhyankar } 220ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 221ec7775f6SShri Abhyankar } 222ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224ec7775f6SShri Abhyankar PetscFunctionReturn(0); 225ec7775f6SShri Abhyankar } 226ec7775f6SShri Abhyankar 227052efed2SBarry Smith /*@ 22887828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 229052efed2SBarry Smith 2303f9fe445SBarry Smith Neighbor-wise Collective on Mat 231fee21e36SBarry Smith 23298a79cdbSBarry Smith Input Parameters: 23398a79cdbSBarry Smith + Y - the matrices 23487828ca2SBarry Smith - a - the PetscScalar 23598a79cdbSBarry Smith 2362860a424SLois Curfman McInnes Level: intermediate 2372860a424SLois Curfman McInnes 23895452b02SPatrick Sanan Notes: 23995452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2406f33a894SBarry 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 2416f33a894SBarry Smith entry). 2426f33a894SBarry Smith 2430c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2440c0fd78eSBarry Smith 2456f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 2466f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 2476f33a894SBarry Smith 248052efed2SBarry Smith .keywords: matrix, add, shift 2496b9ee512SLois Curfman McInnes 2500c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 251052efed2SBarry Smith @*/ 2527087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 253052efed2SBarry Smith { 2546849ba73SBarry Smith PetscErrorCode ierr; 255052efed2SBarry Smith 2563a40ed3dSBarry Smith PetscFunctionBegin; 2570700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 258ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 259ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2604994cf47SJed Brown MatCheckPreallocated(Y,1); 261b50b34bdSBarry Smith 2621c738b47SStefano Zampini if (Y->ops->shift) { 263f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 2641c738b47SStefano Zampini } else { 2651c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 2661c738b47SStefano Zampini } 2677d68702bSBarry Smith 2685528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 269fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 270b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 271b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 272c41cb2e2SAlejandro Lamas Daviña } 273d67ff14aSKarl Rupp #endif 2743a40ed3dSBarry Smith PetscFunctionReturn(0); 275052efed2SBarry Smith } 2766d84be18SBarry Smith 2777087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 27809f38230SBarry Smith { 27909f38230SBarry Smith PetscErrorCode ierr; 28067576b8bSBarry Smith PetscInt i,start,end; 28109f38230SBarry Smith PetscScalar *v; 28209f38230SBarry Smith 28309f38230SBarry Smith PetscFunctionBegin; 28409f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 28509f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 28609f38230SBarry Smith for (i=start; i<end; i++) { 28709f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 28809f38230SBarry Smith } 28909f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 29009f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 29109f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 29209f38230SBarry Smith PetscFunctionReturn(0); 29309f38230SBarry Smith } 29409f38230SBarry Smith 2956d84be18SBarry Smith /*@ 296f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 297f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 298f56f2b3fSBarry Smith INSERT_VALUES. 2996d84be18SBarry Smith 3006d84be18SBarry Smith Input Parameters: 30198a79cdbSBarry Smith + Y - the input matrix 302f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 303f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3046d84be18SBarry Smith 3053f9fe445SBarry Smith Neighbor-wise Collective on Mat and Vec 306fee21e36SBarry Smith 30795452b02SPatrick Sanan Notes: 30895452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3096f33a894SBarry 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 3106f33a894SBarry Smith entry). 3116f33a894SBarry Smith 3122860a424SLois Curfman McInnes Level: intermediate 3132860a424SLois Curfman McInnes 3146b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal 3156b9ee512SLois Curfman McInnes 3160c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 3176d84be18SBarry Smith @*/ 3187087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 3196d84be18SBarry Smith { 3206849ba73SBarry Smith PetscErrorCode ierr; 32167576b8bSBarry Smith PetscInt matlocal,veclocal; 3226d84be18SBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 3240700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3250700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 32667576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 32767576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 32867576b8bSBarry 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); 329f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 330f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 33194d884c6SBarry Smith } else { 33209f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3336d84be18SBarry Smith } 3345528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 3366d84be18SBarry Smith } 337d4bb536fSBarry Smith 338d4bb536fSBarry Smith /*@ 33904aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 340d4bb536fSBarry Smith 3413f9fe445SBarry Smith Logically on Mat 342fee21e36SBarry Smith 34398a79cdbSBarry Smith Input Parameters: 34404aac2b0SHong Zhang + a - the PetscScalar multiplier 34504aac2b0SHong Zhang . Y - the first matrix 34604aac2b0SHong Zhang . X - the second matrix 34704aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 348d4bb536fSBarry Smith 3492860a424SLois Curfman McInnes Level: intermediate 3502860a424SLois Curfman McInnes 351d4bb536fSBarry Smith .keywords: matrix, add 352d4bb536fSBarry Smith 3532860a424SLois Curfman McInnes .seealso: MatAXPY() 354d4bb536fSBarry Smith @*/ 3557087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 356d4bb536fSBarry Smith { 35787828ca2SBarry Smith PetscScalar one = 1.0; 3586849ba73SBarry Smith PetscErrorCode ierr; 35938baddfdSBarry Smith PetscInt mX,mY,nX,nY; 360d4bb536fSBarry Smith 3613a40ed3dSBarry Smith PetscFunctionBegin; 362c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 3630700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 364c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 365329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 366329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 367e32f2f54SBarry 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); 368d4bb536fSBarry Smith 369f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 370cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 3713a40ed3dSBarry Smith PetscFunctionReturn(0); 372d4bb536fSBarry Smith } 373b0a32e0cSBarry Smith 374b0a32e0cSBarry Smith /*@ 375b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 376b0a32e0cSBarry Smith 377b0a32e0cSBarry Smith Collective on Mat 378b0a32e0cSBarry Smith 379b0a32e0cSBarry Smith Input Parameter: 380b0a32e0cSBarry Smith . inmat - the matrix 381b0a32e0cSBarry Smith 382b0a32e0cSBarry Smith Output Parameter: 383f3b1f45cSBarry Smith . mat - the explict operator 384b0a32e0cSBarry Smith 385b0a32e0cSBarry Smith Notes: 386b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 387b0a32e0cSBarry Smith identity matrix. 388b0a32e0cSBarry Smith 389b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 390b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 391b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 392b0a32e0cSBarry Smith 393b0a32e0cSBarry Smith Level: advanced 394b0a32e0cSBarry Smith 395b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 396b0a32e0cSBarry Smith @*/ 3977087cfbeSBarry Smith PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 398b0a32e0cSBarry Smith { 399dfbe8321SBarry Smith PetscErrorCode ierr; 400b0a32e0cSBarry Smith MPI_Comm comm; 40138baddfdSBarry Smith PetscMPIInt size; 402b0a32e0cSBarry Smith 403b0a32e0cSBarry Smith PetscFunctionBegin; 4040700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 4054482741eSBarry Smith PetscValidPointer(mat,2); 406b0a32e0cSBarry Smith 407ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 408b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 409b3d09e86SJed Brown ierr = MatConvert_Shell(inmat,size == 1 ? MATSEQDENSE : MATAIJ,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 41024f910e3SHong Zhang PetscFunctionReturn(0); 41124f910e3SHong Zhang } 4124325cce7SMatthew G Knepley 4134325cce7SMatthew G Knepley /*@ 414f3b1f45cSBarry Smith MatComputeExplicitOperatorTranspose - 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: 426f3b1f45cSBarry Smith This computation is done by applying the transpose of the operator to columns of the 427f3b1f45cSBarry Smith identity matrix. 428f3b1f45cSBarry Smith 429f3b1f45cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 430f3b1f45cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 431f3b1f45cSBarry Smith and is recommended for use only with relatively small systems. 432f3b1f45cSBarry Smith 433f3b1f45cSBarry Smith Level: advanced 434f3b1f45cSBarry Smith 435f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator 436f3b1f45cSBarry Smith @*/ 437f3b1f45cSBarry Smith PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 438f3b1f45cSBarry Smith { 439f3b1f45cSBarry Smith Vec in,out; 440f3b1f45cSBarry Smith PetscErrorCode ierr; 441f3b1f45cSBarry Smith PetscInt i,m,n,M,N,*rows,start,end; 442f3b1f45cSBarry Smith MPI_Comm comm; 443f3b1f45cSBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 444f3b1f45cSBarry Smith PetscMPIInt size; 445f3b1f45cSBarry Smith 446f3b1f45cSBarry Smith PetscFunctionBegin; 447f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 448f3b1f45cSBarry Smith PetscValidPointer(mat,2); 449f3b1f45cSBarry Smith 450f3b1f45cSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 451f3b1f45cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 452f3b1f45cSBarry Smith 453f3b1f45cSBarry Smith ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 454f3b1f45cSBarry Smith ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 455f3b1f45cSBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 456f3b1f45cSBarry Smith ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 457f3b1f45cSBarry Smith ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 458f3b1f45cSBarry Smith ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 459f3b1f45cSBarry Smith for (i=0; i<m; i++) rows[i] = start + i; 460f3b1f45cSBarry Smith 461f3b1f45cSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 462f3b1f45cSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 463f3b1f45cSBarry Smith if (size == 1) { 464f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 465f3b1f45cSBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 466f3b1f45cSBarry Smith } else { 467f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 468f3b1f45cSBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 469f3b1f45cSBarry Smith } 470f3b1f45cSBarry Smith 471f3b1f45cSBarry Smith for (i=0; i<N; i++) { 472f3b1f45cSBarry Smith 473f3b1f45cSBarry Smith ierr = VecSet(in,zero);CHKERRQ(ierr); 474f3b1f45cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 475f3b1f45cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 476f3b1f45cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 477f3b1f45cSBarry Smith 478f3b1f45cSBarry Smith ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 479f3b1f45cSBarry Smith 480f3b1f45cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 481f3b1f45cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 482f3b1f45cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 483f3b1f45cSBarry Smith 484f3b1f45cSBarry Smith } 485f3b1f45cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 486f3b1f45cSBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 487f3b1f45cSBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 488f3b1f45cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 489f3b1f45cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 490f3b1f45cSBarry Smith PetscFunctionReturn(0); 491f3b1f45cSBarry Smith } 492f3b1f45cSBarry Smith 493f3b1f45cSBarry Smith /*@ 4944325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4954325cce7SMatthew G Knepley 4964325cce7SMatthew G Knepley Input Parameters: 4974325cce7SMatthew G Knepley + A - The matrix 4984325cce7SMatthew G Knepley - tol - The zero tolerance 4994325cce7SMatthew G Knepley 5004325cce7SMatthew G Knepley Output Parameters: 5014325cce7SMatthew G Knepley . A - The chopped matrix 5024325cce7SMatthew G Knepley 5034325cce7SMatthew G Knepley Level: intermediate 5044325cce7SMatthew G Knepley 5054325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 5063fc99919SSatish Balay @*/ 5074325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 5084325cce7SMatthew G Knepley { 5094325cce7SMatthew G Knepley PetscScalar *newVals; 5104325cce7SMatthew G Knepley PetscInt *newCols; 5114325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 5124325cce7SMatthew G Knepley PetscErrorCode ierr; 5134325cce7SMatthew G Knepley 5144325cce7SMatthew G Knepley PetscFunctionBegin; 5154325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 5164325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 5174325cce7SMatthew G Knepley PetscInt ncols; 5184325cce7SMatthew G Knepley 5190298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5204325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 5210298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5224325cce7SMatthew G Knepley } 5234325cce7SMatthew G Knepley numRows = rEnd - rStart; 524b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 525dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 5264325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 5274325cce7SMatthew G Knepley const PetscScalar *vals; 5284325cce7SMatthew G Knepley const PetscInt *cols; 529fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 5304325cce7SMatthew G Knepley 5314325cce7SMatthew G Knepley if (r < rEnd) { 5324325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 5334325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 5344325cce7SMatthew G Knepley newCols[c] = cols[c]; 5354325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 5364325cce7SMatthew G Knepley } 537fad45679SMatthew G. Knepley newcols = ncols; 5384325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 539fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 5404325cce7SMatthew G Knepley } 5414325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5424325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5434325cce7SMatthew G Knepley } 5444325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 5454325cce7SMatthew G Knepley PetscFunctionReturn(0); 5464325cce7SMatthew G Knepley } 547