xref: /petsc/src/mat/utils/axpy.c (revision 1b71c1a785b8a166f8741a9698fe876d03e19e81)
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 
1266eab9c35SStefano Zampini     ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr);
1276eab9c35SStefano Zampini     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
128646531bbSStefano Zampini     ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr);
129646531bbSStefano Zampini     ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr);
130646531bbSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr);
131646531bbSStefano Zampini     ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr);
132646531bbSStefano Zampini     ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr);
133646531bbSStefano Zampini     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
134646531bbSStefano Zampini     ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr);
135646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
136646531bbSStefano Zampini       PetscInt          ncols;
137646531bbSStefano Zampini       const PetscInt    *row;
138646531bbSStefano Zampini       const PetscScalar *vals;
139646531bbSStefano Zampini 
140646531bbSStefano Zampini       ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
141646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
142646531bbSStefano Zampini       ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
143646531bbSStefano Zampini       ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
144646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
145646531bbSStefano Zampini       ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
146646531bbSStefano Zampini     }
147646531bbSStefano Zampini     ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148646531bbSStefano Zampini     ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1496eab9c35SStefano Zampini     ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr);
1506eab9c35SStefano Zampini     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
151646531bbSStefano Zampini 
152646531bbSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr);
153646531bbSStefano Zampini     ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr);
154646531bbSStefano Zampini     ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr);
155646531bbSStefano Zampini     ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr);
156646531bbSStefano Zampini     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
157646531bbSStefano Zampini   }
158646531bbSStefano Zampini   PetscFunctionReturn(0);
159646531bbSStefano Zampini }
160646531bbSStefano Zampini 
161f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
162607cd303SBarry Smith {
16338baddfdSBarry Smith   PetscInt          i,start,end,j,ncols,m,n;
1646849ba73SBarry Smith   PetscErrorCode    ierr;
16538baddfdSBarry Smith   const PetscInt    *row;
166b3cc6726SBarry Smith   PetscScalar       *val;
167b3cc6726SBarry Smith   const PetscScalar *vals;
168607cd303SBarry Smith 
169607cd303SBarry Smith   PetscFunctionBegin;
1708dadbd76SSatish Balay   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
17190f02eecSBarry Smith   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
1726eab9c35SStefano Zampini   ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
173f4df32b1SMatthew Knepley   if (a == 1.0) {
174d4bb536fSBarry Smith     for (i = start; i < end; i++) {
175d4bb536fSBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
176d4bb536fSBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
177d4bb536fSBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
178d4bb536fSBarry Smith     }
179d4bb536fSBarry Smith   } else {
18021a3365eSStefano Zampini     PetscInt vs = 100;
18121a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
18221a3365eSStefano Zampini     ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
18306be10caSBarry Smith     for (i=start; i<end; i++) {
184b3cc6726SBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
18521a3365eSStefano Zampini       if (vs < ncols) {
18621a3365eSStefano Zampini         vs   = PetscMin(2*ncols,n);
18721a3365eSStefano Zampini         ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr);
1886f79c3a4SBarry Smith       }
18921a3365eSStefano Zampini       for (j=0; j<ncols; j++) val[j] = a*vals[j];
190b3cc6726SBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
191b3cc6726SBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
1926f79c3a4SBarry Smith     }
193b3cc6726SBarry Smith     ierr = PetscFree(val);CHKERRQ(ierr);
194d4bb536fSBarry Smith   }
1956eab9c35SStefano Zampini   ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
1966d4a8577SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1976d4a8577SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1996f79c3a4SBarry Smith }
200052efed2SBarry Smith 
201ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
202ec7775f6SShri Abhyankar {
203ec7775f6SShri Abhyankar   PetscInt          i,start,end,j,ncols,m,n;
204ec7775f6SShri Abhyankar   PetscErrorCode    ierr;
205ec7775f6SShri Abhyankar   const PetscInt    *row;
206ec7775f6SShri Abhyankar   PetscScalar       *val;
207ec7775f6SShri Abhyankar   const PetscScalar *vals;
208ec7775f6SShri Abhyankar 
209ec7775f6SShri Abhyankar   PetscFunctionBegin;
210ec7775f6SShri Abhyankar   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
211ec7775f6SShri Abhyankar   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
2126eab9c35SStefano Zampini   ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr);
2136eab9c35SStefano Zampini   ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
214ec7775f6SShri Abhyankar   if (a == 1.0) {
215ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
216ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
217ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
218ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
219ec7775f6SShri Abhyankar 
220ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
221ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
222ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
223ec7775f6SShri Abhyankar     }
224ec7775f6SShri Abhyankar   } else {
22521a3365eSStefano Zampini     PetscInt vs = 100;
22621a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
22721a3365eSStefano Zampini     ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
228ec7775f6SShri Abhyankar     for (i=start; i<end; i++) {
229ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
230ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
231ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
232ec7775f6SShri Abhyankar 
233ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
23421a3365eSStefano Zampini       if (vs < ncols) {
23521a3365eSStefano Zampini         vs   = PetscMin(2*ncols,n);
23621a3365eSStefano Zampini         ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr);
237ec7775f6SShri Abhyankar       }
23821a3365eSStefano Zampini       for (j=0; j<ncols; j++) val[j] = a*vals[j];
239ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
240ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
241ec7775f6SShri Abhyankar     }
242ec7775f6SShri Abhyankar     ierr = PetscFree(val);CHKERRQ(ierr);
243ec7775f6SShri Abhyankar   }
2446eab9c35SStefano Zampini   ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr);
2456eab9c35SStefano Zampini   ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
246ec7775f6SShri Abhyankar   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247ec7775f6SShri Abhyankar   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248ec7775f6SShri Abhyankar   PetscFunctionReturn(0);
249ec7775f6SShri Abhyankar }
250ec7775f6SShri Abhyankar 
251052efed2SBarry Smith /*@
25287828ca2SBarry Smith    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
253052efed2SBarry Smith 
2543f9fe445SBarry Smith    Neighbor-wise Collective on Mat
255fee21e36SBarry Smith 
25698a79cdbSBarry Smith    Input Parameters:
25798a79cdbSBarry Smith +  Y - the matrices
25887828ca2SBarry Smith -  a - the PetscScalar
25998a79cdbSBarry Smith 
2602860a424SLois Curfman McInnes    Level: intermediate
2612860a424SLois Curfman McInnes 
26295452b02SPatrick Sanan    Notes:
26395452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2646f33a894SBarry 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
2656f33a894SBarry Smith    entry).
2666f33a894SBarry Smith 
2670c0fd78eSBarry Smith    To form Y = Y + diag(V) use MatDiagonalSet()
2680c0fd78eSBarry Smith 
2696f33a894SBarry Smith    Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is
2706f33a894SBarry Smith     preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix.
2716f33a894SBarry Smith 
272052efed2SBarry Smith .keywords: matrix, add, shift
2736b9ee512SLois Curfman McInnes 
2740c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale()
275052efed2SBarry Smith  @*/
2767087cfbeSBarry Smith PetscErrorCode  MatShift(Mat Y,PetscScalar a)
277052efed2SBarry Smith {
2786849ba73SBarry Smith   PetscErrorCode ierr;
279052efed2SBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
2810700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
282ce94432eSBarry Smith   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
283ce94432eSBarry Smith   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2844994cf47SJed Brown   MatCheckPreallocated(Y,1);
285*1b71c1a7SBarry Smith   if (a == 0.0) PetscFunctionReturn(0);
286b50b34bdSBarry Smith 
2871c738b47SStefano Zampini   if (Y->ops->shift) {
288f4df32b1SMatthew Knepley     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
2891c738b47SStefano Zampini   } else {
2901c738b47SStefano Zampini     ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
2911c738b47SStefano Zampini   }
2927d68702bSBarry Smith 
2935528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
294fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
295b8ced49eSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
296b8ced49eSKarl Rupp     Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
297c41cb2e2SAlejandro Lamas Daviña   }
298d67ff14aSKarl Rupp #endif
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
300052efed2SBarry Smith }
3016d84be18SBarry Smith 
3027087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
30309f38230SBarry Smith {
30409f38230SBarry Smith   PetscErrorCode ierr;
30567576b8bSBarry Smith   PetscInt       i,start,end;
30609f38230SBarry Smith   PetscScalar    *v;
30709f38230SBarry Smith 
30809f38230SBarry Smith   PetscFunctionBegin;
30909f38230SBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
31009f38230SBarry Smith   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
31109f38230SBarry Smith   for (i=start; i<end; i++) {
31209f38230SBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
31309f38230SBarry Smith   }
31409f38230SBarry Smith   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
31509f38230SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
31609f38230SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
31709f38230SBarry Smith   PetscFunctionReturn(0);
31809f38230SBarry Smith }
31909f38230SBarry Smith 
3206d84be18SBarry Smith /*@
321f56f2b3fSBarry Smith    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
322f56f2b3fSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
323f56f2b3fSBarry Smith    INSERT_VALUES.
3246d84be18SBarry Smith 
3256d84be18SBarry Smith    Input Parameters:
32698a79cdbSBarry Smith +  Y - the input matrix
327f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
328f56f2b3fSBarry Smith -  i - INSERT_VALUES or ADD_VALUES
3296d84be18SBarry Smith 
3303f9fe445SBarry Smith    Neighbor-wise Collective on Mat and Vec
331fee21e36SBarry Smith 
33295452b02SPatrick Sanan    Notes:
33395452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3346f33a894SBarry 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
3356f33a894SBarry Smith    entry).
3366f33a894SBarry Smith 
3372860a424SLois Curfman McInnes    Level: intermediate
3382860a424SLois Curfman McInnes 
3396b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal
3406b9ee512SLois Curfman McInnes 
3410c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale()
3426d84be18SBarry Smith @*/
3437087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
3446d84be18SBarry Smith {
3456849ba73SBarry Smith   PetscErrorCode ierr;
34667576b8bSBarry Smith   PetscInt       matlocal,veclocal;
3476d84be18SBarry Smith 
3483a40ed3dSBarry Smith   PetscFunctionBegin;
3490700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
3500700a824SBarry Smith   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
35167576b8bSBarry Smith   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
35267576b8bSBarry Smith   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
35367576b8bSBarry 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);
354f56f2b3fSBarry Smith   if (Y->ops->diagonalset) {
355f56f2b3fSBarry Smith     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
35694d884c6SBarry Smith   } else {
35709f38230SBarry Smith     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
3586d84be18SBarry Smith   }
3595528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
3603a40ed3dSBarry Smith   PetscFunctionReturn(0);
3616d84be18SBarry Smith }
362d4bb536fSBarry Smith 
363d4bb536fSBarry Smith /*@
36404aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
365d4bb536fSBarry Smith 
3663f9fe445SBarry Smith    Logically on Mat
367fee21e36SBarry Smith 
36898a79cdbSBarry Smith    Input Parameters:
36904aac2b0SHong Zhang +  a - the PetscScalar multiplier
37004aac2b0SHong Zhang .  Y - the first matrix
37104aac2b0SHong Zhang .  X - the second matrix
37204aac2b0SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
373d4bb536fSBarry Smith 
3742860a424SLois Curfman McInnes    Level: intermediate
3752860a424SLois Curfman McInnes 
376d4bb536fSBarry Smith .keywords: matrix, add
377d4bb536fSBarry Smith 
3782860a424SLois Curfman McInnes .seealso: MatAXPY()
379d4bb536fSBarry Smith  @*/
3807087cfbeSBarry Smith PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
381d4bb536fSBarry Smith {
38287828ca2SBarry Smith   PetscScalar    one = 1.0;
3836849ba73SBarry Smith   PetscErrorCode ierr;
38438baddfdSBarry Smith   PetscInt       mX,mY,nX,nY;
385d4bb536fSBarry Smith 
3863a40ed3dSBarry Smith   PetscFunctionBegin;
387c5eb9154SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
3880700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
389c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
390329f5518SBarry Smith   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
391329f5518SBarry Smith   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
392e32f2f54SBarry 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);
393f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
394cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
3953a40ed3dSBarry Smith   PetscFunctionReturn(0);
396d4bb536fSBarry Smith }
397b0a32e0cSBarry Smith 
398b0a32e0cSBarry Smith /*@
3990bacdadaSStefano Zampini     MatComputeOperator - Computes the explicit matrix
400b0a32e0cSBarry Smith 
401b0a32e0cSBarry Smith     Collective on Mat
402b0a32e0cSBarry Smith 
403b0a32e0cSBarry Smith     Input Parameter:
404186a3e20SStefano Zampini +   inmat - the matrix
405186a3e20SStefano Zampini -   mattype - the matrix type for the explicit operator
406b0a32e0cSBarry Smith 
407b0a32e0cSBarry Smith     Output Parameter:
408f3b1f45cSBarry Smith .   mat - the explict  operator
409b0a32e0cSBarry Smith 
410b0a32e0cSBarry Smith     Notes:
411186a3e20SStefano Zampini     This computation is done by applying the operators to columns of the identity matrix.
412186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
413186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
414b0a32e0cSBarry Smith 
415b0a32e0cSBarry Smith     Level: advanced
416b0a32e0cSBarry Smith 
417b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator
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 
436f3b1f45cSBarry Smith     Input Parameter:
437f3b1f45cSBarry Smith .   inmat - the matrix
438f3b1f45cSBarry Smith 
439f3b1f45cSBarry Smith     Output Parameter:
440f3b1f45cSBarry Smith .   mat - the explict  operator transposed
441f3b1f45cSBarry Smith 
442f3b1f45cSBarry Smith     Notes:
443186a3e20SStefano Zampini     This computation is done by applying the transpose of the operator to columns of the identity matrix.
444186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
445186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
446f3b1f45cSBarry Smith 
447f3b1f45cSBarry Smith     Level: advanced
448f3b1f45cSBarry Smith 
449f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator
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 {
4814325cce7SMatthew G Knepley   PetscScalar    *newVals;
4824325cce7SMatthew G Knepley   PetscInt       *newCols;
4834325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
4844325cce7SMatthew G Knepley   PetscErrorCode ierr;
4854325cce7SMatthew G Knepley 
4864325cce7SMatthew G Knepley   PetscFunctionBegin;
4874325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4884325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4894325cce7SMatthew G Knepley     PetscInt ncols;
4904325cce7SMatthew G Knepley 
4910298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4924325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4930298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4944325cce7SMatthew G Knepley   }
4954325cce7SMatthew G Knepley   numRows = rEnd - rStart;
496b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
497dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4984325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4994325cce7SMatthew G Knepley     const PetscScalar *vals;
5004325cce7SMatthew G Knepley     const PetscInt    *cols;
501fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
5024325cce7SMatthew G Knepley 
5034325cce7SMatthew G Knepley     if (r < rEnd) {
5044325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
5054325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
5064325cce7SMatthew G Knepley         newCols[c] = cols[c];
5074325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
5084325cce7SMatthew G Knepley       }
509fad45679SMatthew G. Knepley       newcols = ncols;
5104325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
511fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
5124325cce7SMatthew G Knepley     }
5134325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5144325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5154325cce7SMatthew G Knepley   }
5164325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
5174325cce7SMatthew G Knepley   PetscFunctionReturn(0);
5184325cce7SMatthew G Knepley }
519