xref: /petsc/src/mat/utils/axpy.c (revision edf9a46c28dade19307bcf4eebf8bfc3d44ed2f7)
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 
49*edf9a46cSStefano Zampini    Notes: No operation is performed when a is zero.
50*edf9a46cSStefano Zampini 
512860a424SLois Curfman McInnes    Level: intermediate
522860a424SLois Curfman McInnes 
532860a424SLois Curfman McInnes .seealso: MatAYPX()
5406be10caSBarry Smith  @*/
557087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
566f79c3a4SBarry Smith {
57d54c938cSPierre Jolivet   PetscErrorCode ierr;
58646531bbSStefano Zampini   PetscInt       M1,M2,N1,N2;
59c1ac3661SBarry Smith   PetscInt       m1,m2,n1,n2;
60a683ea4aSPierre Jolivet   MatType        t1,t2;
61a683ea4aSPierre Jolivet   PetscBool      sametype,transpose;
626f79c3a4SBarry Smith 
633a40ed3dSBarry Smith   PetscFunctionBegin;
640700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
65c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
66646531bbSStefano Zampini   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
67646531bbSStefano Zampini   PetscCheckSameComm(Y,1,X,3);
68646531bbSStefano Zampini   ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr);
69646531bbSStefano Zampini   ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr);
70646531bbSStefano Zampini   ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr);
71646531bbSStefano Zampini   ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr);
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);
741ccb7edbSStefano Zampini   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (Y)");
751ccb7edbSStefano Zampini   if (!X->assembled) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix (X)");
76*edf9a46cSStefano Zampini   if (a == (PetscScalar)0.0) PetscFunctionReturn(0);
77*edf9a46cSStefano Zampini   if (Y == X) {
78*edf9a46cSStefano Zampini     ierr = MatScale(Y,1.0 + a);CHKERRQ(ierr);
79*edf9a46cSStefano Zampini     PetscFunctionReturn(0);
80*edf9a46cSStefano Zampini   }
81a683ea4aSPierre Jolivet   ierr = MatGetType(X,&t1);CHKERRQ(ierr);
82a683ea4aSPierre Jolivet   ierr = MatGetType(Y,&t2);CHKERRQ(ierr);
83a683ea4aSPierre Jolivet   ierr = PetscStrcmp(t1,t2,&sametype);CHKERRQ(ierr);
84e8136da8SHong Zhang   ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
850ff8bee4SStefano Zampini   if (Y->ops->axpy && sametype) {
86f4df32b1SMatthew Knepley     ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr);
87d4bb536fSBarry Smith   } else {
88a683ea4aSPierre Jolivet     ierr = PetscStrcmp(t1,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr);
89a683ea4aSPierre Jolivet     if (transpose) {
90d54c938cSPierre Jolivet         ierr = MatTransposeAXPY_Private(Y,a,X,str,X);CHKERRQ(ierr);
91a683ea4aSPierre Jolivet     } else {
92a683ea4aSPierre Jolivet       ierr = PetscStrcmp(t2,MATTRANSPOSEMAT,&transpose);CHKERRQ(ierr);
93a683ea4aSPierre Jolivet       if (transpose) {
94d54c938cSPierre Jolivet         ierr = MatTransposeAXPY_Private(Y,a,X,str,Y);CHKERRQ(ierr);
95a683ea4aSPierre Jolivet       } else {
96f4df32b1SMatthew Knepley         ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
97607cd303SBarry Smith       }
98a683ea4aSPierre Jolivet     }
99a683ea4aSPierre Jolivet   }
100e8136da8SHong Zhang   ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
101607cd303SBarry Smith   PetscFunctionReturn(0);
102607cd303SBarry Smith }
103607cd303SBarry Smith 
104646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
105646531bbSStefano Zampini {
106646531bbSStefano Zampini   PetscErrorCode ierr;
10727afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL;
108646531bbSStefano Zampini 
109646531bbSStefano Zampini   PetscFunctionBegin;
110646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
111646531bbSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr);
112646531bbSStefano Zampini   if (preall) {
113646531bbSStefano Zampini     ierr = (*preall)(Y,X,B);CHKERRQ(ierr);
114646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
115646531bbSStefano Zampini     Mat      preallocator;
116646531bbSStefano Zampini     PetscInt r,rstart,rend;
117646531bbSStefano Zampini     PetscInt m,n,M,N;
118646531bbSStefano Zampini 
1196eab9c35SStefano Zampini     ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr);
1206eab9c35SStefano Zampini     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
121646531bbSStefano Zampini     ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr);
122646531bbSStefano Zampini     ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr);
123646531bbSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr);
124646531bbSStefano Zampini     ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr);
125646531bbSStefano Zampini     ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr);
126646531bbSStefano Zampini     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
127646531bbSStefano Zampini     ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr);
128646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
129646531bbSStefano Zampini       PetscInt          ncols;
130646531bbSStefano Zampini       const PetscInt    *row;
131646531bbSStefano Zampini       const PetscScalar *vals;
132646531bbSStefano Zampini 
133646531bbSStefano Zampini       ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
134646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
135646531bbSStefano Zampini       ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
136646531bbSStefano Zampini       ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
137646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
138646531bbSStefano Zampini       ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
139646531bbSStefano Zampini     }
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);
147646531bbSStefano Zampini     ierr = MatSetSizes(*B,m,n,M,N);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;
157*edf9a46cSStefano Zampini   PetscBool      isshell;
158*edf9a46cSStefano Zampini 
159*edf9a46cSStefano Zampini   PetscFunctionBegin;
160*edf9a46cSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)Y,MATSHELL,&isshell);CHKERRQ(ierr);
161*edf9a46cSStefano Zampini   if (isshell) { /* MatShell has special support for AXPY */
162*edf9a46cSStefano Zampini     PetscErrorCode (*f)(Mat,PetscScalar,Mat,MatStructure);
163*edf9a46cSStefano Zampini 
164*edf9a46cSStefano Zampini     ierr = MatGetOperation(Y,MATOP_AXPY,(void (**)(void))&f);CHKERRQ(ierr);
165*edf9a46cSStefano Zampini     if (f) {
166*edf9a46cSStefano Zampini       ierr = (*f)(Y,a,X,str);CHKERRQ(ierr);
167*edf9a46cSStefano Zampini       PetscFunctionReturn(0);
168*edf9a46cSStefano Zampini     }
169*edf9a46cSStefano Zampini   }
170*edf9a46cSStefano Zampini   if (str != DIFFERENT_NONZERO_PATTERN) {
171*edf9a46cSStefano Zampini     PetscInt          i,start,end,j,ncols,m,n;
17238baddfdSBarry Smith     const PetscInt    *row;
173b3cc6726SBarry Smith     PetscScalar       *val;
174b3cc6726SBarry Smith     const PetscScalar *vals;
175607cd303SBarry Smith 
1768dadbd76SSatish Balay     ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
17790f02eecSBarry Smith     ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
1786eab9c35SStefano Zampini     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
179f4df32b1SMatthew Knepley     if (a == 1.0) {
180d4bb536fSBarry Smith       for (i = start; i < end; i++) {
181d4bb536fSBarry Smith         ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
182d4bb536fSBarry Smith         ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
183d4bb536fSBarry Smith         ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
184d4bb536fSBarry Smith       }
185d4bb536fSBarry Smith     } else {
18621a3365eSStefano Zampini       PetscInt vs = 100;
18721a3365eSStefano Zampini       /* realloc if needed, as this function may be used in parallel */
18821a3365eSStefano Zampini       ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
18906be10caSBarry Smith       for (i=start; i<end; i++) {
190b3cc6726SBarry Smith         ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
19121a3365eSStefano Zampini         if (vs < ncols) {
19221a3365eSStefano Zampini           vs   = PetscMin(2*ncols,n);
19321a3365eSStefano Zampini           ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr);
1946f79c3a4SBarry Smith         }
19521a3365eSStefano Zampini         for (j=0; j<ncols; j++) val[j] = a*vals[j];
196b3cc6726SBarry Smith         ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
197b3cc6726SBarry Smith         ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
1986f79c3a4SBarry Smith       }
199b3cc6726SBarry Smith       ierr = PetscFree(val);CHKERRQ(ierr);
200d4bb536fSBarry Smith     }
2016eab9c35SStefano Zampini     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
2026d4a8577SBarry Smith     ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2036d4a8577SBarry Smith     ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
204*edf9a46cSStefano Zampini   } else {
205*edf9a46cSStefano Zampini     Mat B;
206*edf9a46cSStefano Zampini 
207*edf9a46cSStefano Zampini     ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr);
208*edf9a46cSStefano Zampini     ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr);
209*edf9a46cSStefano Zampini     ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr);
210*edf9a46cSStefano Zampini   }
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
2126f79c3a4SBarry Smith }
213052efed2SBarry Smith 
214ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
215ec7775f6SShri Abhyankar {
216ec7775f6SShri Abhyankar   PetscInt          i,start,end,j,ncols,m,n;
217ec7775f6SShri Abhyankar   PetscErrorCode    ierr;
218ec7775f6SShri Abhyankar   const PetscInt    *row;
219ec7775f6SShri Abhyankar   PetscScalar       *val;
220ec7775f6SShri Abhyankar   const PetscScalar *vals;
221ec7775f6SShri Abhyankar 
222ec7775f6SShri Abhyankar   PetscFunctionBegin;
223ec7775f6SShri Abhyankar   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
224ec7775f6SShri Abhyankar   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
2256eab9c35SStefano Zampini   ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr);
2266eab9c35SStefano Zampini   ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
227ec7775f6SShri Abhyankar   if (a == 1.0) {
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);
234ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
235ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
236ec7775f6SShri Abhyankar     }
237ec7775f6SShri Abhyankar   } else {
23821a3365eSStefano Zampini     PetscInt vs = 100;
23921a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
24021a3365eSStefano Zampini     ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
241ec7775f6SShri Abhyankar     for (i=start; i<end; i++) {
242ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
243ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
244ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
245ec7775f6SShri Abhyankar 
246ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
24721a3365eSStefano Zampini       if (vs < ncols) {
24821a3365eSStefano Zampini         vs   = PetscMin(2*ncols,n);
24921a3365eSStefano Zampini         ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr);
250ec7775f6SShri Abhyankar       }
25121a3365eSStefano Zampini       for (j=0; j<ncols; j++) val[j] = a*vals[j];
252ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
253ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
254ec7775f6SShri Abhyankar     }
255ec7775f6SShri Abhyankar     ierr = PetscFree(val);CHKERRQ(ierr);
256ec7775f6SShri Abhyankar   }
2576eab9c35SStefano Zampini   ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr);
2586eab9c35SStefano Zampini   ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
259ec7775f6SShri Abhyankar   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260ec7775f6SShri Abhyankar   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261ec7775f6SShri Abhyankar   PetscFunctionReturn(0);
262ec7775f6SShri Abhyankar }
263ec7775f6SShri Abhyankar 
264052efed2SBarry Smith /*@
26587828ca2SBarry Smith    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
266052efed2SBarry Smith 
2673f9fe445SBarry Smith    Neighbor-wise Collective on Mat
268fee21e36SBarry Smith 
26998a79cdbSBarry Smith    Input Parameters:
27098a79cdbSBarry Smith +  Y - the matrices
27187828ca2SBarry Smith -  a - the PetscScalar
27298a79cdbSBarry Smith 
2732860a424SLois Curfman McInnes    Level: intermediate
2742860a424SLois Curfman McInnes 
27595452b02SPatrick Sanan    Notes:
27695452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2776f33a894SBarry 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
278*edf9a46cSStefano Zampini    entry). No operation is performed when a is zero.
2796f33a894SBarry Smith 
2800c0fd78eSBarry Smith    To form Y = Y + diag(V) use MatDiagonalSet()
2810c0fd78eSBarry Smith 
2820c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale()
283052efed2SBarry Smith  @*/
2847087cfbeSBarry Smith PetscErrorCode  MatShift(Mat Y,PetscScalar a)
285052efed2SBarry Smith {
2866849ba73SBarry Smith   PetscErrorCode ierr;
287052efed2SBarry Smith 
2883a40ed3dSBarry Smith   PetscFunctionBegin;
2890700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
290ce94432eSBarry Smith   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
291ce94432eSBarry Smith   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2924994cf47SJed Brown   MatCheckPreallocated(Y,1);
2931b71c1a7SBarry Smith   if (a == 0.0) PetscFunctionReturn(0);
294b50b34bdSBarry Smith 
2951c738b47SStefano Zampini   if (Y->ops->shift) {
296f4df32b1SMatthew Knepley     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
2971c738b47SStefano Zampini   } else {
2981c738b47SStefano Zampini     ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
2991c738b47SStefano Zampini   }
3007d68702bSBarry Smith 
3015528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
303052efed2SBarry Smith }
3046d84be18SBarry Smith 
3057087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
30609f38230SBarry Smith {
30709f38230SBarry Smith   PetscErrorCode ierr;
30867576b8bSBarry Smith   PetscInt       i,start,end;
30909f38230SBarry Smith   PetscScalar    *v;
31009f38230SBarry Smith 
31109f38230SBarry Smith   PetscFunctionBegin;
31209f38230SBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
31309f38230SBarry Smith   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
31409f38230SBarry Smith   for (i=start; i<end; i++) {
31509f38230SBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
31609f38230SBarry Smith   }
31709f38230SBarry Smith   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
31809f38230SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
31909f38230SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
32009f38230SBarry Smith   PetscFunctionReturn(0);
32109f38230SBarry Smith }
32209f38230SBarry Smith 
3236d84be18SBarry Smith /*@
324f56f2b3fSBarry Smith    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
325f56f2b3fSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
326f56f2b3fSBarry Smith    INSERT_VALUES.
3276d84be18SBarry Smith 
3286d84be18SBarry Smith    Input Parameters:
32998a79cdbSBarry Smith +  Y - the input matrix
330f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
331f56f2b3fSBarry Smith -  i - INSERT_VALUES or ADD_VALUES
3326d84be18SBarry Smith 
333d083f849SBarry Smith    Neighbor-wise Collective on Mat
334fee21e36SBarry Smith 
33595452b02SPatrick Sanan    Notes:
33695452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
3376f33a894SBarry 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
3386f33a894SBarry Smith    entry).
3396f33a894SBarry Smith 
3402860a424SLois Curfman McInnes    Level: intermediate
3412860a424SLois Curfman McInnes 
3420c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale()
3436d84be18SBarry Smith @*/
3447087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
3456d84be18SBarry Smith {
3466849ba73SBarry Smith   PetscErrorCode ierr;
34767576b8bSBarry Smith   PetscInt       matlocal,veclocal;
3486d84be18SBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
3500700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
3510700a824SBarry Smith   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
35267576b8bSBarry Smith   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
35367576b8bSBarry Smith   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
35467576b8bSBarry 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);
355f56f2b3fSBarry Smith   if (Y->ops->diagonalset) {
356f56f2b3fSBarry Smith     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
35794d884c6SBarry Smith   } else {
35809f38230SBarry Smith     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
3596d84be18SBarry Smith   }
3605528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
3626d84be18SBarry Smith }
363d4bb536fSBarry Smith 
364d4bb536fSBarry Smith /*@
36504aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
366d4bb536fSBarry Smith 
3673f9fe445SBarry Smith    Logically on Mat
368fee21e36SBarry Smith 
36998a79cdbSBarry Smith    Input Parameters:
37004aac2b0SHong Zhang +  a - the PetscScalar multiplier
37104aac2b0SHong Zhang .  Y - the first matrix
37204aac2b0SHong Zhang .  X - the second matrix
37304aac2b0SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
374d4bb536fSBarry Smith 
3752860a424SLois Curfman McInnes    Level: intermediate
3762860a424SLois Curfman McInnes 
3772860a424SLois Curfman McInnes .seealso: MatAXPY()
378d4bb536fSBarry Smith  @*/
3797087cfbeSBarry Smith PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
380d4bb536fSBarry Smith {
38187828ca2SBarry Smith   PetscScalar    one = 1.0;
3826849ba73SBarry Smith   PetscErrorCode ierr;
38338baddfdSBarry Smith   PetscInt       mX,mY,nX,nY;
384d4bb536fSBarry Smith 
3853a40ed3dSBarry Smith   PetscFunctionBegin;
386c5eb9154SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
3870700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
388c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
389329f5518SBarry Smith   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
390329f5518SBarry Smith   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
391e32f2f54SBarry 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);
392f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
393cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
3943a40ed3dSBarry Smith   PetscFunctionReturn(0);
395d4bb536fSBarry Smith }
396b0a32e0cSBarry Smith 
397b0a32e0cSBarry Smith /*@
3980bacdadaSStefano Zampini     MatComputeOperator - Computes the explicit matrix
399b0a32e0cSBarry Smith 
400b0a32e0cSBarry Smith     Collective on Mat
401b0a32e0cSBarry Smith 
402b0a32e0cSBarry Smith     Input Parameter:
403186a3e20SStefano Zampini +   inmat - the matrix
404186a3e20SStefano Zampini -   mattype - the matrix type for the explicit operator
405b0a32e0cSBarry Smith 
406b0a32e0cSBarry Smith     Output Parameter:
407f3b1f45cSBarry Smith .   mat - the explict  operator
408b0a32e0cSBarry Smith 
409b0a32e0cSBarry Smith     Notes:
410186a3e20SStefano Zampini     This computation is done by applying the operators to columns of the identity matrix.
411186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
412186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
413b0a32e0cSBarry Smith 
414b0a32e0cSBarry Smith     Level: advanced
415b0a32e0cSBarry Smith 
416b0a32e0cSBarry Smith @*/
4170bacdadaSStefano Zampini PetscErrorCode  MatComputeOperator(Mat inmat,MatType mattype,Mat *mat)
418b0a32e0cSBarry Smith {
419dfbe8321SBarry Smith   PetscErrorCode ierr;
420b0a32e0cSBarry Smith 
421b0a32e0cSBarry Smith   PetscFunctionBegin;
4220700a824SBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
423186a3e20SStefano Zampini   PetscValidPointer(mat,3);
424186a3e20SStefano Zampini   ierr = MatConvert_Shell(inmat,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
42524f910e3SHong Zhang   PetscFunctionReturn(0);
42624f910e3SHong Zhang }
4274325cce7SMatthew G Knepley 
4284325cce7SMatthew G Knepley /*@
4290bacdadaSStefano Zampini     MatComputeOperatorTranspose - Computes the explicit matrix representation of
430f3b1f45cSBarry Smith         a give matrix that can apply MatMultTranspose()
431f3b1f45cSBarry Smith 
432f3b1f45cSBarry Smith     Collective on Mat
433f3b1f45cSBarry Smith 
434f3b1f45cSBarry Smith     Input Parameter:
435f3b1f45cSBarry Smith .   inmat - the matrix
436f3b1f45cSBarry Smith 
437f3b1f45cSBarry Smith     Output Parameter:
438f3b1f45cSBarry Smith .   mat - the explict  operator transposed
439f3b1f45cSBarry Smith 
440f3b1f45cSBarry Smith     Notes:
441186a3e20SStefano Zampini     This computation is done by applying the transpose of the operator to columns of the identity matrix.
442186a3e20SStefano Zampini     This routine is costly in general, and is recommended for use only with relatively small systems.
443186a3e20SStefano Zampini     Currently, this routine uses a dense matrix format if mattype == NULL.
444f3b1f45cSBarry Smith 
445f3b1f45cSBarry Smith     Level: advanced
446f3b1f45cSBarry Smith 
447f3b1f45cSBarry Smith @*/
4480bacdadaSStefano Zampini PetscErrorCode  MatComputeOperatorTranspose(Mat inmat,MatType mattype,Mat *mat)
449f3b1f45cSBarry Smith {
450186a3e20SStefano Zampini   Mat            A;
451f3b1f45cSBarry Smith   PetscErrorCode ierr;
452f3b1f45cSBarry Smith 
453f3b1f45cSBarry Smith   PetscFunctionBegin;
454f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
455186a3e20SStefano Zampini   PetscValidPointer(mat,3);
456186a3e20SStefano Zampini   ierr = MatCreateTranspose(inmat,&A);CHKERRQ(ierr);
457186a3e20SStefano Zampini   ierr = MatConvert_Shell(A,mattype ? mattype : MATDENSE,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
458186a3e20SStefano Zampini   ierr = MatDestroy(&A);CHKERRQ(ierr);
459f3b1f45cSBarry Smith   PetscFunctionReturn(0);
460f3b1f45cSBarry Smith }
461f3b1f45cSBarry Smith 
462f3b1f45cSBarry Smith /*@
4634325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
4644325cce7SMatthew G Knepley 
4654325cce7SMatthew G Knepley   Input Parameters:
4664325cce7SMatthew G Knepley + A   - The matrix
4674325cce7SMatthew G Knepley - tol - The zero tolerance
4684325cce7SMatthew G Knepley 
4694325cce7SMatthew G Knepley   Output Parameters:
4704325cce7SMatthew G Knepley . A - The chopped matrix
4714325cce7SMatthew G Knepley 
4724325cce7SMatthew G Knepley   Level: intermediate
4734325cce7SMatthew G Knepley 
4744325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries()
4753fc99919SSatish Balay  @*/
4764325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol)
4774325cce7SMatthew G Knepley {
4784325cce7SMatthew G Knepley   PetscScalar    *newVals;
4794325cce7SMatthew G Knepley   PetscInt       *newCols;
4804325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
4814325cce7SMatthew G Knepley   PetscErrorCode ierr;
4824325cce7SMatthew G Knepley 
4834325cce7SMatthew G Knepley   PetscFunctionBegin;
4844325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4854325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4864325cce7SMatthew G Knepley     PetscInt ncols;
4874325cce7SMatthew G Knepley 
4880298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4894325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4900298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4914325cce7SMatthew G Knepley   }
4924325cce7SMatthew G Knepley   numRows = rEnd - rStart;
493b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
494dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4954325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4964325cce7SMatthew G Knepley     const PetscScalar *vals;
4974325cce7SMatthew G Knepley     const PetscInt    *cols;
498fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
4994325cce7SMatthew G Knepley 
5004325cce7SMatthew G Knepley     if (r < rEnd) {
5014325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
5024325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
5034325cce7SMatthew G Knepley         newCols[c] = cols[c];
5044325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
5054325cce7SMatthew G Knepley       }
506fad45679SMatthew G. Knepley       newcols = ncols;
5074325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
508fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
5094325cce7SMatthew G Knepley     }
5104325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5114325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5124325cce7SMatthew G Knepley   }
5134325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
5144325cce7SMatthew G Knepley   PetscFunctionReturn(0);
5154325cce7SMatthew G Knepley }
516