xref: /petsc/src/mat/utils/axpy.c (revision 186a3e20b6c0156a6236f097650d994d1bcd0475)
16f79c3a4SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/
36f79c3a4SBarry Smith 
4d54c938cSPierre Jolivet static PetscErrorCode MatTransposeAXPY_Private(Mat Y,PetscScalar a,Mat X,MatStructure str,Mat T)
5d54c938cSPierre Jolivet {
6d54c938cSPierre Jolivet   PetscErrorCode ierr,(*f)(Mat,Mat*);
7d54c938cSPierre Jolivet   Mat            A,F;
8d54c938cSPierre Jolivet 
9d54c938cSPierre Jolivet   PetscFunctionBegin;
10d54c938cSPierre Jolivet   ierr = PetscObjectQueryFunction((PetscObject)T,"MatTransposeGetMat_C",&f);CHKERRQ(ierr);
11d54c938cSPierre Jolivet   if (f) {
12d54c938cSPierre Jolivet     ierr = MatTransposeGetMat(T,&A);CHKERRQ(ierr);
13d54c938cSPierre Jolivet     if (T == X) {
14d54c938cSPierre Jolivet       ierr = PetscInfo(NULL,"Explicitly transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr);
15d54c938cSPierre Jolivet       ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr);
16d54c938cSPierre Jolivet       A = Y;
17d54c938cSPierre Jolivet     } else {
18d54c938cSPierre Jolivet       ierr = PetscInfo(NULL,"Transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr);
19d54c938cSPierre Jolivet       ierr = MatTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr);
20d54c938cSPierre Jolivet     }
21d54c938cSPierre Jolivet   } else {
22d54c938cSPierre Jolivet     ierr = MatHermitianTransposeGetMat(T,&A);CHKERRQ(ierr);
23d54c938cSPierre Jolivet     if (T == X) {
24d54c938cSPierre Jolivet       ierr = PetscInfo(NULL,"Explicitly Hermitian transposing X of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr);
25d54c938cSPierre Jolivet       ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr);
26d54c938cSPierre Jolivet       A = Y;
27d54c938cSPierre Jolivet     } else {
28d54c938cSPierre Jolivet       ierr = PetscInfo(NULL,"Hermitian transposing X because Y of type MATTRANSPOSEMAT to perform MatAXPY()\n");CHKERRQ(ierr);
29d54c938cSPierre Jolivet       ierr = MatHermitianTranspose(X,MAT_INITIAL_MATRIX,&F);CHKERRQ(ierr);
30d54c938cSPierre Jolivet     }
31d54c938cSPierre Jolivet   }
32d54c938cSPierre Jolivet   ierr = MatAXPY(A,a,F,str);CHKERRQ(ierr);
33d54c938cSPierre Jolivet   ierr = MatDestroy(&F);CHKERRQ(ierr);
34d54c938cSPierre Jolivet   PetscFunctionReturn(0);
35d54c938cSPierre Jolivet }
36d54c938cSPierre Jolivet 
3706be10caSBarry Smith /*@
3821c89e3eSBarry Smith    MatAXPY - Computes Y = a*X + Y.
396f79c3a4SBarry Smith 
403f9fe445SBarry Smith    Logically  Collective on Mat
41fee21e36SBarry Smith 
4298a79cdbSBarry Smith    Input Parameters:
43607cd303SBarry Smith +  a - the scalar multiplier
44607cd303SBarry Smith .  X - the first matrix
45607cd303SBarry Smith .  Y - the second matrix
46407f6b05SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN
47407f6b05SHong Zhang          or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
4898a79cdbSBarry Smith 
492860a424SLois Curfman McInnes    Level: intermediate
502860a424SLois Curfman McInnes 
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 {
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);
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) {
84d54c938cSPierre Jolivet         ierr = MatTransposeAXPY_Private(Y,a,X,str,X);CHKERRQ(ierr);
85a683ea4aSPierre Jolivet     } else {
86a683ea4aSPierre Jolivet       ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr);
87a683ea4aSPierre Jolivet       if (transpose) {
88d54c938cSPierre Jolivet         ierr = MatTransposeAXPY_Private(Y,a,X,str,Y);CHKERRQ(ierr);
89a683ea4aSPierre Jolivet       } else {
90646531bbSStefano Zampini         if (str != DIFFERENT_NONZERO_PATTERN) {
91f4df32b1SMatthew Knepley           ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
92646531bbSStefano Zampini         } else {
93646531bbSStefano Zampini           Mat B;
94646531bbSStefano Zampini 
95646531bbSStefano Zampini           ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr);
96646531bbSStefano Zampini           ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr);
97646531bbSStefano Zampini           ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr);
98646531bbSStefano Zampini         }
99607cd303SBarry Smith       }
100a683ea4aSPierre Jolivet     }
101a683ea4aSPierre Jolivet   }
102e8136da8SHong Zhang   ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
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);
382f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
383cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
3843a40ed3dSBarry Smith   PetscFunctionReturn(0);
385d4bb536fSBarry Smith }
386b0a32e0cSBarry Smith 
387b0a32e0cSBarry Smith /*@
388b0a32e0cSBarry Smith     MatComputeExplicitOperator - Computes the explicit matrix
389b0a32e0cSBarry Smith 
390b0a32e0cSBarry Smith     Collective on Mat
391b0a32e0cSBarry Smith 
392b0a32e0cSBarry Smith     Input Parameter:
393*186a3e20SStefano Zampini +   inmat - the matrix
394*186a3e20SStefano Zampini -   mattype - the matrix type for the explicit operator
395b0a32e0cSBarry Smith 
396b0a32e0cSBarry Smith     Output Parameter:
397f3b1f45cSBarry Smith .   mat - the explict  operator
398b0a32e0cSBarry Smith 
399b0a32e0cSBarry Smith     Notes:
400*186a3e20SStefano Zampini     This computation is done by applying the operators to columns of the identity matrix.
401*186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
402*186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
403b0a32e0cSBarry Smith 
404b0a32e0cSBarry Smith     Level: advanced
405b0a32e0cSBarry Smith 
406b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator
407b0a32e0cSBarry Smith @*/
408*186a3e20SStefano Zampini PetscErrorCode  MatComputeExplicitOperator(Mat inmat,MatType mattype,Mat *mat)
409b0a32e0cSBarry Smith {
410dfbe8321SBarry Smith   PetscErrorCode ierr;
411b0a32e0cSBarry Smith 
412b0a32e0cSBarry Smith   PetscFunctionBegin;
4130700a824SBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
414*186a3e20SStefano Zampini   PetscValidPointer(mat,3);
415*186a3e20SStefano Zampini   ierr = MatConvert_Shell(inmat,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
41624f910e3SHong Zhang   PetscFunctionReturn(0);
41724f910e3SHong Zhang }
4184325cce7SMatthew G Knepley 
4194325cce7SMatthew G Knepley /*@
420f3b1f45cSBarry Smith     MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of
421f3b1f45cSBarry Smith         a give matrix that can apply MatMultTranspose()
422f3b1f45cSBarry Smith 
423f3b1f45cSBarry Smith     Collective on Mat
424f3b1f45cSBarry Smith 
425f3b1f45cSBarry Smith     Input Parameter:
426f3b1f45cSBarry Smith .   inmat - the matrix
427f3b1f45cSBarry Smith 
428f3b1f45cSBarry Smith     Output Parameter:
429f3b1f45cSBarry Smith .   mat - the explict  operator transposed
430f3b1f45cSBarry Smith 
431f3b1f45cSBarry Smith     Notes:
432*186a3e20SStefano Zampini     This computation is done by applying the transpose of the operator to columns of the identity matrix.
433*186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
434*186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
435f3b1f45cSBarry Smith 
436f3b1f45cSBarry Smith     Level: advanced
437f3b1f45cSBarry Smith 
438f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator
439f3b1f45cSBarry Smith @*/
440*186a3e20SStefano Zampini PetscErrorCode  MatComputeExplicitOperatorTranspose(Mat inmat,MatType mattype,Mat *mat)
441f3b1f45cSBarry Smith {
442*186a3e20SStefano Zampini   Mat            A;
443f3b1f45cSBarry Smith   PetscErrorCode ierr;
444f3b1f45cSBarry Smith 
445f3b1f45cSBarry Smith   PetscFunctionBegin;
446f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
447*186a3e20SStefano Zampini   PetscValidPointer(mat,3);
448*186a3e20SStefano Zampini   ierr = MatCreateTranspose(inmat,&A);CHKERRQ(ierr);
449*186a3e20SStefano Zampini   ierr = MatConvert_Shell(A,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
450*186a3e20SStefano Zampini   ierr = MatDestroy(&A);CHKERRQ(ierr);
451f3b1f45cSBarry Smith   PetscFunctionReturn(0);
452f3b1f45cSBarry Smith }
453f3b1f45cSBarry Smith 
454f3b1f45cSBarry Smith /*@
4554325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
4564325cce7SMatthew G Knepley 
4574325cce7SMatthew G Knepley   Input Parameters:
4584325cce7SMatthew G Knepley + A   - The matrix
4594325cce7SMatthew G Knepley - tol - The zero tolerance
4604325cce7SMatthew G Knepley 
4614325cce7SMatthew G Knepley   Output Parameters:
4624325cce7SMatthew G Knepley . A - The chopped matrix
4634325cce7SMatthew G Knepley 
4644325cce7SMatthew G Knepley   Level: intermediate
4654325cce7SMatthew G Knepley 
4664325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries()
4673fc99919SSatish Balay  @*/
4684325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol)
4694325cce7SMatthew G Knepley {
4704325cce7SMatthew G Knepley   PetscScalar    *newVals;
4714325cce7SMatthew G Knepley   PetscInt       *newCols;
4724325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
4734325cce7SMatthew G Knepley   PetscErrorCode ierr;
4744325cce7SMatthew G Knepley 
4754325cce7SMatthew G Knepley   PetscFunctionBegin;
4764325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4774325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4784325cce7SMatthew G Knepley     PetscInt ncols;
4794325cce7SMatthew G Knepley 
4800298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4814325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4820298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4834325cce7SMatthew G Knepley   }
4844325cce7SMatthew G Knepley   numRows = rEnd - rStart;
485b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
486dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4874325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4884325cce7SMatthew G Knepley     const PetscScalar *vals;
4894325cce7SMatthew G Knepley     const PetscInt    *cols;
490fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
4914325cce7SMatthew G Knepley 
4924325cce7SMatthew G Knepley     if (r < rEnd) {
4934325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
4944325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
4954325cce7SMatthew G Knepley         newCols[c] = cols[c];
4964325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
4974325cce7SMatthew G Knepley       }
498fad45679SMatthew G. Knepley       newcols = ncols;
4994325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
500fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
5014325cce7SMatthew G Knepley     }
5024325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5034325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5044325cce7SMatthew G Knepley   }
5054325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
5064325cce7SMatthew G Knepley   PetscFunctionReturn(0);
5074325cce7SMatthew G Knepley }
508