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 49edf9a46cSStefano Zampini Notes: No operation is performed when a is zero. 50edf9a46cSStefano Zampini 512860a424SLois Curfman McInnes Level: intermediate 522860a424SLois Curfman McInnes 532860a424SLois Curfman McInnes .seealso: MatAYPX() 5406be10caSBarry Smith @*/ 557087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 566f79c3a4SBarry Smith { 57d54c938cSPierre 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); 7247415954SStefano Zampini if (M1 != M2 || N1 != N2) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_SIZ,"Non conforming matrix add: global sizes X %D x %D, Y %D x %D",M1,N1,M2,N2); 7347415954SStefano Zampini if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: local sizes X %D x %D, Y %D x %D",m1,n1,m2,n2); 741ccb7edbSStefano Zampini if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (Y)"); 751ccb7edbSStefano Zampini if (!X->assembled) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (X)"); 76edf9a46cSStefano Zampini if (a == (PetscScalar)0.0) PetscFunctionReturn(0); 77edf9a46cSStefano Zampini if (Y == X) { 78edf9a46cSStefano Zampini ierr = MatScale(Y,1.0 + a);CHKERRQ(ierr); 79edf9a46cSStefano Zampini PetscFunctionReturn(0); 80edf9a46cSStefano Zampini } 81a683ea4aSPierre Jolivet ierr = MatGetType(X,&t1);CHKERRQ(ierr); 82a683ea4aSPierre Jolivet ierr = MatGetType(Y,&t2);CHKERRQ(ierr); 83a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,t2,&sametype);CHKERRQ(ierr); 84e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 8593c87e65SStefano Zampini if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) { 86f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 87d4bb536fSBarry Smith } else { 88a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 89a683ea4aSPierre Jolivet if (transpose) { 90d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,X);CHKERRQ(ierr); 91a683ea4aSPierre Jolivet } else { 92a683ea4aSPierre Jolivet ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 93a683ea4aSPierre Jolivet if (transpose) { 94d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,Y);CHKERRQ(ierr); 95a683ea4aSPierre Jolivet } else { 96f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 97607cd303SBarry Smith } 98a683ea4aSPierre Jolivet } 99a683ea4aSPierre Jolivet } 100e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 101607cd303SBarry Smith PetscFunctionReturn(0); 102607cd303SBarry Smith } 103607cd303SBarry Smith 104646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 105646531bbSStefano Zampini { 106646531bbSStefano Zampini PetscErrorCode ierr; 10727afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 108646531bbSStefano Zampini 109646531bbSStefano Zampini PetscFunctionBegin; 110646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 111646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 112646531bbSStefano Zampini if (preall) { 113646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 114646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 115646531bbSStefano Zampini Mat preallocator; 116646531bbSStefano Zampini PetscInt r,rstart,rend; 117646531bbSStefano Zampini PetscInt m,n,M,N; 118646531bbSStefano Zampini 1196eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 1206eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 121646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 122646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 123646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 124646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 125a6ab3590SBarry Smith ierr = MatSetLayouts(preallocator,Y->rmap,Y->cmap);CHKERRQ(ierr); 126646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 127646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 128646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 129646531bbSStefano Zampini PetscInt ncols; 130646531bbSStefano Zampini const PetscInt *row; 131646531bbSStefano Zampini const PetscScalar *vals; 132646531bbSStefano Zampini 133646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 134646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 135646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 136646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 137646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 138646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 139646531bbSStefano Zampini } 140b8550948SPierre Jolivet ierr = MatSetOption(preallocator,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 141646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 142646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1436eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 1446eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 145646531bbSStefano Zampini 146646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 147646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 148a6ab3590SBarry Smith ierr = MatSetLayouts(*B,Y->rmap,Y->cmap);CHKERRQ(ierr); 149646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 150646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 151646531bbSStefano Zampini } 152646531bbSStefano Zampini PetscFunctionReturn(0); 153646531bbSStefano Zampini } 154646531bbSStefano Zampini 155f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 156607cd303SBarry Smith { 1576849ba73SBarry Smith PetscErrorCode ierr; 158be705e3aSPierre Jolivet PetscBool isshell,isdense,isnest; 159edf9a46cSStefano Zampini 160edf9a46cSStefano Zampini PetscFunctionBegin; 161e6aea9c0SStefano Zampini ierr = MatIsShell(Y,&isshell);CHKERRQ(ierr); 162edf9a46cSStefano Zampini if (isshell) { /* MatShell has special support for AXPY */ 163edf9a46cSStefano Zampini PetscErrorCode (*f)(Mat,PetscScalar,Mat,MatStructure); 164edf9a46cSStefano Zampini 165edf9a46cSStefano Zampini ierr = MatGetOperation(Y,MATOP_AXPY,(void (**)(void))&f);CHKERRQ(ierr); 166edf9a46cSStefano Zampini if (f) { 167edf9a46cSStefano Zampini ierr = (*f)(Y,a,X,str);CHKERRQ(ierr); 168edf9a46cSStefano Zampini PetscFunctionReturn(0); 169edf9a46cSStefano Zampini } 170edf9a46cSStefano Zampini } 17193c87e65SStefano Zampini /* no need to preallocate if Y is dense */ 17293c87e65SStefano Zampini ierr = PetscObjectBaseTypeCompareAny((PetscObject)Y,&isdense,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 173be705e3aSPierre Jolivet if (isdense) { 174be705e3aSPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)X,MATNEST,&isnest);CHKERRQ(ierr); 175be705e3aSPierre Jolivet if (isnest) { 176be705e3aSPierre Jolivet ierr = MatAXPY_Dense_Nest(Y,a,X);CHKERRQ(ierr); 177be705e3aSPierre Jolivet PetscFunctionReturn(0); 178be705e3aSPierre Jolivet } 179be705e3aSPierre Jolivet if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN; 180be705e3aSPierre Jolivet } 181a6ab3590SBarry Smith if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) { 182edf9a46cSStefano Zampini PetscInt i,start,end,j,ncols,m,n; 18338baddfdSBarry Smith const PetscInt *row; 184b3cc6726SBarry Smith PetscScalar *val; 185b3cc6726SBarry Smith const PetscScalar *vals; 186607cd303SBarry Smith 1878dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 18890f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 1896eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 190f4df32b1SMatthew Knepley if (a == 1.0) { 191d4bb536fSBarry Smith for (i = start; i < end; i++) { 192d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 193d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 194d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 195d4bb536fSBarry Smith } 196d4bb536fSBarry Smith } else { 19721a3365eSStefano Zampini PetscInt vs = 100; 19821a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 19921a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 20006be10caSBarry Smith for (i=start; i<end; i++) { 201b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 20221a3365eSStefano Zampini if (vs < ncols) { 20321a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 20421a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 2056f79c3a4SBarry Smith } 20621a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 207b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 208b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 2096f79c3a4SBarry Smith } 210b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 211d4bb536fSBarry Smith } 2126eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 2136d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2146d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 215edf9a46cSStefano Zampini } else { 216edf9a46cSStefano Zampini Mat B; 217edf9a46cSStefano Zampini 218edf9a46cSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 219edf9a46cSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 22093c87e65SStefano Zampini /* TODO mat_tests-ex37_nsize-1_mat_type-baij_mat_block_size-2 fails with MatHeaderMerge */ 221edf9a46cSStefano Zampini ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 222edf9a46cSStefano Zampini } 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 2246f79c3a4SBarry Smith } 225052efed2SBarry Smith 226ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 227ec7775f6SShri Abhyankar { 228ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 229ec7775f6SShri Abhyankar PetscErrorCode ierr; 230ec7775f6SShri Abhyankar const PetscInt *row; 231ec7775f6SShri Abhyankar PetscScalar *val; 232ec7775f6SShri Abhyankar const PetscScalar *vals; 233ec7775f6SShri Abhyankar 234ec7775f6SShri Abhyankar PetscFunctionBegin; 235ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 236ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 2376eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 2386eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 239ec7775f6SShri Abhyankar if (a == 1.0) { 240ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 241ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 242ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 243ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 244ec7775f6SShri Abhyankar 245ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 246ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 247ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 248ec7775f6SShri Abhyankar } 249ec7775f6SShri Abhyankar } else { 25021a3365eSStefano Zampini PetscInt vs = 100; 25121a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 25221a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 253ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 254ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 255ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 256ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 257ec7775f6SShri Abhyankar 258ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 25921a3365eSStefano Zampini if (vs < ncols) { 26021a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 26121a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 262ec7775f6SShri Abhyankar } 26321a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 264ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 265ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 266ec7775f6SShri Abhyankar } 267ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 268ec7775f6SShri Abhyankar } 2696eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 2706eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 271ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 272ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 273ec7775f6SShri Abhyankar PetscFunctionReturn(0); 274ec7775f6SShri Abhyankar } 275ec7775f6SShri Abhyankar 276052efed2SBarry Smith /*@ 27787828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 278052efed2SBarry Smith 2793f9fe445SBarry Smith Neighbor-wise Collective on Mat 280fee21e36SBarry Smith 28198a79cdbSBarry Smith Input Parameters: 28298a79cdbSBarry Smith + Y - the matrices 28387828ca2SBarry Smith - a - the PetscScalar 28498a79cdbSBarry Smith 2852860a424SLois Curfman McInnes Level: intermediate 2862860a424SLois Curfman McInnes 28795452b02SPatrick Sanan Notes: 28895452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2896f33a894SBarry 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 290edf9a46cSStefano Zampini entry). No operation is performed when a is zero. 2916f33a894SBarry Smith 2920c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2930c0fd78eSBarry Smith 2940c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 295052efed2SBarry Smith @*/ 2967087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 297052efed2SBarry Smith { 2986849ba73SBarry Smith PetscErrorCode ierr; 299052efed2SBarry Smith 3003a40ed3dSBarry Smith PetscFunctionBegin; 3010700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 302ce94432eSBarry Smith if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 303ce94432eSBarry Smith if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 3044994cf47SJed Brown MatCheckPreallocated(Y,1); 3051b71c1a7SBarry Smith if (a == 0.0) PetscFunctionReturn(0); 306b50b34bdSBarry Smith 3071c738b47SStefano Zampini if (Y->ops->shift) { 308f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 3091c738b47SStefano Zampini } else { 3101c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 3111c738b47SStefano Zampini } 3127d68702bSBarry Smith 3135528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3143a40ed3dSBarry Smith PetscFunctionReturn(0); 315052efed2SBarry Smith } 3166d84be18SBarry Smith 3177087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 31809f38230SBarry Smith { 31909f38230SBarry Smith PetscErrorCode ierr; 32067576b8bSBarry Smith PetscInt i,start,end; 321fa2e578bSStefano Zampini const PetscScalar *v; 32209f38230SBarry Smith 32309f38230SBarry Smith PetscFunctionBegin; 32409f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 325fa2e578bSStefano Zampini ierr = VecGetArrayRead(D,&v);CHKERRQ(ierr); 32609f38230SBarry Smith for (i=start; i<end; i++) { 32709f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 32809f38230SBarry Smith } 329fa2e578bSStefano Zampini ierr = VecRestoreArrayRead(D,&v);CHKERRQ(ierr); 33009f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 33109f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 33209f38230SBarry Smith PetscFunctionReturn(0); 33309f38230SBarry Smith } 33409f38230SBarry Smith 3356d84be18SBarry Smith /*@ 336f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 337f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 338f56f2b3fSBarry Smith INSERT_VALUES. 3396d84be18SBarry Smith 34048e586daSJose E. Roman Neighbor-wise Collective on Mat 34148e586daSJose E. Roman 3426d84be18SBarry Smith Input Parameters: 34398a79cdbSBarry Smith + Y - the input matrix 344f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 345f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3466d84be18SBarry Smith 34795452b02SPatrick Sanan Notes: 34895452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3496f33a894SBarry 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 3506f33a894SBarry Smith entry). 3516f33a894SBarry Smith 3522860a424SLois Curfman McInnes Level: intermediate 3532860a424SLois Curfman McInnes 3540c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 3556d84be18SBarry Smith @*/ 3567087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 3576d84be18SBarry Smith { 3586849ba73SBarry Smith PetscErrorCode ierr; 35967576b8bSBarry Smith PetscInt matlocal,veclocal; 3606d84be18SBarry Smith 3613a40ed3dSBarry Smith PetscFunctionBegin; 3620700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3630700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 36467576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 36567576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 36667576b8bSBarry 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); 367f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 368f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 36994d884c6SBarry Smith } else { 37009f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3716d84be18SBarry Smith } 3725528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3733a40ed3dSBarry Smith PetscFunctionReturn(0); 3746d84be18SBarry Smith } 375d4bb536fSBarry Smith 376d4bb536fSBarry Smith /*@ 37704aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 378d4bb536fSBarry Smith 3793f9fe445SBarry Smith Logically on Mat 380fee21e36SBarry Smith 38198a79cdbSBarry Smith Input Parameters: 38204aac2b0SHong Zhang + a - the PetscScalar multiplier 38304aac2b0SHong Zhang . Y - the first matrix 38404aac2b0SHong Zhang . X - the second matrix 38504aac2b0SHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 386d4bb536fSBarry Smith 3872860a424SLois Curfman McInnes Level: intermediate 3882860a424SLois Curfman McInnes 3892860a424SLois Curfman McInnes .seealso: MatAXPY() 390d4bb536fSBarry Smith @*/ 3917087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 392d4bb536fSBarry Smith { 3936849ba73SBarry Smith PetscErrorCode ierr; 394d4bb536fSBarry Smith 3953a40ed3dSBarry Smith PetscFunctionBegin; 396f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 3972ab6f6a8SStefano Zampini ierr = MatAXPY(Y,1.0,X,str);CHKERRQ(ierr); 3983a40ed3dSBarry Smith PetscFunctionReturn(0); 399d4bb536fSBarry Smith } 400b0a32e0cSBarry Smith 401b0a32e0cSBarry Smith /*@ 4020bacdadaSStefano Zampini MatComputeOperator - Computes the explicit matrix 403b0a32e0cSBarry Smith 404b0a32e0cSBarry Smith Collective on Mat 405b0a32e0cSBarry Smith 406b0a32e0cSBarry Smith Input Parameter: 407186a3e20SStefano Zampini + inmat - the matrix 408186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator 409b0a32e0cSBarry Smith 410b0a32e0cSBarry Smith Output Parameter: 411f3b1f45cSBarry Smith . mat - the explict operator 412b0a32e0cSBarry Smith 413b0a32e0cSBarry Smith Notes: 414186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 415186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 416186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 417b0a32e0cSBarry Smith 418b0a32e0cSBarry Smith Level: advanced 419b0a32e0cSBarry Smith 420b0a32e0cSBarry Smith @*/ 4210bacdadaSStefano Zampini PetscErrorCode MatComputeOperator(Mat inmat,MatType mattype,Mat *mat) 422b0a32e0cSBarry Smith { 423dfbe8321SBarry Smith PetscErrorCode ierr; 424b0a32e0cSBarry Smith 425b0a32e0cSBarry Smith PetscFunctionBegin; 4260700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 427186a3e20SStefano Zampini PetscValidPointer(mat,3); 428186a3e20SStefano Zampini ierr = MatConvert_Shell(inmat,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 42924f910e3SHong Zhang PetscFunctionReturn(0); 43024f910e3SHong Zhang } 4314325cce7SMatthew G Knepley 4324325cce7SMatthew G Knepley /*@ 4330bacdadaSStefano Zampini MatComputeOperatorTranspose - Computes the explicit matrix representation of 434f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 435f3b1f45cSBarry Smith 436f3b1f45cSBarry Smith Collective on Mat 437f3b1f45cSBarry Smith 438f3b1f45cSBarry Smith Input Parameter: 439f3b1f45cSBarry Smith . inmat - the matrix 440f3b1f45cSBarry Smith 441f3b1f45cSBarry Smith Output Parameter: 442f3b1f45cSBarry Smith . mat - the explict operator transposed 443f3b1f45cSBarry Smith 444f3b1f45cSBarry Smith Notes: 445186a3e20SStefano Zampini This computation is done by applying the transpose of the operator to columns of the identity matrix. 446186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 447186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 448f3b1f45cSBarry Smith 449f3b1f45cSBarry Smith Level: advanced 450f3b1f45cSBarry Smith 451f3b1f45cSBarry Smith @*/ 4520bacdadaSStefano Zampini PetscErrorCode MatComputeOperatorTranspose(Mat inmat,MatType mattype,Mat *mat) 453f3b1f45cSBarry Smith { 454186a3e20SStefano Zampini Mat A; 455f3b1f45cSBarry Smith PetscErrorCode ierr; 456f3b1f45cSBarry Smith 457f3b1f45cSBarry Smith PetscFunctionBegin; 458f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 459186a3e20SStefano Zampini PetscValidPointer(mat,3); 460186a3e20SStefano Zampini ierr = MatCreateTranspose(inmat,&A);CHKERRQ(ierr); 461186a3e20SStefano Zampini ierr = MatConvert_Shell(A,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 462186a3e20SStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 463f3b1f45cSBarry Smith PetscFunctionReturn(0); 464f3b1f45cSBarry Smith } 465f3b1f45cSBarry Smith 466f3b1f45cSBarry Smith /*@ 4674325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4684325cce7SMatthew G Knepley 4694325cce7SMatthew G Knepley Input Parameters: 4704325cce7SMatthew G Knepley + A - The matrix 4714325cce7SMatthew G Knepley - tol - The zero tolerance 4724325cce7SMatthew G Knepley 4734325cce7SMatthew G Knepley Output Parameters: 4744325cce7SMatthew G Knepley . A - The chopped matrix 4754325cce7SMatthew G Knepley 4764325cce7SMatthew G Knepley Level: intermediate 4774325cce7SMatthew G Knepley 4784325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 4793fc99919SSatish Balay @*/ 4804325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 4814325cce7SMatthew G Knepley { 482038df967SPierre Jolivet Mat a; 4834325cce7SMatthew G Knepley PetscScalar *newVals; 4844325cce7SMatthew G Knepley PetscInt *newCols; 4854325cce7SMatthew G Knepley PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 486038df967SPierre Jolivet PetscBool flg; 4874325cce7SMatthew G Knepley PetscErrorCode ierr; 4884325cce7SMatthew G Knepley 4894325cce7SMatthew G Knepley PetscFunctionBegin; 490038df967SPierre Jolivet ierr = PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, "");CHKERRQ(ierr); 491038df967SPierre Jolivet if (flg) { 492038df967SPierre Jolivet ierr = MatDenseGetLocalMatrix(A, &a);CHKERRQ(ierr); 493038df967SPierre Jolivet ierr = MatDenseGetLDA(a, &r);CHKERRQ(ierr); 494038df967SPierre Jolivet ierr = MatGetSize(a, &rStart, &rEnd);CHKERRQ(ierr); 495038df967SPierre Jolivet ierr = MatDenseGetArray(a, &newVals);CHKERRQ(ierr); 496038df967SPierre Jolivet for (; colMax < rEnd; ++colMax) { 497038df967SPierre Jolivet for (maxRows = 0; maxRows < rStart; ++maxRows) { 498038df967SPierre Jolivet newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) < tol ? 0.0 : newVals[maxRows + colMax * r]; 499038df967SPierre Jolivet } 500038df967SPierre Jolivet } 501038df967SPierre Jolivet ierr = MatDenseRestoreArray(a, &newVals);CHKERRQ(ierr); 502038df967SPierre Jolivet } else { 5034325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 504a4c05331SJacob Faibussowitsch ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 5054325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 5064325cce7SMatthew G Knepley PetscInt ncols; 5074325cce7SMatthew G Knepley 5080298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5094325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 5100298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5114325cce7SMatthew G Knepley } 5124325cce7SMatthew G Knepley numRows = rEnd - rStart; 513*820f2d46SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 514dcca6d9dSJed Brown ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 515038df967SPierre Jolivet ierr = MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);CHKERRQ(ierr); 5164325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 5174325cce7SMatthew G Knepley const PetscScalar *vals; 5184325cce7SMatthew G Knepley const PetscInt *cols; 519fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 5204325cce7SMatthew G Knepley 5214325cce7SMatthew G Knepley if (r < rEnd) { 5224325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 5234325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 5244325cce7SMatthew G Knepley newCols[c] = cols[c]; 5254325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 5264325cce7SMatthew G Knepley } 527fad45679SMatthew G. Knepley newcols = ncols; 5284325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 529fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 5304325cce7SMatthew G Knepley } 5314325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5324325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5334325cce7SMatthew G Knepley } 534a4c05331SJacob Faibussowitsch ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 5354325cce7SMatthew G Knepley ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 536038df967SPierre Jolivet } 5374325cce7SMatthew G Knepley PetscFunctionReturn(0); 5384325cce7SMatthew G Knepley } 539