xref: /petsc/src/mat/utils/axpy.c (revision 27afe9e8a864f3dad535868ca085e9fa806dbfc3)
16f79c3a4SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/
36f79c3a4SBarry Smith 
406be10caSBarry Smith /*@
521c89e3eSBarry Smith    MatAXPY - Computes Y = a*X + Y.
66f79c3a4SBarry Smith 
73f9fe445SBarry Smith    Logically  Collective on Mat
8fee21e36SBarry Smith 
998a79cdbSBarry Smith    Input Parameters:
10607cd303SBarry Smith +  a - the scalar multiplier
11607cd303SBarry Smith .  X - the first matrix
12607cd303SBarry Smith .  Y - the second matrix
13407f6b05SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN
14407f6b05SHong Zhang          or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
1598a79cdbSBarry Smith 
162860a424SLois Curfman McInnes    Level: intermediate
172860a424SLois Curfman McInnes 
189cf4f1e8SLois Curfman McInnes .keywords: matrix, add
19d4bb536fSBarry Smith 
202860a424SLois Curfman McInnes .seealso: MatAYPX()
2106be10caSBarry Smith  @*/
227087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
236f79c3a4SBarry Smith {
246849ba73SBarry Smith   PetscErrorCode ierr;
25646531bbSStefano Zampini   PetscInt       M1,M2,N1,N2;
26c1ac3661SBarry Smith   PetscInt       m1,m2,n1,n2;
270ff8bee4SStefano Zampini   PetscBool      sametype;
286f79c3a4SBarry Smith 
293a40ed3dSBarry Smith   PetscFunctionBegin;
300700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
31c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
32646531bbSStefano Zampini   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
33646531bbSStefano Zampini   PetscCheckSameComm(Y,1,X,3);
34646531bbSStefano Zampini   ierr = MatGetSize(X,&M1,&N1);CHKERRQ(ierr);
35646531bbSStefano Zampini   ierr = MatGetSize(Y,&M2,&N2);CHKERRQ(ierr);
36646531bbSStefano Zampini   ierr = MatGetLocalSize(X,&m1,&n1);CHKERRQ(ierr);
37646531bbSStefano Zampini   ierr = MatGetLocalSize(Y,&m2,&n2);CHKERRQ(ierr);
38646531bbSStefano 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);
39646531bbSStefano 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);
401987afe7SBarry Smith 
410ff8bee4SStefano Zampini   ierr = PetscStrcmp(((PetscObject)X)->type_name,((PetscObject)Y)->type_name,&sametype);CHKERRQ(ierr);
42e8136da8SHong Zhang   ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
430ff8bee4SStefano Zampini   if (Y->ops->axpy && sametype) {
44f4df32b1SMatthew Knepley     ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr);
45d4bb536fSBarry Smith   } else {
46646531bbSStefano Zampini     if (str != DIFFERENT_NONZERO_PATTERN) {
47f4df32b1SMatthew Knepley       ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
48646531bbSStefano Zampini     } else {
49646531bbSStefano Zampini       Mat B;
50646531bbSStefano Zampini 
51646531bbSStefano Zampini       ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr);
52646531bbSStefano Zampini       ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr);
53646531bbSStefano Zampini       ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr);
54646531bbSStefano Zampini     }
55607cd303SBarry Smith   }
56e8136da8SHong Zhang   ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
57fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
58b8ced49eSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
59b8ced49eSKarl Rupp     Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
60c41cb2e2SAlejandro Lamas Daviña   }
61d67ff14aSKarl Rupp #endif
62607cd303SBarry Smith   PetscFunctionReturn(0);
63607cd303SBarry Smith }
64607cd303SBarry Smith 
65646531bbSStefano Zampini PetscErrorCode MatAXPY_Basic_Preallocate(Mat Y, Mat X, Mat *B)
66646531bbSStefano Zampini {
67646531bbSStefano Zampini   PetscErrorCode ierr;
68*27afe9e8SStefano Zampini   PetscErrorCode (*preall)(Mat,Mat,Mat*) = NULL;
69646531bbSStefano Zampini 
70646531bbSStefano Zampini   PetscFunctionBegin;
71646531bbSStefano Zampini   /* look for any available faster alternative to the general preallocator */
72646531bbSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)Y,"MatAXPYGetPreallocation_C",&preall);CHKERRQ(ierr);
73646531bbSStefano Zampini   if (preall) {
74646531bbSStefano Zampini     ierr = (*preall)(Y,X,B);CHKERRQ(ierr);
75646531bbSStefano Zampini   } else { /* Use MatPrellocator, assumes same row-col distribution */
76646531bbSStefano Zampini     Mat      preallocator;
77646531bbSStefano Zampini     PetscInt r,rstart,rend;
78646531bbSStefano Zampini     PetscInt m,n,M,N;
79646531bbSStefano Zampini 
80646531bbSStefano Zampini     ierr = MatGetSize(Y,&M,&N);CHKERRQ(ierr);
81646531bbSStefano Zampini     ierr = MatGetLocalSize(Y,&m,&n);CHKERRQ(ierr);
82646531bbSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)Y),&preallocator);CHKERRQ(ierr);
83646531bbSStefano Zampini     ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr);
84646531bbSStefano Zampini     ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr);
85646531bbSStefano Zampini     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
86646531bbSStefano Zampini     ierr = MatGetOwnershipRange(preallocator,&rstart,&rend);CHKERRQ(ierr);
87646531bbSStefano Zampini     for (r = rstart; r < rend; ++r) {
88646531bbSStefano Zampini       PetscInt          ncols;
89646531bbSStefano Zampini       const PetscInt    *row;
90646531bbSStefano Zampini       const PetscScalar *vals;
91646531bbSStefano Zampini 
92646531bbSStefano Zampini       ierr = MatGetRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
93646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
94646531bbSStefano Zampini       ierr = MatRestoreRow(Y,r,&ncols,&row,&vals);CHKERRQ(ierr);
95646531bbSStefano Zampini       ierr = MatGetRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
96646531bbSStefano Zampini       ierr = MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);CHKERRQ(ierr);
97646531bbSStefano Zampini       ierr = MatRestoreRow(X,r,&ncols,&row,&vals);CHKERRQ(ierr);
98646531bbSStefano Zampini     }
99646531bbSStefano Zampini     ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100646531bbSStefano Zampini     ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101646531bbSStefano Zampini 
102646531bbSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)Y),B);CHKERRQ(ierr);
103646531bbSStefano Zampini     ierr = MatSetType(*B,((PetscObject)Y)->type_name);CHKERRQ(ierr);
104646531bbSStefano Zampini     ierr = MatSetSizes(*B,m,n,M,N);CHKERRQ(ierr);
105646531bbSStefano Zampini     ierr = MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*B);CHKERRQ(ierr);
106646531bbSStefano Zampini     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
107646531bbSStefano Zampini   }
108646531bbSStefano Zampini   PetscFunctionReturn(0);
109646531bbSStefano Zampini }
110646531bbSStefano Zampini 
111f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
112607cd303SBarry Smith {
11338baddfdSBarry Smith   PetscInt          i,start,end,j,ncols,m,n;
1146849ba73SBarry Smith   PetscErrorCode    ierr;
11538baddfdSBarry Smith   const PetscInt    *row;
116b3cc6726SBarry Smith   PetscScalar       *val;
117b3cc6726SBarry Smith   const PetscScalar *vals;
118607cd303SBarry Smith 
119607cd303SBarry Smith   PetscFunctionBegin;
1208dadbd76SSatish Balay   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
12190f02eecSBarry Smith   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
122f4df32b1SMatthew Knepley   if (a == 1.0) {
123d4bb536fSBarry Smith     for (i = start; i < end; i++) {
124d4bb536fSBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
125d4bb536fSBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
126d4bb536fSBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
127d4bb536fSBarry Smith     }
128d4bb536fSBarry Smith   } else {
12921a3365eSStefano Zampini     PetscInt vs = 100;
13021a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
13121a3365eSStefano Zampini     ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
13206be10caSBarry Smith     for (i=start; i<end; i++) {
133b3cc6726SBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
13421a3365eSStefano Zampini       if (vs < ncols) {
13521a3365eSStefano Zampini         vs   = PetscMin(2*ncols,n);
13621a3365eSStefano Zampini         ierr = PetscRealloc(vs*sizeof(*val),&val);CHKERRQ(ierr);
1376f79c3a4SBarry Smith       }
13821a3365eSStefano Zampini       for (j=0; j<ncols; j++) val[j] = a*vals[j];
139b3cc6726SBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
140b3cc6726SBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
1416f79c3a4SBarry Smith     }
142b3cc6726SBarry Smith     ierr = PetscFree(val);CHKERRQ(ierr);
143d4bb536fSBarry Smith   }
1446d4a8577SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1456d4a8577SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
1476f79c3a4SBarry Smith }
148052efed2SBarry Smith 
149ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
150ec7775f6SShri Abhyankar {
151ec7775f6SShri Abhyankar   PetscInt          i,start,end,j,ncols,m,n;
152ec7775f6SShri Abhyankar   PetscErrorCode    ierr;
153ec7775f6SShri Abhyankar   const PetscInt    *row;
154ec7775f6SShri Abhyankar   PetscScalar       *val;
155ec7775f6SShri Abhyankar   const PetscScalar *vals;
156ec7775f6SShri Abhyankar 
157ec7775f6SShri Abhyankar   PetscFunctionBegin;
158ec7775f6SShri Abhyankar   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
159ec7775f6SShri Abhyankar   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
160ec7775f6SShri Abhyankar   if (a == 1.0) {
161ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
162ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
163ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
164ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
165ec7775f6SShri Abhyankar 
166ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
167ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
168ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
169ec7775f6SShri Abhyankar     }
170ec7775f6SShri Abhyankar   } else {
17121a3365eSStefano Zampini     PetscInt vs = 100;
17221a3365eSStefano Zampini     /* realloc if needed, as this function may be used in parallel */
17321a3365eSStefano Zampini     ierr = PetscMalloc1(vs,&val);CHKERRQ(ierr);
174ec7775f6SShri Abhyankar     for (i=start; i<end; i++) {
175ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
176ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
177ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
178ec7775f6SShri Abhyankar 
179ec7775f6SShri Abhyankar       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);
183ec7775f6SShri Abhyankar       }
18421a3365eSStefano Zampini       for (j=0; j<ncols; j++) val[j] = a*vals[j];
185ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
186ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
187ec7775f6SShri Abhyankar     }
188ec7775f6SShri Abhyankar     ierr = PetscFree(val);CHKERRQ(ierr);
189ec7775f6SShri Abhyankar   }
190ec7775f6SShri Abhyankar   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191ec7775f6SShri Abhyankar   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192ec7775f6SShri Abhyankar   PetscFunctionReturn(0);
193ec7775f6SShri Abhyankar }
194ec7775f6SShri Abhyankar 
195052efed2SBarry Smith /*@
19687828ca2SBarry Smith    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
197052efed2SBarry Smith 
1983f9fe445SBarry Smith    Neighbor-wise Collective on Mat
199fee21e36SBarry Smith 
20098a79cdbSBarry Smith    Input Parameters:
20198a79cdbSBarry Smith +  Y - the matrices
20287828ca2SBarry Smith -  a - the PetscScalar
20398a79cdbSBarry Smith 
2042860a424SLois Curfman McInnes    Level: intermediate
2052860a424SLois Curfman McInnes 
20695452b02SPatrick Sanan    Notes:
20795452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2086f33a894SBarry 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
2096f33a894SBarry Smith    entry).
2106f33a894SBarry Smith 
2110c0fd78eSBarry Smith    To form Y = Y + diag(V) use MatDiagonalSet()
2120c0fd78eSBarry Smith 
2136f33a894SBarry Smith    Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is
2146f33a894SBarry Smith     preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix.
2156f33a894SBarry Smith 
216052efed2SBarry Smith .keywords: matrix, add, shift
2176b9ee512SLois Curfman McInnes 
2180c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale()
219052efed2SBarry Smith  @*/
2207087cfbeSBarry Smith PetscErrorCode  MatShift(Mat Y,PetscScalar a)
221052efed2SBarry Smith {
2226849ba73SBarry Smith   PetscErrorCode ierr;
223052efed2SBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
2250700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
226ce94432eSBarry Smith   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
227ce94432eSBarry Smith   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2284994cf47SJed Brown   MatCheckPreallocated(Y,1);
229b50b34bdSBarry Smith 
2301c738b47SStefano Zampini   if (Y->ops->shift) {
231f4df32b1SMatthew Knepley     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
2321c738b47SStefano Zampini   } else {
2331c738b47SStefano Zampini     ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
2341c738b47SStefano Zampini   }
2357d68702bSBarry Smith 
2365528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
237fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
238b8ced49eSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
239b8ced49eSKarl Rupp     Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
240c41cb2e2SAlejandro Lamas Daviña   }
241d67ff14aSKarl Rupp #endif
2423a40ed3dSBarry Smith   PetscFunctionReturn(0);
243052efed2SBarry Smith }
2446d84be18SBarry Smith 
2457087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
24609f38230SBarry Smith {
24709f38230SBarry Smith   PetscErrorCode ierr;
24867576b8bSBarry Smith   PetscInt       i,start,end;
24909f38230SBarry Smith   PetscScalar    *v;
25009f38230SBarry Smith 
25109f38230SBarry Smith   PetscFunctionBegin;
25209f38230SBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
25309f38230SBarry Smith   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
25409f38230SBarry Smith   for (i=start; i<end; i++) {
25509f38230SBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
25609f38230SBarry Smith   }
25709f38230SBarry Smith   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
25809f38230SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25909f38230SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
26009f38230SBarry Smith   PetscFunctionReturn(0);
26109f38230SBarry Smith }
26209f38230SBarry Smith 
2636d84be18SBarry Smith /*@
264f56f2b3fSBarry Smith    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
265f56f2b3fSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
266f56f2b3fSBarry Smith    INSERT_VALUES.
2676d84be18SBarry Smith 
2686d84be18SBarry Smith    Input Parameters:
26998a79cdbSBarry Smith +  Y - the input matrix
270f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
271f56f2b3fSBarry Smith -  i - INSERT_VALUES or ADD_VALUES
2726d84be18SBarry Smith 
2733f9fe445SBarry Smith    Neighbor-wise Collective on Mat and Vec
274fee21e36SBarry Smith 
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
2786f33a894SBarry Smith    entry).
2796f33a894SBarry Smith 
2802860a424SLois Curfman McInnes    Level: intermediate
2812860a424SLois Curfman McInnes 
2826b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal
2836b9ee512SLois Curfman McInnes 
2840c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale()
2856d84be18SBarry Smith @*/
2867087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
2876d84be18SBarry Smith {
2886849ba73SBarry Smith   PetscErrorCode ierr;
28967576b8bSBarry Smith   PetscInt       matlocal,veclocal;
2906d84be18SBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
2920700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
2930700a824SBarry Smith   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
29467576b8bSBarry Smith   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
29567576b8bSBarry Smith   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
29667576b8bSBarry 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);
297f56f2b3fSBarry Smith   if (Y->ops->diagonalset) {
298f56f2b3fSBarry Smith     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
29994d884c6SBarry Smith   } else {
30009f38230SBarry Smith     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
3016d84be18SBarry Smith   }
3025528ad4fSTristan Konolige   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
3033a40ed3dSBarry Smith   PetscFunctionReturn(0);
3046d84be18SBarry Smith }
305d4bb536fSBarry Smith 
306d4bb536fSBarry Smith /*@
30704aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
308d4bb536fSBarry Smith 
3093f9fe445SBarry Smith    Logically on Mat
310fee21e36SBarry Smith 
31198a79cdbSBarry Smith    Input Parameters:
31204aac2b0SHong Zhang +  a - the PetscScalar multiplier
31304aac2b0SHong Zhang .  Y - the first matrix
31404aac2b0SHong Zhang .  X - the second matrix
31504aac2b0SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
316d4bb536fSBarry Smith 
3172860a424SLois Curfman McInnes    Level: intermediate
3182860a424SLois Curfman McInnes 
319d4bb536fSBarry Smith .keywords: matrix, add
320d4bb536fSBarry Smith 
3212860a424SLois Curfman McInnes .seealso: MatAXPY()
322d4bb536fSBarry Smith  @*/
3237087cfbeSBarry Smith PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
324d4bb536fSBarry Smith {
32587828ca2SBarry Smith   PetscScalar    one = 1.0;
3266849ba73SBarry Smith   PetscErrorCode ierr;
32738baddfdSBarry Smith   PetscInt       mX,mY,nX,nY;
328d4bb536fSBarry Smith 
3293a40ed3dSBarry Smith   PetscFunctionBegin;
330c5eb9154SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
3310700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
332c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
333329f5518SBarry Smith   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
334329f5518SBarry Smith   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
335e32f2f54SBarry 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);
336d4bb536fSBarry Smith 
337f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
338cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
340d4bb536fSBarry Smith }
341b0a32e0cSBarry Smith 
342b0a32e0cSBarry Smith /*@
343b0a32e0cSBarry Smith     MatComputeExplicitOperator - Computes the explicit matrix
344b0a32e0cSBarry Smith 
345b0a32e0cSBarry Smith     Collective on Mat
346b0a32e0cSBarry Smith 
347b0a32e0cSBarry Smith     Input Parameter:
348b0a32e0cSBarry Smith .   inmat - the matrix
349b0a32e0cSBarry Smith 
350b0a32e0cSBarry Smith     Output Parameter:
351f3b1f45cSBarry Smith .   mat - the explict  operator
352b0a32e0cSBarry Smith 
353b0a32e0cSBarry Smith     Notes:
354b0a32e0cSBarry Smith     This computation is done by applying the operators to columns of the
355b0a32e0cSBarry Smith     identity matrix.
356b0a32e0cSBarry Smith 
357b0a32e0cSBarry Smith     Currently, this routine uses a dense matrix format when 1 processor
358b0a32e0cSBarry Smith     is used and a sparse format otherwise.  This routine is costly in general,
359b0a32e0cSBarry Smith     and is recommended for use only with relatively small systems.
360b0a32e0cSBarry Smith 
361b0a32e0cSBarry Smith     Level: advanced
362b0a32e0cSBarry Smith 
363b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator
364b0a32e0cSBarry Smith @*/
3657087cfbeSBarry Smith PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
366b0a32e0cSBarry Smith {
367dfbe8321SBarry Smith   PetscErrorCode ierr;
368b0a32e0cSBarry Smith   MPI_Comm       comm;
36938baddfdSBarry Smith   PetscMPIInt    size;
370b0a32e0cSBarry Smith 
371b0a32e0cSBarry Smith   PetscFunctionBegin;
3720700a824SBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
3734482741eSBarry Smith   PetscValidPointer(mat,2);
374b0a32e0cSBarry Smith 
375ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
376b0a32e0cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
377b3d09e86SJed Brown   ierr = MatConvert_Shell(inmat,size == 1 ? MATSEQDENSE : MATAIJ,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
37824f910e3SHong Zhang   PetscFunctionReturn(0);
37924f910e3SHong Zhang }
3804325cce7SMatthew G Knepley 
3814325cce7SMatthew G Knepley /*@
382f3b1f45cSBarry Smith     MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of
383f3b1f45cSBarry Smith         a give matrix that can apply MatMultTranspose()
384f3b1f45cSBarry Smith 
385f3b1f45cSBarry Smith     Collective on Mat
386f3b1f45cSBarry Smith 
387f3b1f45cSBarry Smith     Input Parameter:
388f3b1f45cSBarry Smith .   inmat - the matrix
389f3b1f45cSBarry Smith 
390f3b1f45cSBarry Smith     Output Parameter:
391f3b1f45cSBarry Smith .   mat - the explict  operator transposed
392f3b1f45cSBarry Smith 
393f3b1f45cSBarry Smith     Notes:
394f3b1f45cSBarry Smith     This computation is done by applying the transpose of the operator to columns of the
395f3b1f45cSBarry Smith     identity matrix.
396f3b1f45cSBarry Smith 
397f3b1f45cSBarry Smith     Currently, this routine uses a dense matrix format when 1 processor
398f3b1f45cSBarry Smith     is used and a sparse format otherwise.  This routine is costly in general,
399f3b1f45cSBarry Smith     and is recommended for use only with relatively small systems.
400f3b1f45cSBarry Smith 
401f3b1f45cSBarry Smith     Level: advanced
402f3b1f45cSBarry Smith 
403f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator
404f3b1f45cSBarry Smith @*/
405f3b1f45cSBarry Smith PetscErrorCode  MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat)
406f3b1f45cSBarry Smith {
407f3b1f45cSBarry Smith   Vec            in,out;
408f3b1f45cSBarry Smith   PetscErrorCode ierr;
409f3b1f45cSBarry Smith   PetscInt       i,m,n,M,N,*rows,start,end;
410f3b1f45cSBarry Smith   MPI_Comm       comm;
411f3b1f45cSBarry Smith   PetscScalar    *array,zero = 0.0,one = 1.0;
412f3b1f45cSBarry Smith   PetscMPIInt    size;
413f3b1f45cSBarry Smith 
414f3b1f45cSBarry Smith   PetscFunctionBegin;
415f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
416f3b1f45cSBarry Smith   PetscValidPointer(mat,2);
417f3b1f45cSBarry Smith 
418f3b1f45cSBarry Smith   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
419f3b1f45cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
420f3b1f45cSBarry Smith 
421f3b1f45cSBarry Smith   ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr);
422f3b1f45cSBarry Smith   ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr);
423f3b1f45cSBarry Smith   ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr);
424f3b1f45cSBarry Smith   ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
425f3b1f45cSBarry Smith   ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr);
426f3b1f45cSBarry Smith   ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr);
427f3b1f45cSBarry Smith   for (i=0; i<m; i++) rows[i] = start + i;
428f3b1f45cSBarry Smith 
429f3b1f45cSBarry Smith   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
430f3b1f45cSBarry Smith   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
431f3b1f45cSBarry Smith   if (size == 1) {
432f3b1f45cSBarry Smith     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
433f3b1f45cSBarry Smith     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
434f3b1f45cSBarry Smith   } else {
435f3b1f45cSBarry Smith     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
436f3b1f45cSBarry Smith     ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr);
437f3b1f45cSBarry Smith   }
438f3b1f45cSBarry Smith 
439f3b1f45cSBarry Smith   for (i=0; i<N; i++) {
440f3b1f45cSBarry Smith 
441f3b1f45cSBarry Smith     ierr = VecSet(in,zero);CHKERRQ(ierr);
442f3b1f45cSBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
443f3b1f45cSBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
444f3b1f45cSBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
445f3b1f45cSBarry Smith 
446f3b1f45cSBarry Smith     ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr);
447f3b1f45cSBarry Smith 
448f3b1f45cSBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
449f3b1f45cSBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
450f3b1f45cSBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
451f3b1f45cSBarry Smith 
452f3b1f45cSBarry Smith   }
453f3b1f45cSBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
454f3b1f45cSBarry Smith   ierr = VecDestroy(&out);CHKERRQ(ierr);
455f3b1f45cSBarry Smith   ierr = VecDestroy(&in);CHKERRQ(ierr);
456f3b1f45cSBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
457f3b1f45cSBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
458f3b1f45cSBarry Smith   PetscFunctionReturn(0);
459f3b1f45cSBarry Smith }
460f3b1f45cSBarry Smith 
461f3b1f45cSBarry Smith /*@
4624325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
4634325cce7SMatthew G Knepley 
4644325cce7SMatthew G Knepley   Input Parameters:
4654325cce7SMatthew G Knepley + A   - The matrix
4664325cce7SMatthew G Knepley - tol - The zero tolerance
4674325cce7SMatthew G Knepley 
4684325cce7SMatthew G Knepley   Output Parameters:
4694325cce7SMatthew G Knepley . A - The chopped matrix
4704325cce7SMatthew G Knepley 
4714325cce7SMatthew G Knepley   Level: intermediate
4724325cce7SMatthew G Knepley 
4734325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries()
4743fc99919SSatish Balay  @*/
4754325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol)
4764325cce7SMatthew G Knepley {
4774325cce7SMatthew G Knepley   PetscScalar    *newVals;
4784325cce7SMatthew G Knepley   PetscInt       *newCols;
4794325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
4804325cce7SMatthew G Knepley   PetscErrorCode ierr;
4814325cce7SMatthew G Knepley 
4824325cce7SMatthew G Knepley   PetscFunctionBegin;
4834325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4844325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4854325cce7SMatthew G Knepley     PetscInt ncols;
4864325cce7SMatthew G Knepley 
4870298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4884325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4890298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4904325cce7SMatthew G Knepley   }
4914325cce7SMatthew G Knepley   numRows = rEnd - rStart;
492b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
493dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4944325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4954325cce7SMatthew G Knepley     const PetscScalar *vals;
4964325cce7SMatthew G Knepley     const PetscInt    *cols;
497fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
4984325cce7SMatthew G Knepley 
4994325cce7SMatthew G Knepley     if (r < rEnd) {
5004325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
5014325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
5024325cce7SMatthew G Knepley         newCols[c] = cols[c];
5034325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
5044325cce7SMatthew G Knepley       }
505fad45679SMatthew G. Knepley       newcols = ncols;
5064325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
507fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
5084325cce7SMatthew G Knepley     }
5094325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5104325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5114325cce7SMatthew G Knepley   }
5124325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
5134325cce7SMatthew G Knepley   PetscFunctionReturn(0);
5144325cce7SMatthew G Knepley }
515