16f79c3a4SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 4*d54c938cSPierre Jolivet static PetscErrorCode MatTransposeAXPY_Private(Mat Y,PetscScalar a,Mat X,MatStructure str,Mat T) 5*d54c938cSPierre Jolivet { 6*d54c938cSPierre Jolivet PetscErrorCode ierr,(*f)(Mat,Mat*); 7*d54c938cSPierre Jolivet Mat A,F; 8*d54c938cSPierre Jolivet 9*d54c938cSPierre Jolivet PetscFunctionBegin; 10*d54c938cSPierre Jolivet ierr = PetscObjectQueryFunction((PetscObject)T,"MatTransposeGetMat_C",&f);CHKERRQ(ierr); 11*d54c938cSPierre Jolivet if (f) { 12*d54c938cSPierre Jolivet ierr = MatTransposeGetMat(T,&A);CHKERRQ(ierr); 13*d54c938cSPierre Jolivet if (T == X) { 14*d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Explicitly transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 15*d54c938cSPierre Jolivet ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 16*d54c938cSPierre Jolivet A = Y; 17*d54c938cSPierre Jolivet } else { 18*d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 19*d54c938cSPierre Jolivet ierr = MatTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 20*d54c938cSPierre Jolivet } 21*d54c938cSPierre Jolivet } else { 22*d54c938cSPierre Jolivet ierr = MatHermitianTransposeGetMat(T,&A);CHKERRQ(ierr); 23*d54c938cSPierre Jolivet if (T == X) { 24*d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Explicitly Hermitian transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 25*d54c938cSPierre Jolivet ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 26*d54c938cSPierre Jolivet A = Y; 27*d54c938cSPierre Jolivet } else { 28*d54c938cSPierre Jolivet ierr = PetscInfo(NULL,"Hermitian transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr); 29*d54c938cSPierre Jolivet ierr = MatHermitianTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr); 30*d54c938cSPierre Jolivet } 31*d54c938cSPierre Jolivet } 32*d54c938cSPierre Jolivet ierr = MatAXPY(A,a,F,str);CHKERRQ(ierr); 33*d54c938cSPierre Jolivet ierr = MatDestroy(&F);CHKERRQ(ierr); 34*d54c938cSPierre Jolivet PetscFunctionReturn(0); 35*d54c938cSPierre Jolivet } 36*d54c938cSPierre 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 519cf4f1e8SLois Curfman McInnes .keywords: matrix, add 52d4bb536fSBarry Smith 532860a424SLois Curfman McInnes .seealso: MatAYPX() 5406be10caSBarry Smith @*/ 557087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 566f79c3a4SBarry Smith { 57*d54c938cSPierre Jolivet PetscErrorCode ierr; 58646531bbSStefano Zampini PetscInt M1,M2,N1,N2; 59c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 60a683ea4aSPierre Jolivet MatType t1,t2; 61a683ea4aSPierre Jolivet PetscBool sametype,transpose; 626f79c3a4SBarry Smith 633a40ed3dSBarry Smith PetscFunctionBegin; 640700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 65c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 66646531bbSStefano Zampini PetscValidHeaderSpecific(X,MAT_CLASSID,3); 67646531bbSStefano Zampini PetscCheckSameComm(Y,1,X,3); 68646531bbSStefano Zampini ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr); 69646531bbSStefano Zampini ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr); 70646531bbSStefano Zampini ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr); 71646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr); 72646531bbSStefano 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); 73646531bbSStefano 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); 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) { 84*d54c938cSPierre 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) { 88*d54c938cSPierre 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); 103fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 104b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 105b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 106c41cb2e2SAlejandro Lamas Daviña } 107d67ff14aSKarl Rupp #endif 108607cd303SBarry Smith PetscFunctionReturn(0); 109607cd303SBarry Smith } 110607cd303SBarry Smith 111646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 112646531bbSStefano Zampini { 113646531bbSStefano Zampini PetscErrorCode ierr; 11427afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 115646531bbSStefano Zampini 116646531bbSStefano Zampini PetscFunctionBegin; 117646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 118646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 119646531bbSStefano Zampini if (preall) { 120646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 121646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 122646531bbSStefano Zampini Mat preallocator; 123646531bbSStefano Zampini PetscInt r,rstart,rend; 124646531bbSStefano Zampini PetscInt m,n,M,N; 125646531bbSStefano Zampini 126646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 127646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 128646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 129646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 130646531bbSStefano Zampini ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 131646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 132646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 133646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 134646531bbSStefano Zampini PetscInt ncols; 135646531bbSStefano Zampini const PetscInt *row; 136646531bbSStefano Zampini const PetscScalar *vals; 137646531bbSStefano Zampini 138646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 139646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 140646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 141646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 142646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 143646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 144646531bbSStefano Zampini } 145646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147646531bbSStefano Zampini 148646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 149646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 150646531bbSStefano Zampini ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr); 151646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 152646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 153646531bbSStefano Zampini } 154646531bbSStefano Zampini PetscFunctionReturn(0); 155646531bbSStefano Zampini } 156646531bbSStefano Zampini 157f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 158607cd303SBarry Smith { 15938baddfdSBarry Smith PetscInt i,start,end,j,ncols,m,n; 1606849ba73SBarry Smith PetscErrorCode ierr; 16138baddfdSBarry Smith const PetscInt *row; 162b3cc6726SBarry Smith PetscScalar *val; 163b3cc6726SBarry Smith const PetscScalar *vals; 164607cd303SBarry Smith 165607cd303SBarry Smith PetscFunctionBegin; 1668dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 16790f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);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 } 1906d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1916d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 1936f79c3a4SBarry Smith } 194052efed2SBarry Smith 195ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 196ec7775f6SShri Abhyankar { 197ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 198ec7775f6SShri Abhyankar PetscErrorCode ierr; 199ec7775f6SShri Abhyankar const PetscInt *row; 200ec7775f6SShri Abhyankar PetscScalar *val; 201ec7775f6SShri Abhyankar const PetscScalar *vals; 202ec7775f6SShri Abhyankar 203ec7775f6SShri Abhyankar PetscFunctionBegin; 204ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 205ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 206ec7775f6SShri Abhyankar if (a == 1.0) { 207ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 208ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 209ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 210ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 211ec7775f6SShri Abhyankar 212ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 213ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 214ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 215ec7775f6SShri Abhyankar } 216ec7775f6SShri Abhyankar } else { 21721a3365eSStefano Zampini PetscInt vs = 100; 21821a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 21921a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 220ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 221ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 222ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 223ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 224ec7775f6SShri Abhyankar 225ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 22621a3365eSStefano Zampini if (vs < ncols) { 22721a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 22821a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 229ec7775f6SShri Abhyankar } 23021a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 231ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 232ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 233ec7775f6SShri Abhyankar } 234ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 235ec7775f6SShri Abhyankar } 236ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 237ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 238ec7775f6SShri Abhyankar PetscFunctionReturn(0); 239ec7775f6SShri Abhyankar } 240ec7775f6SShri Abhyankar 241052efed2SBarry Smith /*@ 24287828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 243052efed2SBarry Smith 2443f9fe445SBarry Smith Neighbor-wise Collective on Mat 245fee21e36SBarry Smith 24698a79cdbSBarry Smith Input Parameters: 24798a79cdbSBarry Smith + Y - the matrices 24887828ca2SBarry Smith - a - the PetscScalar 24998a79cdbSBarry Smith 2502860a424SLois Curfman McInnes Level: intermediate 2512860a424SLois Curfman McInnes 25295452b02SPatrick Sanan Notes: 25395452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2546f33a894SBarry 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 2556f33a894SBarry Smith entry). 2566f33a894SBarry Smith 2570c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2580c0fd78eSBarry Smith 2596f33a894SBarry Smith Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 2606f33a894SBarry Smith preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 2616f33a894SBarry Smith 262052efed2SBarry Smith .keywords: matrix, add, shift 2636b9ee512SLois Curfman McInnes 2640c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 265052efed2SBarry Smith @*/ 2667087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 267052efed2SBarry Smith { 2686849ba73SBarry Smith PetscErrorCode ierr; 269052efed2SBarry Smith 2703a40ed3dSBarry Smith PetscFunctionBegin; 2710700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 272ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 273ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2744994cf47SJed Brown MatCheckPreallocated(Y,1); 275b50b34bdSBarry Smith 2761c738b47SStefano Zampini if (Y->ops->shift) { 277f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 2781c738b47SStefano Zampini } else { 2791c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 2801c738b47SStefano Zampini } 2817d68702bSBarry Smith 2825528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 283fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 284b8ced49eSKarl Rupp if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 285b8ced49eSKarl Rupp Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU; 286c41cb2e2SAlejandro Lamas Daviña } 287d67ff14aSKarl Rupp #endif 2883a40ed3dSBarry Smith PetscFunctionReturn(0); 289052efed2SBarry Smith } 2906d84be18SBarry Smith 2917087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 29209f38230SBarry Smith { 29309f38230SBarry Smith PetscErrorCode ierr; 29467576b8bSBarry Smith PetscInt i,start,end; 29509f38230SBarry Smith PetscScalar *v; 29609f38230SBarry Smith 29709f38230SBarry Smith PetscFunctionBegin; 29809f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 29909f38230SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 30009f38230SBarry Smith for (i=start; i<end; i++) { 30109f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 30209f38230SBarry Smith } 30309f38230SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 30409f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30509f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30609f38230SBarry Smith PetscFunctionReturn(0); 30709f38230SBarry Smith } 30809f38230SBarry Smith 3096d84be18SBarry Smith /*@ 310f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 311f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 312f56f2b3fSBarry Smith INSERT_VALUES. 3136d84be18SBarry Smith 3146d84be18SBarry Smith Input Parameters: 31598a79cdbSBarry Smith + Y - the input matrix 316f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 317f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3186d84be18SBarry Smith 3193f9fe445SBarry Smith Neighbor-wise Collective on Mat and Vec 320fee21e36SBarry Smith 32195452b02SPatrick Sanan Notes: 32295452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3236f33a894SBarry 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 3246f33a894SBarry Smith entry). 3256f33a894SBarry Smith 3262860a424SLois Curfman McInnes Level: intermediate 3272860a424SLois Curfman McInnes 3286b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal 3296b9ee512SLois Curfman McInnes 3300c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 3316d84be18SBarry Smith @*/ 3327087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 3336d84be18SBarry Smith { 3346849ba73SBarry Smith PetscErrorCode ierr; 33567576b8bSBarry Smith PetscInt matlocal,veclocal; 3366d84be18SBarry Smith 3373a40ed3dSBarry Smith PetscFunctionBegin; 3380700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3390700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 34067576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 34167576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 34267576b8bSBarry 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); 343f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 344f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 34594d884c6SBarry Smith } else { 34609f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3476d84be18SBarry Smith } 3485528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3493a40ed3dSBarry Smith PetscFunctionReturn(0); 3506d84be18SBarry Smith } 351d4bb536fSBarry Smith 352d4bb536fSBarry Smith /*@ 35304aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 354d4bb536fSBarry Smith 3553f9fe445SBarry Smith Logically on Mat 356fee21e36SBarry Smith 35798a79cdbSBarry Smith Input Parameters: 35804aac2b0SHong Zhang + a - the PetscScalar multiplier 35904aac2b0SHong Zhang . Y - the first matrix 36004aac2b0SHong Zhang . X - the second matrix 36104aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 362d4bb536fSBarry Smith 3632860a424SLois Curfman McInnes Level: intermediate 3642860a424SLois Curfman McInnes 365d4bb536fSBarry Smith .keywords: matrix, add 366d4bb536fSBarry Smith 3672860a424SLois Curfman McInnes .seealso: MatAXPY() 368d4bb536fSBarry Smith @*/ 3697087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 370d4bb536fSBarry Smith { 37187828ca2SBarry Smith PetscScalar one = 1.0; 3726849ba73SBarry Smith PetscErrorCode ierr; 37338baddfdSBarry Smith PetscInt mX,mY,nX,nY; 374d4bb536fSBarry Smith 3753a40ed3dSBarry Smith PetscFunctionBegin; 376c5eb9154SBarry Smith PetscValidHeaderSpecific(X,MAT_CLASSID,3); 3770700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 378c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 379329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 380329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 381e32f2f54SBarry 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); 382d4bb536fSBarry Smith 383f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 384cb9801acSJed Brown ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 3853a40ed3dSBarry Smith PetscFunctionReturn(0); 386d4bb536fSBarry Smith } 387b0a32e0cSBarry Smith 388b0a32e0cSBarry Smith /*@ 389b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 390b0a32e0cSBarry Smith 391b0a32e0cSBarry Smith Collective on Mat 392b0a32e0cSBarry Smith 393b0a32e0cSBarry Smith Input Parameter: 394b0a32e0cSBarry Smith . inmat - the matrix 395b0a32e0cSBarry Smith 396b0a32e0cSBarry Smith Output Parameter: 397f3b1f45cSBarry Smith . mat - the explict operator 398b0a32e0cSBarry Smith 399b0a32e0cSBarry Smith Notes: 400b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 401b0a32e0cSBarry Smith identity matrix. 402b0a32e0cSBarry Smith 403b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 404b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 405b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 406b0a32e0cSBarry Smith 407b0a32e0cSBarry Smith Level: advanced 408b0a32e0cSBarry Smith 409b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 410b0a32e0cSBarry Smith @*/ 4117087cfbeSBarry Smith PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 412b0a32e0cSBarry Smith { 413dfbe8321SBarry Smith PetscErrorCode ierr; 414b0a32e0cSBarry Smith MPI_Comm comm; 41538baddfdSBarry Smith PetscMPIInt size; 416b0a32e0cSBarry Smith 417b0a32e0cSBarry Smith PetscFunctionBegin; 4180700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 4194482741eSBarry Smith PetscValidPointer(mat,2); 420b0a32e0cSBarry Smith 421ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 422b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 423b3d09e86SJed Brown ierr = MatConvert_Shell(inmat,size == 1 ? MATSEQDENSE : MATAIJ,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 42424f910e3SHong Zhang PetscFunctionReturn(0); 42524f910e3SHong Zhang } 4264325cce7SMatthew G Knepley 4274325cce7SMatthew G Knepley /*@ 428f3b1f45cSBarry Smith MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of 429f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 430f3b1f45cSBarry Smith 431f3b1f45cSBarry Smith Collective on Mat 432f3b1f45cSBarry Smith 433f3b1f45cSBarry Smith Input Parameter: 434f3b1f45cSBarry Smith . inmat - the matrix 435f3b1f45cSBarry Smith 436f3b1f45cSBarry Smith Output Parameter: 437f3b1f45cSBarry Smith . mat - the explict operator transposed 438f3b1f45cSBarry Smith 439f3b1f45cSBarry Smith Notes: 440f3b1f45cSBarry Smith This computation is done by applying the transpose of the operator to columns of the 441f3b1f45cSBarry Smith identity matrix. 442f3b1f45cSBarry Smith 443f3b1f45cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 444f3b1f45cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 445f3b1f45cSBarry Smith and is recommended for use only with relatively small systems. 446f3b1f45cSBarry Smith 447f3b1f45cSBarry Smith Level: advanced 448f3b1f45cSBarry Smith 449f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator 450f3b1f45cSBarry Smith @*/ 451f3b1f45cSBarry Smith PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 452f3b1f45cSBarry Smith { 453f3b1f45cSBarry Smith Vec in,out; 454f3b1f45cSBarry Smith PetscErrorCode ierr; 455f3b1f45cSBarry Smith PetscInt i,m,n,M,N,*rows,start,end; 456f3b1f45cSBarry Smith MPI_Comm comm; 457f3b1f45cSBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 458f3b1f45cSBarry Smith PetscMPIInt size; 459f3b1f45cSBarry Smith 460f3b1f45cSBarry Smith PetscFunctionBegin; 461f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 462f3b1f45cSBarry Smith PetscValidPointer(mat,2); 463f3b1f45cSBarry Smith 464f3b1f45cSBarry Smith ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 465f3b1f45cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 466f3b1f45cSBarry Smith 467f3b1f45cSBarry Smith ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 468f3b1f45cSBarry Smith ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 469f3b1f45cSBarry Smith ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 470f3b1f45cSBarry Smith ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 471f3b1f45cSBarry Smith ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 472f3b1f45cSBarry Smith ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 473f3b1f45cSBarry Smith for (i=0; i<m; i++) rows[i] = start + i; 474f3b1f45cSBarry Smith 475f3b1f45cSBarry Smith ierr = MatCreate(comm,mat);CHKERRQ(ierr); 476f3b1f45cSBarry Smith ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 477f3b1f45cSBarry Smith if (size == 1) { 478f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 479f3b1f45cSBarry Smith ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 480f3b1f45cSBarry Smith } else { 481f3b1f45cSBarry Smith ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 482f3b1f45cSBarry Smith ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 483f3b1f45cSBarry Smith } 484f3b1f45cSBarry Smith 485f3b1f45cSBarry Smith for (i=0; i<N; i++) { 486f3b1f45cSBarry Smith 487f3b1f45cSBarry Smith ierr = VecSet(in,zero);CHKERRQ(ierr); 488f3b1f45cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 489f3b1f45cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 490f3b1f45cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 491f3b1f45cSBarry Smith 492f3b1f45cSBarry Smith ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 493f3b1f45cSBarry Smith 494f3b1f45cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 495f3b1f45cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 496f3b1f45cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 497f3b1f45cSBarry Smith 498f3b1f45cSBarry Smith } 499f3b1f45cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 500f3b1f45cSBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 501f3b1f45cSBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 502f3b1f45cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 503f3b1f45cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 504f3b1f45cSBarry Smith PetscFunctionReturn(0); 505f3b1f45cSBarry Smith } 506f3b1f45cSBarry Smith 507f3b1f45cSBarry Smith /*@ 5084325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 5094325cce7SMatthew G Knepley 5104325cce7SMatthew G Knepley Input Parameters: 5114325cce7SMatthew G Knepley + A - The matrix 5124325cce7SMatthew G Knepley - tol - The zero tolerance 5134325cce7SMatthew G Knepley 5144325cce7SMatthew G Knepley Output Parameters: 5154325cce7SMatthew G Knepley . A - The chopped matrix 5164325cce7SMatthew G Knepley 5174325cce7SMatthew G Knepley Level: intermediate 5184325cce7SMatthew G Knepley 5194325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 5203fc99919SSatish Balay @*/ 5214325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 5224325cce7SMatthew G Knepley { 5234325cce7SMatthew G Knepley PetscScalar *newVals; 5244325cce7SMatthew G Knepley PetscInt *newCols; 5254325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 5264325cce7SMatthew G Knepley PetscErrorCode ierr; 5274325cce7SMatthew G Knepley 5284325cce7SMatthew G Knepley PetscFunctionBegin; 5294325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 5304325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 5314325cce7SMatthew G Knepley PetscInt ncols; 5324325cce7SMatthew G Knepley 5330298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5344325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 5350298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5364325cce7SMatthew G Knepley } 5374325cce7SMatthew G Knepley numRows = rEnd - rStart; 538b2566f29SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 539dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 5404325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 5414325cce7SMatthew G Knepley const PetscScalar *vals; 5424325cce7SMatthew G Knepley const PetscInt *cols; 543fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 5444325cce7SMatthew G Knepley 5454325cce7SMatthew G Knepley if (r < rEnd) { 5464325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 5474325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 5484325cce7SMatthew G Knepley newCols[c] = cols[c]; 5494325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 5504325cce7SMatthew G Knepley } 551fad45679SMatthew G. Knepley newcols = ncols; 5524325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 553fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 5544325cce7SMatthew G Knepley } 5554325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5564325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5574325cce7SMatthew G Knepley } 5584325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 5594325cce7SMatthew G Knepley PetscFunctionReturn(0); 5604325cce7SMatthew G Knepley } 561