xref: /petsc/src/mat/utils/axpy.c (revision 1c738b47add20c2cbe69d056fd94a20e1fa7e005)
16f79c3a4SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/
36f79c3a4SBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "MatAXPY"
606be10caSBarry Smith /*@
721c89e3eSBarry Smith    MatAXPY - Computes Y = a*X + Y.
86f79c3a4SBarry Smith 
93f9fe445SBarry Smith    Logically  Collective on Mat
10fee21e36SBarry Smith 
1198a79cdbSBarry Smith    Input Parameters:
12607cd303SBarry Smith +  a - the scalar multiplier
13607cd303SBarry Smith .  X - the first matrix
14607cd303SBarry Smith .  Y - the second matrix
15407f6b05SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN
16407f6b05SHong Zhang          or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
1798a79cdbSBarry Smith 
182860a424SLois Curfman McInnes    Level: intermediate
192860a424SLois Curfman McInnes 
209cf4f1e8SLois Curfman McInnes .keywords: matrix, add
21d4bb536fSBarry Smith 
222860a424SLois Curfman McInnes .seealso: MatAYPX()
2306be10caSBarry Smith  @*/
247087cfbeSBarry Smith PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
256f79c3a4SBarry Smith {
266849ba73SBarry Smith   PetscErrorCode ierr;
27c1ac3661SBarry Smith   PetscInt       m1,m2,n1,n2;
286f79c3a4SBarry Smith 
293a40ed3dSBarry Smith   PetscFunctionBegin;
300700a824SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
310700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
32c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
33273d9f13SBarry Smith   ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr);
34273d9f13SBarry Smith   ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr);
35e32f2f54SBarry 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);
361987afe7SBarry Smith 
37e8136da8SHong Zhang   ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
38ab784542SHong Zhang   if (Y->ops->axpy) {
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);
442eea6e16SBarry Smith #if defined(PETSC_HAVE_CUSP)
452eea6e16SBarry Smith   if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
462eea6e16SBarry Smith     Y->valid_GPU_matrix = PETSC_CUSP_CPU;
472eea6e16SBarry Smith   }
48c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VIENNACL)
49d67ff14aSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) {
50d67ff14aSKarl Rupp     Y->valid_GPU_matrix = PETSC_VIENNACL_CPU;
51d67ff14aSKarl Rupp   }
52c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VECCUDA)
53c41cb2e2SAlejandro Lamas Daviña   if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) {
54c41cb2e2SAlejandro Lamas Daviña     Y->valid_GPU_matrix = PETSC_CUDA_CPU;
55c41cb2e2SAlejandro Lamas Daviña   }
56d67ff14aSKarl Rupp #endif
57607cd303SBarry Smith   PetscFunctionReturn(0);
58607cd303SBarry Smith }
59607cd303SBarry Smith 
60607cd303SBarry Smith #undef __FUNCT__
61607cd303SBarry Smith #define __FUNCT__ "MatAXPY_Basic"
62f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
63607cd303SBarry Smith {
6438baddfdSBarry Smith   PetscInt          i,start,end,j,ncols,m,n;
656849ba73SBarry Smith   PetscErrorCode    ierr;
6638baddfdSBarry Smith   const PetscInt    *row;
67b3cc6726SBarry Smith   PetscScalar       *val;
68b3cc6726SBarry Smith   const PetscScalar *vals;
69607cd303SBarry Smith 
70607cd303SBarry Smith   PetscFunctionBegin;
718dadbd76SSatish Balay   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
7290f02eecSBarry Smith   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
73f4df32b1SMatthew Knepley   if (a == 1.0) {
74d4bb536fSBarry Smith     for (i = start; i < end; i++) {
75d4bb536fSBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
76d4bb536fSBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
77d4bb536fSBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
78d4bb536fSBarry Smith     }
79d4bb536fSBarry Smith   } else {
80854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr);
8106be10caSBarry Smith     for (i=start; i<end; i++) {
82b3cc6726SBarry Smith       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
8306be10caSBarry Smith       for (j=0; j<ncols; j++) {
84f4df32b1SMatthew Knepley         val[j] = a*vals[j];
856f79c3a4SBarry Smith       }
86b3cc6726SBarry Smith       ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
87b3cc6726SBarry Smith       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
886f79c3a4SBarry Smith     }
89b3cc6726SBarry Smith     ierr = PetscFree(val);CHKERRQ(ierr);
90d4bb536fSBarry Smith   }
916d4a8577SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
926d4a8577SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
933a40ed3dSBarry Smith   PetscFunctionReturn(0);
946f79c3a4SBarry Smith }
95052efed2SBarry Smith 
964a2ae208SSatish Balay #undef __FUNCT__
97ec7775f6SShri Abhyankar #define __FUNCT__ "MatAXPY_BasicWithPreallocation"
98ec7775f6SShri Abhyankar PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
99ec7775f6SShri Abhyankar {
100ec7775f6SShri Abhyankar   PetscInt          i,start,end,j,ncols,m,n;
101ec7775f6SShri Abhyankar   PetscErrorCode    ierr;
102ec7775f6SShri Abhyankar   const PetscInt    *row;
103ec7775f6SShri Abhyankar   PetscScalar       *val;
104ec7775f6SShri Abhyankar   const PetscScalar *vals;
105ec7775f6SShri Abhyankar 
106ec7775f6SShri Abhyankar   PetscFunctionBegin;
107ec7775f6SShri Abhyankar   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
108ec7775f6SShri Abhyankar   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
109ec7775f6SShri Abhyankar   if (a == 1.0) {
110ec7775f6SShri Abhyankar     for (i = start; i < end; i++) {
111ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
112ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
113ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
114ec7775f6SShri Abhyankar 
115ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
116ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
117ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
118ec7775f6SShri Abhyankar     }
119ec7775f6SShri Abhyankar   } else {
120854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr);
121ec7775f6SShri Abhyankar     for (i=start; i<end; i++) {
122ec7775f6SShri Abhyankar       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
123ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
124ec7775f6SShri Abhyankar       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
125ec7775f6SShri Abhyankar 
126ec7775f6SShri Abhyankar       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
127ec7775f6SShri Abhyankar       for (j=0; j<ncols; j++) {
128ec7775f6SShri Abhyankar         val[j] = a*vals[j];
129ec7775f6SShri Abhyankar       }
130ec7775f6SShri Abhyankar       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
131ec7775f6SShri Abhyankar       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
132ec7775f6SShri Abhyankar     }
133ec7775f6SShri Abhyankar     ierr = PetscFree(val);CHKERRQ(ierr);
134ec7775f6SShri Abhyankar   }
135ec7775f6SShri Abhyankar   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136ec7775f6SShri Abhyankar   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137ec7775f6SShri Abhyankar   PetscFunctionReturn(0);
138ec7775f6SShri Abhyankar }
139ec7775f6SShri Abhyankar 
140ec7775f6SShri Abhyankar #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "MatShift"
142052efed2SBarry Smith /*@
14387828ca2SBarry Smith    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
144052efed2SBarry Smith 
1453f9fe445SBarry Smith    Neighbor-wise Collective on Mat
146fee21e36SBarry Smith 
14798a79cdbSBarry Smith    Input Parameters:
14898a79cdbSBarry Smith +  Y - the matrices
14987828ca2SBarry Smith -  a - the PetscScalar
15098a79cdbSBarry Smith 
1512860a424SLois Curfman McInnes    Level: intermediate
1522860a424SLois Curfman McInnes 
1536f33a894SBarry Smith    Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
1546f33a894SBarry 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
1556f33a894SBarry Smith    entry).
1566f33a894SBarry Smith 
1576f33a894SBarry Smith    Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is
1586f33a894SBarry Smith     preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix.
1596f33a894SBarry Smith 
160052efed2SBarry Smith .keywords: matrix, add, shift
1616b9ee512SLois Curfman McInnes 
162f56f2b3fSBarry Smith .seealso: MatDiagonalSet()
163052efed2SBarry Smith  @*/
1647087cfbeSBarry Smith PetscErrorCode  MatShift(Mat Y,PetscScalar a)
165052efed2SBarry Smith {
1666849ba73SBarry Smith   PetscErrorCode ierr;
167052efed2SBarry Smith 
1683a40ed3dSBarry Smith   PetscFunctionBegin;
1690700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
170ce94432eSBarry Smith   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
171ce94432eSBarry Smith   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1724994cf47SJed Brown   MatCheckPreallocated(Y,1);
173b50b34bdSBarry Smith 
174*1c738b47SStefano Zampini   if (Y->ops->shift) {
175f4df32b1SMatthew Knepley     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
176*1c738b47SStefano Zampini   } else {
177*1c738b47SStefano Zampini     ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
178*1c738b47SStefano Zampini   }
1797d68702bSBarry Smith 
1802eea6e16SBarry Smith #if defined(PETSC_HAVE_CUSP)
1812eea6e16SBarry Smith   if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
1822eea6e16SBarry Smith     Y->valid_GPU_matrix = PETSC_CUSP_CPU;
1832eea6e16SBarry Smith   }
184c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VIENNACL)
185d67ff14aSKarl Rupp   if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) {
186d67ff14aSKarl Rupp     Y->valid_GPU_matrix = PETSC_VIENNACL_CPU;
187d67ff14aSKarl Rupp   }
188c41cb2e2SAlejandro Lamas Daviña #elif defined(PETSC_HAVE_VECCUDA)
189c41cb2e2SAlejandro Lamas Daviña   if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) {
190c41cb2e2SAlejandro Lamas Daviña     Y->valid_GPU_matrix = PETSC_CUDA_CPU;
191c41cb2e2SAlejandro Lamas Daviña   }
192d67ff14aSKarl Rupp #endif
1933a40ed3dSBarry Smith   PetscFunctionReturn(0);
194052efed2SBarry Smith }
1956d84be18SBarry Smith 
1964a2ae208SSatish Balay #undef __FUNCT__
19709f38230SBarry Smith #define __FUNCT__ "MatDiagonalSet_Default"
1987087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
19909f38230SBarry Smith {
20009f38230SBarry Smith   PetscErrorCode ierr;
20167576b8bSBarry Smith   PetscInt       i,start,end;
20209f38230SBarry Smith   PetscScalar    *v;
20309f38230SBarry Smith 
20409f38230SBarry Smith   PetscFunctionBegin;
20509f38230SBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
20609f38230SBarry Smith   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
20709f38230SBarry Smith   for (i=start; i<end; i++) {
20809f38230SBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
20909f38230SBarry Smith   }
21009f38230SBarry Smith   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
21109f38230SBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
21209f38230SBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
21309f38230SBarry Smith   PetscFunctionReturn(0);
21409f38230SBarry Smith }
21509f38230SBarry Smith 
21609f38230SBarry Smith #undef __FUNCT__
2174a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalSet"
2186d84be18SBarry Smith /*@
219f56f2b3fSBarry Smith    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
220f56f2b3fSBarry Smith    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
221f56f2b3fSBarry Smith    INSERT_VALUES.
2226d84be18SBarry Smith 
2236d84be18SBarry Smith    Input Parameters:
22498a79cdbSBarry Smith +  Y - the input matrix
225f56f2b3fSBarry Smith .  D - the diagonal matrix, represented as a vector
226f56f2b3fSBarry Smith -  i - INSERT_VALUES or ADD_VALUES
2276d84be18SBarry Smith 
2283f9fe445SBarry Smith    Neighbor-wise Collective on Mat and Vec
229fee21e36SBarry Smith 
2306f33a894SBarry Smith    Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
2316f33a894SBarry 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
2326f33a894SBarry Smith    entry).
2336f33a894SBarry Smith 
2342860a424SLois Curfman McInnes    Level: intermediate
2352860a424SLois Curfman McInnes 
2366b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal
2376b9ee512SLois Curfman McInnes 
2386b9ee512SLois Curfman McInnes .seealso: MatShift()
2396d84be18SBarry Smith @*/
2407087cfbeSBarry Smith PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
2416d84be18SBarry Smith {
2426849ba73SBarry Smith   PetscErrorCode ierr;
24367576b8bSBarry Smith   PetscInt       matlocal,veclocal;
2446d84be18SBarry Smith 
2453a40ed3dSBarry Smith   PetscFunctionBegin;
2460700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
2470700a824SBarry Smith   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
24867576b8bSBarry Smith   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
24967576b8bSBarry Smith   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
25067576b8bSBarry 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);
251f56f2b3fSBarry Smith   if (Y->ops->diagonalset) {
252f56f2b3fSBarry Smith     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
25394d884c6SBarry Smith   } else {
25409f38230SBarry Smith     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
2556d84be18SBarry Smith   }
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2576d84be18SBarry Smith }
258d4bb536fSBarry Smith 
2594a2ae208SSatish Balay #undef __FUNCT__
2604a2ae208SSatish Balay #define __FUNCT__ "MatAYPX"
261d4bb536fSBarry Smith /*@
26204aac2b0SHong Zhang    MatAYPX - Computes Y = a*Y + X.
263d4bb536fSBarry Smith 
2643f9fe445SBarry Smith    Logically on Mat
265fee21e36SBarry Smith 
26698a79cdbSBarry Smith    Input Parameters:
26704aac2b0SHong Zhang +  a - the PetscScalar multiplier
26804aac2b0SHong Zhang .  Y - the first matrix
26904aac2b0SHong Zhang .  X - the second matrix
27004aac2b0SHong Zhang -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
271d4bb536fSBarry Smith 
2722860a424SLois Curfman McInnes    Level: intermediate
2732860a424SLois Curfman McInnes 
274d4bb536fSBarry Smith .keywords: matrix, add
275d4bb536fSBarry Smith 
2762860a424SLois Curfman McInnes .seealso: MatAXPY()
277d4bb536fSBarry Smith  @*/
2787087cfbeSBarry Smith PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
279d4bb536fSBarry Smith {
28087828ca2SBarry Smith   PetscScalar    one = 1.0;
2816849ba73SBarry Smith   PetscErrorCode ierr;
28238baddfdSBarry Smith   PetscInt       mX,mY,nX,nY;
283d4bb536fSBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
285c5eb9154SBarry Smith   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
2860700a824SBarry Smith   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
287c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(Y,a,2);
288329f5518SBarry Smith   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
289329f5518SBarry Smith   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
290e32f2f54SBarry 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);
291d4bb536fSBarry Smith 
292f4df32b1SMatthew Knepley   ierr = MatScale(Y,a);CHKERRQ(ierr);
293cb9801acSJed Brown   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
2943a40ed3dSBarry Smith   PetscFunctionReturn(0);
295d4bb536fSBarry Smith }
296b0a32e0cSBarry Smith 
2974a2ae208SSatish Balay #undef __FUNCT__
2984a2ae208SSatish Balay #define __FUNCT__ "MatComputeExplicitOperator"
299b0a32e0cSBarry Smith /*@
300b0a32e0cSBarry Smith     MatComputeExplicitOperator - Computes the explicit matrix
301b0a32e0cSBarry Smith 
302b0a32e0cSBarry Smith     Collective on Mat
303b0a32e0cSBarry Smith 
304b0a32e0cSBarry Smith     Input Parameter:
305b0a32e0cSBarry Smith .   inmat - the matrix
306b0a32e0cSBarry Smith 
307b0a32e0cSBarry Smith     Output Parameter:
308b0a32e0cSBarry Smith .   mat - the explict preconditioned operator
309b0a32e0cSBarry Smith 
310b0a32e0cSBarry Smith     Notes:
311b0a32e0cSBarry Smith     This computation is done by applying the operators to columns of the
312b0a32e0cSBarry Smith     identity matrix.
313b0a32e0cSBarry Smith 
314b0a32e0cSBarry Smith     Currently, this routine uses a dense matrix format when 1 processor
315b0a32e0cSBarry Smith     is used and a sparse format otherwise.  This routine is costly in general,
316b0a32e0cSBarry Smith     and is recommended for use only with relatively small systems.
317b0a32e0cSBarry Smith 
318b0a32e0cSBarry Smith     Level: advanced
319b0a32e0cSBarry Smith 
320b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator
321b0a32e0cSBarry Smith @*/
3227087cfbeSBarry Smith PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
323b0a32e0cSBarry Smith {
324b0a32e0cSBarry Smith   Vec            in,out;
325dfbe8321SBarry Smith   PetscErrorCode ierr;
3260b12b109SJed Brown   PetscInt       i,m,n,M,N,*rows,start,end;
327b0a32e0cSBarry Smith   MPI_Comm       comm;
32887828ca2SBarry Smith   PetscScalar    *array,zero = 0.0,one = 1.0;
32938baddfdSBarry Smith   PetscMPIInt    size;
330b0a32e0cSBarry Smith 
331b0a32e0cSBarry Smith   PetscFunctionBegin;
3320700a824SBarry Smith   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
3334482741eSBarry Smith   PetscValidPointer(mat,2);
334b0a32e0cSBarry Smith 
335ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
336b0a32e0cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
337b0a32e0cSBarry Smith 
3380b12b109SJed Brown   ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr);
3390b12b109SJed Brown   ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr);
3402a7a6963SBarry Smith   ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr);
3410b12b109SJed Brown   ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
3420b12b109SJed Brown   ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr);
343785e854fSJed Brown   ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr);
3448865f1eaSKarl Rupp   for (i=0; i<m; i++) rows[i] = start + i;
345b0a32e0cSBarry Smith 
346f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
3470b12b109SJed Brown   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
348b0a32e0cSBarry Smith   if (size == 1) {
349be5d1d56SKris Buschelman     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
3500298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
351b0a32e0cSBarry Smith   } else {
352be5d1d56SKris Buschelman     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
3530298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr);
354b0a32e0cSBarry Smith   }
355b0a32e0cSBarry Smith 
3560b12b109SJed Brown   for (i=0; i<N; i++) {
357b0a32e0cSBarry Smith 
3582dcb1b2aSMatthew Knepley     ierr = VecSet(in,zero);CHKERRQ(ierr);
359b0a32e0cSBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
360b0a32e0cSBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
361b0a32e0cSBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
362b0a32e0cSBarry Smith 
363b0a32e0cSBarry Smith     ierr = MatMult(inmat,in,out);CHKERRQ(ierr);
364b0a32e0cSBarry Smith 
365b0a32e0cSBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
366b0a32e0cSBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
367b0a32e0cSBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
368b0a32e0cSBarry Smith 
369b0a32e0cSBarry Smith   }
370b0a32e0cSBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
3716bf464f9SBarry Smith   ierr = VecDestroy(&out);CHKERRQ(ierr);
3726bf464f9SBarry Smith   ierr = VecDestroy(&in);CHKERRQ(ierr);
373b0a32e0cSBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
374b0a32e0cSBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37524f910e3SHong Zhang   PetscFunctionReturn(0);
37624f910e3SHong Zhang }
3774325cce7SMatthew G Knepley 
3784325cce7SMatthew G Knepley #undef __FUNCT__
3794325cce7SMatthew G Knepley #define __FUNCT__ "MatChop"
3804325cce7SMatthew G Knepley /*@
3814325cce7SMatthew G Knepley   MatChop - Set all values in the matrix less than the tolerance to zero
3824325cce7SMatthew G Knepley 
3834325cce7SMatthew G Knepley   Input Parameters:
3844325cce7SMatthew G Knepley + A   - The matrix
3854325cce7SMatthew G Knepley - tol - The zero tolerance
3864325cce7SMatthew G Knepley 
3874325cce7SMatthew G Knepley   Output Parameters:
3884325cce7SMatthew G Knepley . A - The chopped matrix
3894325cce7SMatthew G Knepley 
3904325cce7SMatthew G Knepley   Level: intermediate
3914325cce7SMatthew G Knepley 
3924325cce7SMatthew G Knepley .seealso: MatCreate(), MatZeroEntries()
3933fc99919SSatish Balay  @*/
3944325cce7SMatthew G Knepley PetscErrorCode MatChop(Mat A, PetscReal tol)
3954325cce7SMatthew G Knepley {
3964325cce7SMatthew G Knepley   PetscScalar    *newVals;
3974325cce7SMatthew G Knepley   PetscInt       *newCols;
3984325cce7SMatthew G Knepley   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
3994325cce7SMatthew G Knepley   PetscErrorCode ierr;
4004325cce7SMatthew G Knepley 
4014325cce7SMatthew G Knepley   PetscFunctionBegin;
4024325cce7SMatthew G Knepley   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
4034325cce7SMatthew G Knepley   for (r = rStart; r < rEnd; ++r) {
4044325cce7SMatthew G Knepley     PetscInt ncols;
4054325cce7SMatthew G Knepley 
4060298fd71SBarry Smith     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4074325cce7SMatthew G Knepley     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
4080298fd71SBarry Smith     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
4094325cce7SMatthew G Knepley   }
4104325cce7SMatthew G Knepley   numRows = rEnd - rStart;
411b2566f29SBarry Smith   ierr    = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
412dcca6d9dSJed Brown   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
4134325cce7SMatthew G Knepley   for (r = rStart; r < rStart+maxRows; ++r) {
4144325cce7SMatthew G Knepley     const PetscScalar *vals;
4154325cce7SMatthew G Knepley     const PetscInt    *cols;
416fad45679SMatthew G. Knepley     PetscInt           ncols, newcols, c;
4174325cce7SMatthew G Knepley 
4184325cce7SMatthew G Knepley     if (r < rEnd) {
4194325cce7SMatthew G Knepley       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
4204325cce7SMatthew G Knepley       for (c = 0; c < ncols; ++c) {
4214325cce7SMatthew G Knepley         newCols[c] = cols[c];
4224325cce7SMatthew G Knepley         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
4234325cce7SMatthew G Knepley       }
424fad45679SMatthew G. Knepley       newcols = ncols;
4254325cce7SMatthew G Knepley       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
426fad45679SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
4274325cce7SMatthew G Knepley     }
4284325cce7SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4294325cce7SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4304325cce7SMatthew G Knepley   }
4314325cce7SMatthew G Knepley   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
4324325cce7SMatthew G Knepley   PetscFunctionReturn(0);
4334325cce7SMatthew G Knepley }
434