xref: /petsc/src/mat/utils/axpy.c (revision fd3149344b786979040ec64d65bd274fdf13dcc2)
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;
25c1ac3661SBarry Smith   PetscInt       m1,m2,n1,n2;
260ff8bee4SStefano Zampini   PetscBool      sametype;
276f79c3a4SBarry Smith 
283a40ed3dSBarry Smith   PetscFunctionBegin;
290700a824SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
300700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
31c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
32273d9f13SBarry Smith   ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr);
33273d9f13SBarry Smith   ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr);
34e32f2f54SBarry Smith   if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %D %D %D %D",m1,m2,n1,n2);
351987afe7SBarry Smith 
360ff8bee4SStefano Zampini   ierr = PetscStrcmp(((PetscObject)X)->type_name,((PetscObject)Y)->type_name,&sametype);CHKERRQ(ierr);
37e8136da8SHong Zhang   ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
380ff8bee4SStefano Zampini   if (Y->ops->axpy && sametype) {
39f4df32b1SMatthew Knepley     ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr);
40d4bb536fSBarry Smith   } else {
41f4df32b1SMatthew Knepley     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
42607cd303SBarry Smith   }
43e8136da8SHong Zhang   ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
44*fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
45b8ced49eSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
46b8ced49eSKarl Rupp     Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
47c41cb2e2SAlejandro Lamas Daviña   }
48d67ff14aSKarl Rupp #endif
49607cd303SBarry Smith   PetscFunctionReturn(0);
50607cd303SBarry Smith }
51607cd303SBarry Smith 
52f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
53607cd303SBarry Smith {
5438baddfdSBarry Smith   PetscInt          i,start,end,j,ncols,m,n;
556849ba73SBarry Smith   PetscErrorCode    ierr;
5638baddfdSBarry Smith   const PetscInt    *row;
57b3cc6726SBarry Smith   PetscScalar       *val;
58b3cc6726SBarry Smith   const PetscScalar *vals;
59607cd303SBarry Smith 
60607cd303SBarry Smith   PetscFunctionBegin;
618dadbd76SSatish Balay   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
6290f02eecSBarry Smith   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
63f4df32b1SMatthew Knepley   if (a == 1.0) {
64d4bb536fSBarry Smith     for (i = start; i < end; i++) {
65d4bb536fSBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
66d4bb536fSBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
67d4bb536fSBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
68d4bb536fSBarry Smith     }
69d4bb536fSBarry Smith   } else {
70854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr);
7106be10caSBarry Smith     for (i=start; i<end; i++) {
72b3cc6726SBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
7306be10caSBarry Smith       for (j=0; j<ncols; j++) {
74f4df32b1SMatthew Knepley         val[j] = a*vals[j];
756f79c3a4SBarry Smith       }
76b3cc6726SBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
77b3cc6726SBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
786f79c3a4SBarry Smith     }
79b3cc6726SBarry Smith     ierr = PetscFree(val);CHKERRQ(ierr);
80d4bb536fSBarry Smith   }
816d4a8577SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
826d4a8577SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
833a40ed3dSBarry Smith   PetscFunctionReturn(0);
846f79c3a4SBarry Smith }
85052efed2SBarry Smith 
86ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
87ec7775f6SShri Abhyankar {
88ec7775f6SShri Abhyankar   PetscInt          i,start,end,j,ncols,m,n;
89ec7775f6SShri Abhyankar   PetscErrorCode    ierr;
90ec7775f6SShri Abhyankar   const PetscInt    *row;
91ec7775f6SShri Abhyankar   PetscScalar       *val;
92ec7775f6SShri Abhyankar   const PetscScalar *vals;
93ec7775f6SShri Abhyankar 
94ec7775f6SShri Abhyankar   PetscFunctionBegin;
95ec7775f6SShri Abhyankar   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
96ec7775f6SShri Abhyankar   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
97ec7775f6SShri Abhyankar   if (a == 1.0) {
98ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
99ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
100ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
101ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
102ec7775f6SShri Abhyankar 
103ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
104ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
105ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
106ec7775f6SShri Abhyankar     }
107ec7775f6SShri Abhyankar   } else {
108854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr);
109ec7775f6SShri Abhyankar     for (i=start; i<end; i++) {
110ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
111ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
112ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
113ec7775f6SShri Abhyankar 
114ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
115ec7775f6SShri Abhyankar       for (j=0; j<ncols; j++) {
116ec7775f6SShri Abhyankar         val[j] = a*vals[j];
117ec7775f6SShri Abhyankar       }
118ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
119ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
120ec7775f6SShri Abhyankar     }
121ec7775f6SShri Abhyankar     ierr = PetscFree(val);CHKERRQ(ierr);
122ec7775f6SShri Abhyankar   }
123ec7775f6SShri Abhyankar   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
124ec7775f6SShri Abhyankar   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
125ec7775f6SShri Abhyankar   PetscFunctionReturn(0);
126ec7775f6SShri Abhyankar }
127ec7775f6SShri Abhyankar 
128052efed2SBarry Smith /*@
12987828ca2SBarry Smith    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
130052efed2SBarry Smith 
1313f9fe445SBarry Smith    Neighbor-wise Collective on Mat
132fee21e36SBarry Smith 
13398a79cdbSBarry Smith    Input Parameters:
13498a79cdbSBarry Smith +  Y - the matrices
13587828ca2SBarry Smith -  a - the PetscScalar
13698a79cdbSBarry Smith 
1372860a424SLois Curfman McInnes    Level: intermediate
1382860a424SLois Curfman McInnes 
13995452b02SPatrick Sanan    Notes:
14095452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
1416f33a894SBarry 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
1426f33a894SBarry Smith    entry).
1436f33a894SBarry Smith 
1440c0fd78eSBarry Smith    To form Y = Y + diag(V) use MatDiagonalSet()
1450c0fd78eSBarry Smith 
1466f33a894SBarry Smith    Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is
1476f33a894SBarry Smith     preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix.
1486f33a894SBarry Smith 
149052efed2SBarry Smith .keywords: matrix, add, shift
1506b9ee512SLois Curfman McInnes 
1510c0fd78eSBarry Smith .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale()
152052efed2SBarry Smith  @*/
1537087cfbeSBarry Smith PetscErrorCode  MatShift(Mat Y,PetscScalar a)
154052efed2SBarry Smith {
1556849ba73SBarry Smith   PetscErrorCode ierr;
156052efed2SBarry Smith 
1573a40ed3dSBarry Smith   PetscFunctionBegin;
1580700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
159ce94432eSBarry Smith   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
160ce94432eSBarry Smith   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1614994cf47SJed Brown   MatCheckPreallocated(Y,1);
162b50b34bdSBarry Smith 
1631c738b47SStefano Zampini   if (Y->ops->shift) {
164f4df32b1SMatthew Knepley     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
1651c738b47SStefano Zampini   } else {
1661c738b47SStefano Zampini     ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
1671c738b47SStefano Zampini   }
1687d68702bSBarry Smith 
169*fd314934SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
170b8ced49eSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
171b8ced49eSKarl Rupp     Y->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
172c41cb2e2SAlejandro Lamas Daviña   }
173d67ff14aSKarl Rupp #endif
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175052efed2SBarry Smith }
1766d84be18SBarry Smith 
1777087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
17809f38230SBarry Smith {
17909f38230SBarry Smith   PetscErrorCode ierr;
18067576b8bSBarry Smith   PetscInt       i,start,end;
18109f38230SBarry Smith   PetscScalar    *v;
18209f38230SBarry Smith 
18309f38230SBarry Smith   PetscFunctionBegin;
18409f38230SBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
18509f38230SBarry Smith   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
18609f38230SBarry Smith   for (i=start; i<end; i++) {
18709f38230SBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
18809f38230SBarry Smith   }
18909f38230SBarry Smith   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
19009f38230SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19109f38230SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19209f38230SBarry Smith   PetscFunctionReturn(0);
19309f38230SBarry Smith }
19409f38230SBarry Smith 
1956d84be18SBarry Smith /*@
196f56f2b3fSBarry Smith    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
197f56f2b3fSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
198f56f2b3fSBarry Smith    INSERT_VALUES.
1996d84be18SBarry Smith 
2006d84be18SBarry Smith    Input Parameters:
20198a79cdbSBarry Smith +  Y - the input matrix
202f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
203f56f2b3fSBarry Smith -  i - INSERT_VALUES or ADD_VALUES
2046d84be18SBarry Smith 
2053f9fe445SBarry Smith    Neighbor-wise Collective on Mat and Vec
206fee21e36SBarry Smith 
20795452b02SPatrick Sanan    Notes:
20895452b02SPatrick Sanan     If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2096f33a894SBarry 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
2106f33a894SBarry Smith    entry).
2116f33a894SBarry Smith 
2122860a424SLois Curfman McInnes    Level: intermediate
2132860a424SLois Curfman McInnes 
2146b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal
2156b9ee512SLois Curfman McInnes 
2160c0fd78eSBarry Smith .seealso: MatShift(), MatScale(), MatDiagonalScale()
2176d84be18SBarry Smith @*/
2187087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
2196d84be18SBarry Smith {
2206849ba73SBarry Smith   PetscErrorCode ierr;
22167576b8bSBarry Smith   PetscInt       matlocal,veclocal;
2226d84be18SBarry Smith 
2233a40ed3dSBarry Smith   PetscFunctionBegin;
2240700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
2250700a824SBarry Smith   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
22667576b8bSBarry Smith   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
22767576b8bSBarry Smith   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
22867576b8bSBarry 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);
229f56f2b3fSBarry Smith   if (Y->ops->diagonalset) {
230f56f2b3fSBarry Smith     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
23194d884c6SBarry Smith   } else {
23209f38230SBarry Smith     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
2336d84be18SBarry Smith   }
2343a40ed3dSBarry Smith   PetscFunctionReturn(0);
2356d84be18SBarry Smith }
236d4bb536fSBarry Smith 
237d4bb536fSBarry Smith /*@
23804aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
239d4bb536fSBarry Smith 
2403f9fe445SBarry Smith    Logically on Mat
241fee21e36SBarry Smith 
24298a79cdbSBarry Smith    Input Parameters:
24304aac2b0SHong Zhang +  a - the PetscScalar multiplier
24404aac2b0SHong Zhang .  Y - the first matrix
24504aac2b0SHong Zhang .  X - the second matrix
24604aac2b0SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
247d4bb536fSBarry Smith 
2482860a424SLois Curfman McInnes    Level: intermediate
2492860a424SLois Curfman McInnes 
250d4bb536fSBarry Smith .keywords: matrix, add
251d4bb536fSBarry Smith 
2522860a424SLois Curfman McInnes .seealso: MatAXPY()
253d4bb536fSBarry Smith  @*/
2547087cfbeSBarry Smith PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
255d4bb536fSBarry Smith {
25687828ca2SBarry Smith   PetscScalar    one = 1.0;
2576849ba73SBarry Smith   PetscErrorCode ierr;
25838baddfdSBarry Smith   PetscInt       mX,mY,nX,nY;
259d4bb536fSBarry Smith 
2603a40ed3dSBarry Smith   PetscFunctionBegin;
261c5eb9154SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
2620700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
263c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
264329f5518SBarry Smith   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
265329f5518SBarry Smith   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
266e32f2f54SBarry 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);
267d4bb536fSBarry Smith 
268f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
269cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
2703a40ed3dSBarry Smith   PetscFunctionReturn(0);
271d4bb536fSBarry Smith }
272b0a32e0cSBarry Smith 
273b0a32e0cSBarry Smith /*@
274b0a32e0cSBarry Smith     MatComputeExplicitOperator - Computes the explicit matrix
275b0a32e0cSBarry Smith 
276b0a32e0cSBarry Smith     Collective on Mat
277b0a32e0cSBarry Smith 
278b0a32e0cSBarry Smith     Input Parameter:
279b0a32e0cSBarry Smith .   inmat - the matrix
280b0a32e0cSBarry Smith 
281b0a32e0cSBarry Smith     Output Parameter:
282f3b1f45cSBarry Smith .   mat - the explict  operator
283b0a32e0cSBarry Smith 
284b0a32e0cSBarry Smith     Notes:
285b0a32e0cSBarry Smith     This computation is done by applying the operators to columns of the
286b0a32e0cSBarry Smith     identity matrix.
287b0a32e0cSBarry Smith 
288b0a32e0cSBarry Smith     Currently, this routine uses a dense matrix format when 1 processor
289b0a32e0cSBarry Smith     is used and a sparse format otherwise.  This routine is costly in general,
290b0a32e0cSBarry Smith     and is recommended for use only with relatively small systems.
291b0a32e0cSBarry Smith 
292b0a32e0cSBarry Smith     Level: advanced
293b0a32e0cSBarry Smith 
294b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator
295b0a32e0cSBarry Smith @*/
2967087cfbeSBarry Smith PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
297b0a32e0cSBarry Smith {
298b0a32e0cSBarry Smith   Vec            in,out;
299dfbe8321SBarry Smith   PetscErrorCode ierr;
3000b12b109SJed Brown   PetscInt       i,m,n,M,N,*rows,start,end;
301b0a32e0cSBarry Smith   MPI_Comm       comm;
30287828ca2SBarry Smith   PetscScalar    *array,zero = 0.0,one = 1.0;
30338baddfdSBarry Smith   PetscMPIInt    size;
304b0a32e0cSBarry Smith 
305b0a32e0cSBarry Smith   PetscFunctionBegin;
3060700a824SBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
3074482741eSBarry Smith   PetscValidPointer(mat,2);
308b0a32e0cSBarry Smith 
309ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
310b0a32e0cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
311b0a32e0cSBarry Smith 
3120b12b109SJed Brown   ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr);
3130b12b109SJed Brown   ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr);
3142a7a6963SBarry Smith   ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr);
3150b12b109SJed Brown   ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
3160b12b109SJed Brown   ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr);
317785e854fSJed Brown   ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr);
3188865f1eaSKarl Rupp   for (i=0; i<m; i++) rows[i] = start + i;
319b0a32e0cSBarry Smith 
320f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
3210b12b109SJed Brown   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
322b0a32e0cSBarry Smith   if (size == 1) {
323be5d1d56SKris Buschelman     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
3240298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
325b0a32e0cSBarry Smith   } else {
326be5d1d56SKris Buschelman     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
3270298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr);
328b0a32e0cSBarry Smith   }
329b0a32e0cSBarry Smith 
3300b12b109SJed Brown   for (i=0; i<N; i++) {
331b0a32e0cSBarry Smith 
3322dcb1b2aSMatthew Knepley     ierr = VecSet(in,zero);CHKERRQ(ierr);
333b0a32e0cSBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
334b0a32e0cSBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
335b0a32e0cSBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
336b0a32e0cSBarry Smith 
337b0a32e0cSBarry Smith     ierr = MatMult(inmat,in,out);CHKERRQ(ierr);
338b0a32e0cSBarry Smith 
339b0a32e0cSBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
340b0a32e0cSBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
341b0a32e0cSBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
342b0a32e0cSBarry Smith 
343b0a32e0cSBarry Smith   }
344b0a32e0cSBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
3456bf464f9SBarry Smith   ierr = VecDestroy(&out);CHKERRQ(ierr);
3466bf464f9SBarry Smith   ierr = VecDestroy(&in);CHKERRQ(ierr);
347b0a32e0cSBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
348b0a32e0cSBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34924f910e3SHong Zhang   PetscFunctionReturn(0);
35024f910e3SHong Zhang }
3514325cce7SMatthew G Knepley 
3524325cce7SMatthew G Knepley /*@
353f3b1f45cSBarry Smith     MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of
354f3b1f45cSBarry Smith         a give matrix that can apply MatMultTranspose()
355f3b1f45cSBarry Smith 
356f3b1f45cSBarry Smith     Collective on Mat
357f3b1f45cSBarry Smith 
358f3b1f45cSBarry Smith     Input Parameter:
359f3b1f45cSBarry Smith .   inmat - the matrix
360f3b1f45cSBarry Smith 
361f3b1f45cSBarry Smith     Output Parameter:
362f3b1f45cSBarry Smith .   mat - the explict  operator transposed
363f3b1f45cSBarry Smith 
364f3b1f45cSBarry Smith     Notes:
365f3b1f45cSBarry Smith     This computation is done by applying the transpose of the operator to columns of the
366f3b1f45cSBarry Smith     identity matrix.
367f3b1f45cSBarry Smith 
368f3b1f45cSBarry Smith     Currently, this routine uses a dense matrix format when 1 processor
369f3b1f45cSBarry Smith     is used and a sparse format otherwise.  This routine is costly in general,
370f3b1f45cSBarry Smith     and is recommended for use only with relatively small systems.
371f3b1f45cSBarry Smith 
372f3b1f45cSBarry Smith     Level: advanced
373f3b1f45cSBarry Smith 
374f3b1f45cSBarry Smith .keywords: Mat, compute, explicit, operator
375f3b1f45cSBarry Smith @*/
376f3b1f45cSBarry Smith PetscErrorCode  MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat)
377f3b1f45cSBarry Smith {
378f3b1f45cSBarry Smith   Vec            in,out;
379f3b1f45cSBarry Smith   PetscErrorCode ierr;
380f3b1f45cSBarry Smith   PetscInt       i,m,n,M,N,*rows,start,end;
381f3b1f45cSBarry Smith   MPI_Comm       comm;
382f3b1f45cSBarry Smith   PetscScalar    *array,zero = 0.0,one = 1.0;
383f3b1f45cSBarry Smith   PetscMPIInt    size;
384f3b1f45cSBarry Smith 
385f3b1f45cSBarry Smith   PetscFunctionBegin;
386f3b1f45cSBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
387f3b1f45cSBarry Smith   PetscValidPointer(mat,2);
388f3b1f45cSBarry Smith 
389f3b1f45cSBarry Smith   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
390f3b1f45cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
391f3b1f45cSBarry Smith 
392f3b1f45cSBarry Smith   ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr);
393f3b1f45cSBarry Smith   ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr);
394f3b1f45cSBarry Smith   ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr);
395f3b1f45cSBarry Smith   ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
396f3b1f45cSBarry Smith   ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr);
397f3b1f45cSBarry Smith   ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr);
398f3b1f45cSBarry Smith   for (i=0; i<m; i++) rows[i] = start + i;
399f3b1f45cSBarry Smith 
400f3b1f45cSBarry Smith   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
401f3b1f45cSBarry Smith   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
402f3b1f45cSBarry Smith   if (size == 1) {
403f3b1f45cSBarry Smith     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
404f3b1f45cSBarry Smith     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
405f3b1f45cSBarry Smith   } else {
406f3b1f45cSBarry Smith     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
407f3b1f45cSBarry Smith     ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr);
408f3b1f45cSBarry Smith   }
409f3b1f45cSBarry Smith 
410f3b1f45cSBarry Smith   for (i=0; i<N; i++) {
411f3b1f45cSBarry Smith 
412f3b1f45cSBarry Smith     ierr = VecSet(in,zero);CHKERRQ(ierr);
413f3b1f45cSBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
414f3b1f45cSBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
415f3b1f45cSBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
416f3b1f45cSBarry Smith 
417f3b1f45cSBarry Smith     ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr);
418f3b1f45cSBarry Smith 
419f3b1f45cSBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
420f3b1f45cSBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
421f3b1f45cSBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
422f3b1f45cSBarry Smith 
423f3b1f45cSBarry Smith   }
424f3b1f45cSBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
425f3b1f45cSBarry Smith   ierr = VecDestroy(&out);CHKERRQ(ierr);
426f3b1f45cSBarry Smith   ierr = VecDestroy(&in);CHKERRQ(ierr);
427f3b1f45cSBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
428f3b1f45cSBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
429f3b1f45cSBarry Smith   PetscFunctionReturn(0);
430f3b1f45cSBarry Smith }
431f3b1f45cSBarry Smith 
432f3b1f45cSBarry Smith /*@
4334325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
4344325cce7SMatthew G Knepley 
4354325cce7SMatthew G Knepley   Input Parameters:
4364325cce7SMatthew G Knepley + A   - The matrix
4374325cce7SMatthew G Knepley - tol - The zero tolerance
4384325cce7SMatthew G Knepley 
4394325cce7SMatthew G Knepley   Output Parameters:
4404325cce7SMatthew G Knepley . A - The chopped matrix
4414325cce7SMatthew G Knepley 
4424325cce7SMatthew G Knepley   Level: intermediate
4434325cce7SMatthew G Knepley 
4444325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries()
4453fc99919SSatish Balay  @*/
4464325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol)
4474325cce7SMatthew G Knepley {
4484325cce7SMatthew G Knepley   PetscScalar    *newVals;
4494325cce7SMatthew G Knepley   PetscInt       *newCols;
4504325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
4514325cce7SMatthew G Knepley   PetscErrorCode ierr;
4524325cce7SMatthew G Knepley 
4534325cce7SMatthew G Knepley   PetscFunctionBegin;
4544325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4554325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4564325cce7SMatthew G Knepley     PetscInt ncols;
4574325cce7SMatthew G Knepley 
4580298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4594325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4600298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4614325cce7SMatthew G Knepley   }
4624325cce7SMatthew G Knepley   numRows = rEnd - rStart;
463b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
464dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4654325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4664325cce7SMatthew G Knepley     const PetscScalar *vals;
4674325cce7SMatthew G Knepley     const PetscInt    *cols;
468fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
4694325cce7SMatthew G Knepley 
4704325cce7SMatthew G Knepley     if (r < rEnd) {
4714325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
4724325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
4734325cce7SMatthew G Knepley         newCols[c] = cols[c];
4744325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
4754325cce7SMatthew G Knepley       }
476fad45679SMatthew G. Knepley       newcols = ncols;
4774325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
478fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
4794325cce7SMatthew G Knepley     }
4804325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4814325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4824325cce7SMatthew G Knepley   }
4834325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
4844325cce7SMatthew G Knepley   PetscFunctionReturn(0);
4854325cce7SMatthew G Knepley }
486