xref: /petsc/src/mat/utils/axpy.c (revision bcab024551cbaaaa7f49e357634a3584f91cb6d5)
1 
2 #include <petsc-private/matimpl.h>  /*I   "petscmat.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "MatAXPY"
6 /*@
7    MatAXPY - Computes Y = a*X + Y.
8 
9    Logically  Collective on Mat
10 
11    Input Parameters:
12 +  a - the scalar multiplier
13 .  X - the first matrix
14 .  Y - the second matrix
15 -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN
16          or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)
17 
18    Level: intermediate
19 
20 .keywords: matrix, add
21 
22 .seealso: MatAYPX()
23  @*/
24 PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
25 {
26   PetscErrorCode ierr;
27   PetscInt       m1,m2,n1,n2;
28 
29   PetscFunctionBegin;
30   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
31   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
32   PetscValidLogicalCollectiveScalar(Y,a,2);
33   ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr);
34   ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr);
35   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);
36 
37   ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
38   if (Y->ops->axpy) {
39     ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr);
40   } else {
41     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
42   }
43   ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr);
44 #if defined(PETSC_HAVE_CUSP)
45   if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
46     Y->valid_GPU_matrix = PETSC_CUSP_CPU;
47   }
48 #endif
49 #if defined(PETSC_HAVE_VIENNACL)
50   if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) {
51     Y->valid_GPU_matrix = PETSC_VIENNACL_CPU;
52   }
53 #endif
54   PetscFunctionReturn(0);
55 }
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "MatAXPY_Basic"
59 PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
60 {
61   PetscInt          i,start,end,j,ncols,m,n;
62   PetscErrorCode    ierr;
63   const PetscInt    *row;
64   PetscScalar       *val;
65   const PetscScalar *vals;
66 
67   PetscFunctionBegin;
68   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
69   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
70   if (a == 1.0) {
71     for (i = start; i < end; i++) {
72       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
73       ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
74       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
75     }
76   } else {
77     ierr = PetscMalloc1((n+1),&val);CHKERRQ(ierr);
78     for (i=start; i<end; i++) {
79       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
80       for (j=0; j<ncols; j++) {
81         val[j] = a*vals[j];
82       }
83       ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
84       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
85     }
86     ierr = PetscFree(val);CHKERRQ(ierr);
87   }
88   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "MatAXPY_BasicWithPreallocation"
95 PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str)
96 {
97   PetscInt          i,start,end,j,ncols,m,n;
98   PetscErrorCode    ierr;
99   const PetscInt    *row;
100   PetscScalar       *val;
101   const PetscScalar *vals;
102 
103   PetscFunctionBegin;
104   ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr);
105   ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
106   if (a == 1.0) {
107     for (i = start; i < end; i++) {
108       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
109       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
110       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
111 
112       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
113       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
114       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
115     }
116   } else {
117     ierr = PetscMalloc1((n+1),&val);CHKERRQ(ierr);
118     for (i=start; i<end; i++) {
119       ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
120       ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr);
121       ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr);
122 
123       ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
124       for (j=0; j<ncols; j++) {
125         val[j] = a*vals[j];
126       }
127       ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr);
128       ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr);
129     }
130     ierr = PetscFree(val);CHKERRQ(ierr);
131   }
132   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134   PetscFunctionReturn(0);
135 }
136 
137 #undef __FUNCT__
138 #define __FUNCT__ "MatShift"
139 /*@
140    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.
141 
142    Neighbor-wise Collective on Mat
143 
144    Input Parameters:
145 +  Y - the matrices
146 -  a - the PetscScalar
147 
148    Level: intermediate
149 
150 .keywords: matrix, add, shift
151 
152 .seealso: MatDiagonalSet()
153  @*/
154 PetscErrorCode  MatShift(Mat Y,PetscScalar a)
155 {
156   PetscErrorCode ierr;
157   PetscInt       i,start,end;
158 
159   PetscFunctionBegin;
160   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
161   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
162   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
163   MatCheckPreallocated(Y,1);
164 
165   if (Y->ops->shift) {
166     ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr);
167   } else {
168     PetscScalar alpha = a;
169     ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
170     for (i=start; i<end; i++) {
171       ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr);
172     }
173     ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174     ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175   }
176 #if defined(PETSC_HAVE_CUSP)
177   if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
178     Y->valid_GPU_matrix = PETSC_CUSP_CPU;
179   }
180 #endif
181 #if defined(PETSC_HAVE_VIENNACL)
182   if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) {
183     Y->valid_GPU_matrix = PETSC_VIENNACL_CPU;
184   }
185 #endif
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "MatDiagonalSet_Default"
191 PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
192 {
193   PetscErrorCode ierr;
194   PetscInt       i,start,end;
195   PetscScalar    *v;
196 
197   PetscFunctionBegin;
198   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
199   ierr = VecGetArray(D,&v);CHKERRQ(ierr);
200   for (i=start; i<end; i++) {
201     ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr);
202   }
203   ierr = VecRestoreArray(D,&v);CHKERRQ(ierr);
204   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
205   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "MatDiagonalSet"
211 /*@
212    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
213    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
214    INSERT_VALUES.
215 
216    Input Parameters:
217 +  Y - the input matrix
218 .  D - the diagonal matrix, represented as a vector
219 -  i - INSERT_VALUES or ADD_VALUES
220 
221    Neighbor-wise Collective on Mat and Vec
222 
223    Level: intermediate
224 
225 .keywords: matrix, add, shift, diagonal
226 
227 .seealso: MatShift()
228 @*/
229 PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
230 {
231   PetscErrorCode ierr;
232   PetscInt       matlocal,veclocal;
233 
234   PetscFunctionBegin;
235   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
236   PetscValidHeaderSpecific(D,VEC_CLASSID,2);
237   ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr);
238   ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr);
239   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);
240   if (Y->ops->diagonalset) {
241     ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr);
242   } else {
243     ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr);
244   }
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "MatAYPX"
250 /*@
251    MatAYPX - Computes Y = a*Y + X.
252 
253    Logically on Mat
254 
255    Input Parameters:
256 +  a - the PetscScalar multiplier
257 .  Y - the first matrix
258 .  X - the second matrix
259 -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN
260 
261    Level: intermediate
262 
263 .keywords: matrix, add
264 
265 .seealso: MatAXPY()
266  @*/
267 PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
268 {
269   PetscScalar    one = 1.0;
270   PetscErrorCode ierr;
271   PetscInt       mX,mY,nX,nY;
272 
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(X,MAT_CLASSID,3);
275   PetscValidHeaderSpecific(Y,MAT_CLASSID,1);
276   PetscValidLogicalCollectiveScalar(Y,a,2);
277   ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr);
278   ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr);
279   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);
280 
281   ierr = MatScale(Y,a);CHKERRQ(ierr);
282   ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "MatComputeExplicitOperator"
288 /*@
289     MatComputeExplicitOperator - Computes the explicit matrix
290 
291     Collective on Mat
292 
293     Input Parameter:
294 .   inmat - the matrix
295 
296     Output Parameter:
297 .   mat - the explict preconditioned operator
298 
299     Notes:
300     This computation is done by applying the operators to columns of the
301     identity matrix.
302 
303     Currently, this routine uses a dense matrix format when 1 processor
304     is used and a sparse format otherwise.  This routine is costly in general,
305     and is recommended for use only with relatively small systems.
306 
307     Level: advanced
308 
309 .keywords: Mat, compute, explicit, operator
310 @*/
311 PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
312 {
313   Vec            in,out;
314   PetscErrorCode ierr;
315   PetscInt       i,m,n,M,N,*rows,start,end;
316   MPI_Comm       comm;
317   PetscScalar    *array,zero = 0.0,one = 1.0;
318   PetscMPIInt    size;
319 
320   PetscFunctionBegin;
321   PetscValidHeaderSpecific(inmat,MAT_CLASSID,1);
322   PetscValidPointer(mat,2);
323 
324   ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr);
325   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
326 
327   ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr);
328   ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr);
329   ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr);
330   ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
331   ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr);
332   ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr);
333   for (i=0; i<m; i++) rows[i] = start + i;
334 
335   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
336   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
337   if (size == 1) {
338     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
339     ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr);
340   } else {
341     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
342     ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr);
343   }
344 
345   for (i=0; i<N; i++) {
346 
347     ierr = VecSet(in,zero);CHKERRQ(ierr);
348     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
349     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
350     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
351 
352     ierr = MatMult(inmat,in,out);CHKERRQ(ierr);
353 
354     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
355     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
356     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
357 
358   }
359   ierr = PetscFree(rows);CHKERRQ(ierr);
360   ierr = VecDestroy(&out);CHKERRQ(ierr);
361   ierr = VecDestroy(&in);CHKERRQ(ierr);
362   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
363   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "MatChop"
369 /*@
370   MatChop - Set all values in the matrix less than the tolerance to zero
371 
372   Input Parameters:
373 + A   - The matrix
374 - tol - The zero tolerance
375 
376   Output Parameters:
377 . A - The chopped matrix
378 
379   Level: intermediate
380 
381 .seealso: MatCreate(), MatZeroEntries()
382  @*/
383 PetscErrorCode MatChop(Mat A, PetscReal tol)
384 {
385   PetscScalar    *newVals;
386   PetscInt       *newCols;
387   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
392   for (r = rStart; r < rEnd; ++r) {
393     PetscInt ncols;
394 
395     ierr   = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
396     colMax = PetscMax(colMax, ncols);CHKERRQ(ierr);
397     ierr   = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr);
398   }
399   numRows = rEnd - rStart;
400   ierr    = MPI_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
401   ierr    = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr);
402   for (r = rStart; r < rStart+maxRows; ++r) {
403     const PetscScalar *vals;
404     const PetscInt    *cols;
405     PetscInt           ncols, newcols, c;
406 
407     if (r < rEnd) {
408       ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
409       for (c = 0; c < ncols; ++c) {
410         newCols[c] = cols[c];
411         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
412       }
413       newcols = ncols;
414       ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr);
415       ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr);
416     }
417     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
418     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
419   }
420   ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr);
421   PetscFunctionReturn(0);
422 }
423