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 46531d0c6aSBarry Smith - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, UNKNOWN_NONZERO_PATTERN, or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 4798a79cdbSBarry Smith 48edf9a46cSStefano Zampini Notes: No operation is performed when a is zero. 49edf9a46cSStefano Zampini 502860a424SLois Curfman McInnes Level: intermediate 512860a424SLois Curfman McInnes 522860a424SLois Curfman McInnes .seealso: MatAYPX() 5306be10caSBarry Smith @*/ 547087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 556f79c3a4SBarry Smith { 56d54c938cSPierre Jolivet PetscErrorCode ierr; 57646531bbSStefano Zampini PetscInt M1,M2,N1,N2; 58c1ac3661SBarry Smith PetscInt m1,m2,n1,n2; 59a683ea4aSPierre Jolivet MatType t1,t2; 60a683ea4aSPierre Jolivet PetscBool sametype,transpose; 616f79c3a4SBarry Smith 623a40ed3dSBarry Smith PetscFunctionBegin; 630700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 64c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(Y,a,2); 65646531bbSStefano Zampini PetscValidHeaderSpecific(X,MAT_CLASSID,3); 66646531bbSStefano Zampini PetscCheckSameComm(Y,1,X,3); 67646531bbSStefano Zampini ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr); 68646531bbSStefano Zampini ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr); 69646531bbSStefano Zampini ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr); 70646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr); 71*2c71b3e2SJacob Faibussowitsch PetscCheck(M1 == M2 && N1 == N2,PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_SIZ,"Non conforming matrix add: global sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT,M1,N1,M2,N2); 72*2c71b3e2SJacob Faibussowitsch PetscCheck(m1 == m2 && n1 == n2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: local sizes X %" PetscInt_FMT " x %" PetscInt_FMT ", Y %" PetscInt_FMT " x %" PetscInt_FMT,m1,n1,m2,n2); 73*2c71b3e2SJacob Faibussowitsch PetscCheck(Y->assembled,PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (Y)"); 74*2c71b3e2SJacob Faibussowitsch PetscCheck(X->assembled,PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (X)"); 75edf9a46cSStefano Zampini if (a == (PetscScalar)0.0) PetscFunctionReturn(0); 76edf9a46cSStefano Zampini if (Y == X) { 77edf9a46cSStefano Zampini ierr = MatScale(Y,1.0 + a);CHKERRQ(ierr); 78edf9a46cSStefano Zampini PetscFunctionReturn(0); 79edf9a46cSStefano Zampini } 80a683ea4aSPierre Jolivet ierr = MatGetType(X,&t1);CHKERRQ(ierr); 81a683ea4aSPierre Jolivet ierr = MatGetType(Y,&t2);CHKERRQ(ierr); 82a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,t2,&sametype);CHKERRQ(ierr); 83e8136da8SHong Zhang ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 8493c87e65SStefano Zampini if (Y->ops->axpy && (sametype || X->ops->axpy == Y->ops->axpy)) { 85f4df32b1SMatthew Knepley ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 86d4bb536fSBarry Smith } else { 87a683ea4aSPierre Jolivet ierr = PetscStrcmp(t1,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 88a683ea4aSPierre Jolivet if (transpose) { 89d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,X);CHKERRQ(ierr); 90a683ea4aSPierre Jolivet } else { 91a683ea4aSPierre Jolivet ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr); 92a683ea4aSPierre Jolivet if (transpose) { 93d54c938cSPierre Jolivet ierr = MatTransposeAXPY_Private(Y,a,X,str,Y);CHKERRQ(ierr); 94a683ea4aSPierre Jolivet } else { 95f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 96607cd303SBarry Smith } 97a683ea4aSPierre Jolivet } 98a683ea4aSPierre Jolivet } 99e8136da8SHong Zhang ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 100607cd303SBarry Smith PetscFunctionReturn(0); 101607cd303SBarry Smith } 102607cd303SBarry Smith 103646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B) 104646531bbSStefano Zampini { 105646531bbSStefano Zampini PetscErrorCode ierr; 10627afe9e8SStefano Zampini PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL; 107646531bbSStefano Zampini 108646531bbSStefano Zampini PetscFunctionBegin; 109646531bbSStefano Zampini /* look for any available faster alternative to the general preallocator */ 110646531bbSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr); 111646531bbSStefano Zampini if (preall) { 112646531bbSStefano Zampini ierr = (*preall)(Y,X,B);CHKERRQ(ierr); 113646531bbSStefano Zampini } else { /* Use MatPrellocator, assumes same row-col distribution */ 114646531bbSStefano Zampini Mat preallocator; 115646531bbSStefano Zampini PetscInt r,rstart,rend; 116646531bbSStefano Zampini PetscInt m,n,M,N; 117646531bbSStefano Zampini 1186eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 1196eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 120646531bbSStefano Zampini ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr); 121646531bbSStefano Zampini ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr); 122646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr); 123646531bbSStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 124a6ab3590SBarry Smith ierr = MatSetLayouts(preallocator,Y->rmap,Y->cmap);CHKERRQ(ierr); 125646531bbSStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 126646531bbSStefano Zampini ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr); 127646531bbSStefano Zampini for (r = rstart; r < rend; ++r) { 128646531bbSStefano Zampini PetscInt ncols; 129646531bbSStefano Zampini const PetscInt *row; 130646531bbSStefano Zampini const PetscScalar *vals; 131646531bbSStefano Zampini 132646531bbSStefano Zampini ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 133646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 134646531bbSStefano Zampini ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr); 135646531bbSStefano Zampini ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 136646531bbSStefano Zampini ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr); 137646531bbSStefano Zampini ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr); 138646531bbSStefano Zampini } 139b8550948SPierre Jolivet ierr = MatSetOption(preallocator,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 140646531bbSStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 141646531bbSStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1426eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 1436eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 144646531bbSStefano Zampini 145646531bbSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr); 146646531bbSStefano Zampini ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr); 147a6ab3590SBarry Smith ierr = MatSetLayouts(*B,Y->rmap,Y->cmap);CHKERRQ(ierr); 148646531bbSStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr); 149646531bbSStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 150646531bbSStefano Zampini } 151646531bbSStefano Zampini PetscFunctionReturn(0); 152646531bbSStefano Zampini } 153646531bbSStefano Zampini 154f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 155607cd303SBarry Smith { 1566849ba73SBarry Smith PetscErrorCode ierr; 157be705e3aSPierre Jolivet PetscBool isshell,isdense,isnest; 158edf9a46cSStefano Zampini 159edf9a46cSStefano Zampini PetscFunctionBegin; 160e6aea9c0SStefano Zampini ierr = MatIsShell(Y,&isshell);CHKERRQ(ierr); 161edf9a46cSStefano Zampini if (isshell) { /* MatShell has special support for AXPY */ 162edf9a46cSStefano Zampini PetscErrorCode (*f)(Mat,PetscScalar,Mat,MatStructure); 163edf9a46cSStefano Zampini 164edf9a46cSStefano Zampini ierr = MatGetOperation(Y,MATOP_AXPY,(void (**)(void))&f);CHKERRQ(ierr); 165edf9a46cSStefano Zampini if (f) { 166edf9a46cSStefano Zampini ierr = (*f)(Y,a,X,str);CHKERRQ(ierr); 167edf9a46cSStefano Zampini PetscFunctionReturn(0); 168edf9a46cSStefano Zampini } 169edf9a46cSStefano Zampini } 17093c87e65SStefano Zampini /* no need to preallocate if Y is dense */ 17193c87e65SStefano Zampini ierr = PetscObjectBaseTypeCompareAny((PetscObject)Y,&isdense,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 172be705e3aSPierre Jolivet if (isdense) { 173be705e3aSPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)X,MATNEST,&isnest);CHKERRQ(ierr); 174be705e3aSPierre Jolivet if (isnest) { 175be705e3aSPierre Jolivet ierr = MatAXPY_Dense_Nest(Y,a,X);CHKERRQ(ierr); 176be705e3aSPierre Jolivet PetscFunctionReturn(0); 177be705e3aSPierre Jolivet } 178be705e3aSPierre Jolivet if (str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN) str = SUBSET_NONZERO_PATTERN; 179be705e3aSPierre Jolivet } 180a6ab3590SBarry Smith if (str != DIFFERENT_NONZERO_PATTERN && str != UNKNOWN_NONZERO_PATTERN) { 181edf9a46cSStefano Zampini PetscInt i,start,end,j,ncols,m,n; 18238baddfdSBarry Smith const PetscInt *row; 183b3cc6726SBarry Smith PetscScalar *val; 184b3cc6726SBarry Smith const PetscScalar *vals; 185607cd303SBarry Smith 1868dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 18790f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 1886eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 189f4df32b1SMatthew Knepley if (a == 1.0) { 190d4bb536fSBarry Smith for (i = start; i < end; i++) { 191d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 192d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 193d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 194d4bb536fSBarry Smith } 195d4bb536fSBarry Smith } else { 19621a3365eSStefano Zampini PetscInt vs = 100; 19721a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 19821a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 19906be10caSBarry Smith for (i=start; i<end; i++) { 200b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 20121a3365eSStefano Zampini if (vs < ncols) { 20221a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 20321a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 2046f79c3a4SBarry Smith } 20521a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 206b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 207b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 2086f79c3a4SBarry Smith } 209b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 210d4bb536fSBarry Smith } 2116eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 2126d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2136d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 214edf9a46cSStefano Zampini } else { 215edf9a46cSStefano Zampini Mat B; 216edf9a46cSStefano Zampini 217edf9a46cSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 218edf9a46cSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 21979c2fd05SStefano Zampini ierr = MatHeaderMerge(Y,&B);CHKERRQ(ierr); 220edf9a46cSStefano Zampini } 2213a40ed3dSBarry Smith PetscFunctionReturn(0); 2226f79c3a4SBarry Smith } 223052efed2SBarry Smith 224ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 225ec7775f6SShri Abhyankar { 226ec7775f6SShri Abhyankar PetscInt i,start,end,j,ncols,m,n; 227ec7775f6SShri Abhyankar PetscErrorCode ierr; 228ec7775f6SShri Abhyankar const PetscInt *row; 229ec7775f6SShri Abhyankar PetscScalar *val; 230ec7775f6SShri Abhyankar const PetscScalar *vals; 231ec7775f6SShri Abhyankar 232ec7775f6SShri Abhyankar PetscFunctionBegin; 233ec7775f6SShri Abhyankar ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 234ec7775f6SShri Abhyankar ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 2356eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr); 2366eab9c35SStefano Zampini ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 237ec7775f6SShri Abhyankar if (a == 1.0) { 238ec7775f6SShri Abhyankar for (i = start; i < end; i++) { 239ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 240ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 241ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 242ec7775f6SShri Abhyankar 243ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 244ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 245ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 246ec7775f6SShri Abhyankar } 247ec7775f6SShri Abhyankar } else { 24821a3365eSStefano Zampini PetscInt vs = 100; 24921a3365eSStefano Zampini /* realloc if needed, as this function may be used in parallel */ 25021a3365eSStefano Zampini ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr); 251ec7775f6SShri Abhyankar for (i=start; i<end; i++) { 252ec7775f6SShri Abhyankar ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 253ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 254ec7775f6SShri Abhyankar ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 255ec7775f6SShri Abhyankar 256ec7775f6SShri Abhyankar ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 25721a3365eSStefano Zampini if (vs < ncols) { 25821a3365eSStefano Zampini vs = PetscMin(2*ncols,n); 25921a3365eSStefano Zampini ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr); 260ec7775f6SShri Abhyankar } 26121a3365eSStefano Zampini for (j=0; j<ncols; j++) val[j] = a*vals[j]; 262ec7775f6SShri Abhyankar ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 263ec7775f6SShri Abhyankar ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 264ec7775f6SShri Abhyankar } 265ec7775f6SShri Abhyankar ierr = PetscFree(val);CHKERRQ(ierr); 266ec7775f6SShri Abhyankar } 2676eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr); 2686eab9c35SStefano Zampini ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 269ec7775f6SShri Abhyankar ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 270ec7775f6SShri Abhyankar ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 271ec7775f6SShri Abhyankar PetscFunctionReturn(0); 272ec7775f6SShri Abhyankar } 273ec7775f6SShri Abhyankar 274052efed2SBarry Smith /*@ 27587828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 276052efed2SBarry Smith 2773f9fe445SBarry Smith Neighbor-wise Collective on Mat 278fee21e36SBarry Smith 27998a79cdbSBarry Smith Input Parameters: 28098a79cdbSBarry Smith + Y - the matrices 28187828ca2SBarry Smith - a - the PetscScalar 28298a79cdbSBarry Smith 2832860a424SLois Curfman McInnes Level: intermediate 2842860a424SLois Curfman McInnes 28595452b02SPatrick Sanan Notes: 28695452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 2876f33a894SBarry Smith fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 288edf9a46cSStefano Zampini entry). No operation is performed when a is zero. 2896f33a894SBarry Smith 2900c0fd78eSBarry Smith To form Y = Y + diag(V) use MatDiagonalSet() 2910c0fd78eSBarry Smith 2920c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 293052efed2SBarry Smith @*/ 2947087cfbeSBarry Smith PetscErrorCode MatShift(Mat Y,PetscScalar a) 295052efed2SBarry Smith { 2966849ba73SBarry Smith PetscErrorCode ierr; 297052efed2SBarry Smith 2983a40ed3dSBarry Smith PetscFunctionBegin; 2990700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 300*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!Y->assembled,PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 301*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(Y->factortype,PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 3024994cf47SJed Brown MatCheckPreallocated(Y,1); 3031b71c1a7SBarry Smith if (a == 0.0) PetscFunctionReturn(0); 304b50b34bdSBarry Smith 3051c738b47SStefano Zampini if (Y->ops->shift) { 306f4df32b1SMatthew Knepley ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 3071c738b47SStefano Zampini } else { 3081c738b47SStefano Zampini ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 3091c738b47SStefano Zampini } 3107d68702bSBarry Smith 3115528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 313052efed2SBarry Smith } 3146d84be18SBarry Smith 3157087cfbeSBarry Smith PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 31609f38230SBarry Smith { 31709f38230SBarry Smith PetscErrorCode ierr; 31867576b8bSBarry Smith PetscInt i,start,end; 319fa2e578bSStefano Zampini const PetscScalar *v; 32009f38230SBarry Smith 32109f38230SBarry Smith PetscFunctionBegin; 32209f38230SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 323fa2e578bSStefano Zampini ierr = VecGetArrayRead(D,&v);CHKERRQ(ierr); 32409f38230SBarry Smith for (i=start; i<end; i++) { 32509f38230SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 32609f38230SBarry Smith } 327fa2e578bSStefano Zampini ierr = VecRestoreArrayRead(D,&v);CHKERRQ(ierr); 32809f38230SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 32909f38230SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 33009f38230SBarry Smith PetscFunctionReturn(0); 33109f38230SBarry Smith } 33209f38230SBarry Smith 3336d84be18SBarry Smith /*@ 334f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 335f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 336f56f2b3fSBarry Smith INSERT_VALUES. 3376d84be18SBarry Smith 33848e586daSJose E. Roman Neighbor-wise Collective on Mat 33948e586daSJose E. Roman 3406d84be18SBarry Smith Input Parameters: 34198a79cdbSBarry Smith + Y - the input matrix 342f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 343f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 3446d84be18SBarry Smith 34595452b02SPatrick Sanan Notes: 34695452b02SPatrick Sanan If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 3476f33a894SBarry 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 3486f33a894SBarry Smith entry). 3496f33a894SBarry Smith 3502860a424SLois Curfman McInnes Level: intermediate 3512860a424SLois Curfman McInnes 3520c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale() 3536d84be18SBarry Smith @*/ 3547087cfbeSBarry Smith PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 3556d84be18SBarry Smith { 3566849ba73SBarry Smith PetscErrorCode ierr; 35767576b8bSBarry Smith PetscInt matlocal,veclocal; 3586d84be18SBarry Smith 3593a40ed3dSBarry Smith PetscFunctionBegin; 3600700a824SBarry Smith PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 3610700a824SBarry Smith PetscValidHeaderSpecific(D,VEC_CLASSID,2); 36267576b8bSBarry Smith ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 36367576b8bSBarry Smith ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 364*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(matlocal != veclocal,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number local rows of matrix %" PetscInt_FMT " does not match that of vector for diagonal %" PetscInt_FMT,matlocal,veclocal); 365f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 366f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 36794d884c6SBarry Smith } else { 36809f38230SBarry Smith ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 3696d84be18SBarry Smith } 3705528ad4fSTristan Konolige ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 3713a40ed3dSBarry Smith PetscFunctionReturn(0); 3726d84be18SBarry Smith } 373d4bb536fSBarry Smith 374d4bb536fSBarry Smith /*@ 37504aac2b0SHong Zhang MatAYPX - Computes Y = a*Y + X. 376d4bb536fSBarry Smith 3773f9fe445SBarry Smith Logically on Mat 378fee21e36SBarry Smith 37998a79cdbSBarry Smith Input Parameters: 38004aac2b0SHong Zhang + a - the PetscScalar multiplier 38104aac2b0SHong Zhang . Y - the first matrix 38204aac2b0SHong Zhang . X - the second matrix 383531d0c6aSBarry Smith - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, UNKNOWN_NONZERO_PATTERN, or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 384d4bb536fSBarry Smith 3852860a424SLois Curfman McInnes Level: intermediate 3862860a424SLois Curfman McInnes 3872860a424SLois Curfman McInnes .seealso: MatAXPY() 388d4bb536fSBarry Smith @*/ 3897087cfbeSBarry Smith PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 390d4bb536fSBarry Smith { 3916849ba73SBarry Smith PetscErrorCode ierr; 392d4bb536fSBarry Smith 3933a40ed3dSBarry Smith PetscFunctionBegin; 394f4df32b1SMatthew Knepley ierr = MatScale(Y,a);CHKERRQ(ierr); 3952ab6f6a8SStefano Zampini ierr = MatAXPY(Y,1.0,X,str);CHKERRQ(ierr); 3963a40ed3dSBarry Smith PetscFunctionReturn(0); 397d4bb536fSBarry Smith } 398b0a32e0cSBarry Smith 399b0a32e0cSBarry Smith /*@ 4000bacdadaSStefano Zampini MatComputeOperator - Computes the explicit matrix 401b0a32e0cSBarry Smith 402b0a32e0cSBarry Smith Collective on Mat 403b0a32e0cSBarry Smith 404d8d19677SJose E. Roman Input Parameters: 405186a3e20SStefano Zampini + inmat - the matrix 406186a3e20SStefano Zampini - mattype - the matrix type for the explicit operator 407b0a32e0cSBarry Smith 408b0a32e0cSBarry Smith Output Parameter: 409a5b23f4aSJose E. Roman . mat - the explicit operator 410b0a32e0cSBarry Smith 411b0a32e0cSBarry Smith Notes: 412186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 413186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 414186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 415b0a32e0cSBarry Smith 416b0a32e0cSBarry Smith Level: advanced 417b0a32e0cSBarry Smith 418b0a32e0cSBarry Smith @*/ 4190bacdadaSStefano Zampini PetscErrorCode MatComputeOperator(Mat inmat,MatType mattype,Mat *mat) 420b0a32e0cSBarry Smith { 421dfbe8321SBarry Smith PetscErrorCode ierr; 422b0a32e0cSBarry Smith 423b0a32e0cSBarry Smith PetscFunctionBegin; 4240700a824SBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 425186a3e20SStefano Zampini PetscValidPointer(mat,3); 426186a3e20SStefano Zampini ierr = MatConvert_Shell(inmat,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 42724f910e3SHong Zhang PetscFunctionReturn(0); 42824f910e3SHong Zhang } 4294325cce7SMatthew G Knepley 4304325cce7SMatthew G Knepley /*@ 4310bacdadaSStefano Zampini MatComputeOperatorTranspose - Computes the explicit matrix representation of 432f3b1f45cSBarry Smith a give matrix that can apply MatMultTranspose() 433f3b1f45cSBarry Smith 434f3b1f45cSBarry Smith Collective on Mat 435f3b1f45cSBarry Smith 4366b867d5aSJose E. Roman Input Parameters: 4376b867d5aSJose E. Roman + inmat - the matrix 4386b867d5aSJose E. Roman - mattype - the matrix type for the explicit operator 439f3b1f45cSBarry Smith 440f3b1f45cSBarry Smith Output Parameter: 441a5b23f4aSJose E. Roman . mat - the explicit operator transposed 442f3b1f45cSBarry Smith 443f3b1f45cSBarry Smith Notes: 444186a3e20SStefano Zampini This computation is done by applying the transpose of the operator to columns of the identity matrix. 445186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 446186a3e20SStefano Zampini Currently, this routine uses a dense matrix format if mattype == NULL. 447f3b1f45cSBarry Smith 448f3b1f45cSBarry Smith Level: advanced 449f3b1f45cSBarry Smith 450f3b1f45cSBarry Smith @*/ 4510bacdadaSStefano Zampini PetscErrorCode MatComputeOperatorTranspose(Mat inmat,MatType mattype,Mat *mat) 452f3b1f45cSBarry Smith { 453186a3e20SStefano Zampini Mat A; 454f3b1f45cSBarry Smith PetscErrorCode ierr; 455f3b1f45cSBarry Smith 456f3b1f45cSBarry Smith PetscFunctionBegin; 457f3b1f45cSBarry Smith PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 458186a3e20SStefano Zampini PetscValidPointer(mat,3); 459186a3e20SStefano Zampini ierr = MatCreateTranspose(inmat,&A);CHKERRQ(ierr); 460186a3e20SStefano Zampini ierr = MatConvert_Shell(A,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 461186a3e20SStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 462f3b1f45cSBarry Smith PetscFunctionReturn(0); 463f3b1f45cSBarry Smith } 464f3b1f45cSBarry Smith 465f3b1f45cSBarry Smith /*@ 4664325cce7SMatthew G Knepley MatChop - Set all values in the matrix less than the tolerance to zero 4674325cce7SMatthew G Knepley 4684325cce7SMatthew G Knepley Input Parameters: 4694325cce7SMatthew G Knepley + A - The matrix 4704325cce7SMatthew G Knepley - tol - The zero tolerance 4714325cce7SMatthew G Knepley 4724325cce7SMatthew G Knepley Output Parameters: 4734325cce7SMatthew G Knepley . A - The chopped matrix 4744325cce7SMatthew G Knepley 4754325cce7SMatthew G Knepley Level: intermediate 4764325cce7SMatthew G Knepley 4774325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries() 4783fc99919SSatish Balay @*/ 4794325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol) 4804325cce7SMatthew G Knepley { 481038df967SPierre Jolivet Mat a; 4824325cce7SMatthew G Knepley PetscScalar *newVals; 483cf5d3338SPierre Jolivet PetscInt *newCols, rStart, rEnd, numRows, maxRows, r, colMax = 0; 484038df967SPierre Jolivet PetscBool flg; 4854325cce7SMatthew G Knepley PetscErrorCode ierr; 4864325cce7SMatthew G Knepley 4874325cce7SMatthew G Knepley PetscFunctionBegin; 488038df967SPierre Jolivet ierr = PetscObjectBaseTypeCompareAny((PetscObject)A, &flg, MATSEQDENSE, MATMPIDENSE, "");CHKERRQ(ierr); 489038df967SPierre Jolivet if (flg) { 490038df967SPierre Jolivet ierr = MatDenseGetLocalMatrix(A, &a);CHKERRQ(ierr); 491038df967SPierre Jolivet ierr = MatDenseGetLDA(a, &r);CHKERRQ(ierr); 492038df967SPierre Jolivet ierr = MatGetSize(a, &rStart, &rEnd);CHKERRQ(ierr); 493038df967SPierre Jolivet ierr = MatDenseGetArray(a, &newVals);CHKERRQ(ierr); 494038df967SPierre Jolivet for (; colMax < rEnd; ++colMax) { 495038df967SPierre Jolivet for (maxRows = 0; maxRows < rStart; ++maxRows) { 496038df967SPierre Jolivet newVals[maxRows + colMax * r] = PetscAbsScalar(newVals[maxRows + colMax * r]) < tol ? 0.0 : newVals[maxRows + colMax * r]; 497038df967SPierre Jolivet } 498038df967SPierre Jolivet } 499038df967SPierre Jolivet ierr = MatDenseRestoreArray(a, &newVals);CHKERRQ(ierr); 500038df967SPierre Jolivet } else { 5014325cce7SMatthew G Knepley ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 502a4c05331SJacob Faibussowitsch ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 5034325cce7SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 5044325cce7SMatthew G Knepley PetscInt ncols; 5054325cce7SMatthew G Knepley 5060298fd71SBarry Smith ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5074325cce7SMatthew G Knepley colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 5080298fd71SBarry Smith ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 5094325cce7SMatthew G Knepley } 5104325cce7SMatthew G Knepley numRows = rEnd - rStart; 511820f2d46SBarry Smith ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 512dcca6d9dSJed Brown ierr = PetscMalloc2(colMax, &newCols, colMax, &newVals);CHKERRQ(ierr); 513cf5d3338SPierre Jolivet ierr = MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &flg);CHKERRQ(ierr); /* cache user-defined value */ 514038df967SPierre Jolivet ierr = MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);CHKERRQ(ierr); 515cf5d3338SPierre Jolivet /* short-circuit code in MatAssemblyBegin() and MatAssemblyEnd() */ 516cf5d3338SPierre Jolivet /* that are potentially called many times depending on the distribution of A */ 5174325cce7SMatthew G Knepley for (r = rStart; r < rStart+maxRows; ++r) { 5184325cce7SMatthew G Knepley const PetscScalar *vals; 5194325cce7SMatthew G Knepley const PetscInt *cols; 520fad45679SMatthew G. Knepley PetscInt ncols, newcols, c; 5214325cce7SMatthew G Knepley 5224325cce7SMatthew G Knepley if (r < rEnd) { 5234325cce7SMatthew G Knepley ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 5244325cce7SMatthew G Knepley for (c = 0; c < ncols; ++c) { 5254325cce7SMatthew G Knepley newCols[c] = cols[c]; 5264325cce7SMatthew G Knepley newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 5274325cce7SMatthew G Knepley } 528fad45679SMatthew G. Knepley newcols = ncols; 5294325cce7SMatthew G Knepley ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 530fad45679SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 5314325cce7SMatthew G Knepley } 5324325cce7SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5334325cce7SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5344325cce7SMatthew G Knepley } 535a4c05331SJacob Faibussowitsch ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 5364325cce7SMatthew G Knepley ierr = PetscFree2(newCols, newVals);CHKERRQ(ierr); 537cf5d3338SPierre Jolivet ierr = MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, flg);CHKERRQ(ierr); /* reset option to its user-defined value */ 538038df967SPierre Jolivet } 5394325cce7SMatthew G Knepley PetscFunctionReturn(0); 5404325cce7SMatthew G Knepley } 541