Actual source code: axpy.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/

  6: /*@
  7:    MatAXPY - Computes Y = a*X + Y.

  9:    Logically  Collective on Mat

 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)

 18:    Level: intermediate

 20: .keywords: matrix, add

 22: .seealso: MatAYPX()
 23:  @*/
 24: PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
 25: {
 27:   PetscInt       m1,m2,n1,n2;

 33:   MatGetSize(X,&m1,&n1);
 34:   MatGetSize(Y,&m2,&n2);
 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);

 37:   PetscLogEventBegin(MAT_AXPY,Y,0,0,0);
 38:   if (Y->ops->axpy) {
 39:     (*Y->ops->axpy)(Y,a,X,str);
 40:   } else {
 41:     MatAXPY_Basic(Y,a,X,str);
 42:   }
 43:   PetscLogEventEnd(MAT_AXPY,Y,0,0,0);
 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:   return(0);
 55: }

 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;

 68:   MatGetSize(X,&m,&n);
 69:   MatGetOwnershipRange(X,&start,&end);
 70:   if (a == 1.0) {
 71:     for (i = start; i < end; i++) {
 72:       MatGetRow(X,i,&ncols,&row,&vals);
 73:       MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);
 74:       MatRestoreRow(X,i,&ncols,&row,&vals);
 75:     }
 76:   } else {
 77:     PetscMalloc1(n+1,&val);
 78:     for (i=start; i<end; i++) {
 79:       MatGetRow(X,i,&ncols,&row,&vals);
 80:       for (j=0; j<ncols; j++) {
 81:         val[j] = a*vals[j];
 82:       }
 83:       MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);
 84:       MatRestoreRow(X,i,&ncols,&row,&vals);
 85:     }
 86:     PetscFree(val);
 87:   }
 88:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
 89:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
 90:   return(0);
 91: }

 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;

104:   MatGetSize(X,&m,&n);
105:   MatGetOwnershipRange(X,&start,&end);
106:   if (a == 1.0) {
107:     for (i = start; i < end; i++) {
108:       MatGetRow(Y,i,&ncols,&row,&vals);
109:       MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);
110:       MatRestoreRow(Y,i,&ncols,&row,&vals);

112:       MatGetRow(X,i,&ncols,&row,&vals);
113:       MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);
114:       MatRestoreRow(X,i,&ncols,&row,&vals);
115:     }
116:   } else {
117:     PetscMalloc1(n+1,&val);
118:     for (i=start; i<end; i++) {
119:       MatGetRow(Y,i,&ncols,&row,&vals);
120:       MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);
121:       MatRestoreRow(Y,i,&ncols,&row,&vals);

123:       MatGetRow(X,i,&ncols,&row,&vals);
124:       for (j=0; j<ncols; j++) {
125:         val[j] = a*vals[j];
126:       }
127:       MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);
128:       MatRestoreRow(X,i,&ncols,&row,&vals);
129:     }
130:     PetscFree(val);
131:   }
132:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
133:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
134:   return(0);
135: }

139: /*@
140:    MatShift - Computes Y =  Y + a I, where a is a PetscScalar and I is the identity matrix.

142:    Neighbor-wise Collective on Mat

144:    Input Parameters:
145: +  Y - the matrices
146: -  a - the PetscScalar

148:    Level: intermediate

150:    Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
151:    fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
152:    entry).

154:    Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is
155:     preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix.

157: .keywords: matrix, add, shift

159: .seealso: MatDiagonalSet()
160:  @*/
161: PetscErrorCode  MatShift(Mat Y,PetscScalar a)
162: {

167:   if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
168:   if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
169:   MatCheckPreallocated(Y,1);

171:   (*Y->ops->shift)(Y,a);

173: #if defined(PETSC_HAVE_CUSP)
174:   if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
175:     Y->valid_GPU_matrix = PETSC_CUSP_CPU;
176:   }
177: #endif
178: #if defined(PETSC_HAVE_VIENNACL)
179:   if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) {
180:     Y->valid_GPU_matrix = PETSC_VIENNACL_CPU;
181:   }
182: #endif
183:   return(0);
184: }

188: PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
189: {
191:   PetscInt       i,start,end;
192:   PetscScalar    *v;

195:   MatGetOwnershipRange(Y,&start,&end);
196:   VecGetArray(D,&v);
197:   for (i=start; i<end; i++) {
198:     MatSetValues(Y,1,&i,1,&i,v+i-start,is);
199:   }
200:   VecRestoreArray(D,&v);
201:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
202:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
203:   return(0);
204: }

208: /*@
209:    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
210:    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
211:    INSERT_VALUES.

213:    Input Parameters:
214: +  Y - the input matrix
215: .  D - the diagonal matrix, represented as a vector
216: -  i - INSERT_VALUES or ADD_VALUES

218:    Neighbor-wise Collective on Mat and Vec

220:    Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially
221:    fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an
222:    entry).

224:    Level: intermediate

226: .keywords: matrix, add, shift, diagonal

228: .seealso: MatShift()
229: @*/
230: PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
231: {
233:   PetscInt       matlocal,veclocal;

238:   MatGetLocalSize(Y,&matlocal,NULL);
239:   VecGetLocalSize(D,&veclocal);
240:   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);
241:   if (Y->ops->diagonalset) {
242:     (*Y->ops->diagonalset)(Y,D,is);
243:   } else {
244:     MatDiagonalSet_Default(Y,D,is);
245:   }
246:   return(0);
247: }

251: /*@
252:    MatAYPX - Computes Y = a*Y + X.

254:    Logically on Mat

256:    Input Parameters:
257: +  a - the PetscScalar multiplier
258: .  Y - the first matrix
259: .  X - the second matrix
260: -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN

262:    Level: intermediate

264: .keywords: matrix, add

266: .seealso: MatAXPY()
267:  @*/
268: PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
269: {
270:   PetscScalar    one = 1.0;
272:   PetscInt       mX,mY,nX,nY;

278:   MatGetSize(X,&mX,&nX);
279:   MatGetSize(X,&mY,&nY);
280:   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);

282:   MatScale(Y,a);
283:   MatAXPY(Y,one,X,str);
284:   return(0);
285: }

289: /*@
290:     MatComputeExplicitOperator - Computes the explicit matrix

292:     Collective on Mat

294:     Input Parameter:
295: .   inmat - the matrix

297:     Output Parameter:
298: .   mat - the explict preconditioned operator

300:     Notes:
301:     This computation is done by applying the operators to columns of the
302:     identity matrix.

304:     Currently, this routine uses a dense matrix format when 1 processor
305:     is used and a sparse format otherwise.  This routine is costly in general,
306:     and is recommended for use only with relatively small systems.

308:     Level: advanced

310: .keywords: Mat, compute, explicit, operator
311: @*/
312: PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
313: {
314:   Vec            in,out;
316:   PetscInt       i,m,n,M,N,*rows,start,end;
317:   MPI_Comm       comm;
318:   PetscScalar    *array,zero = 0.0,one = 1.0;
319:   PetscMPIInt    size;


325:   PetscObjectGetComm((PetscObject)inmat,&comm);
326:   MPI_Comm_size(comm,&size);

328:   MatGetLocalSize(inmat,&m,&n);
329:   MatGetSize(inmat,&M,&N);
330:   MatCreateVecs(inmat,&in,&out);
331:   VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
332:   VecGetOwnershipRange(out,&start,&end);
333:   PetscMalloc1(m,&rows);
334:   for (i=0; i<m; i++) rows[i] = start + i;

336:   MatCreate(comm,mat);
337:   MatSetSizes(*mat,m,n,M,N);
338:   if (size == 1) {
339:     MatSetType(*mat,MATSEQDENSE);
340:     MatSeqDenseSetPreallocation(*mat,NULL);
341:   } else {
342:     MatSetType(*mat,MATMPIAIJ);
343:     MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);
344:   }

346:   for (i=0; i<N; i++) {

348:     VecSet(in,zero);
349:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
350:     VecAssemblyBegin(in);
351:     VecAssemblyEnd(in);

353:     MatMult(inmat,in,out);

355:     VecGetArray(out,&array);
356:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
357:     VecRestoreArray(out,&array);

359:   }
360:   PetscFree(rows);
361:   VecDestroy(&out);
362:   VecDestroy(&in);
363:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
364:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
365:   return(0);
366: }

370: /*@
371:   MatChop - Set all values in the matrix less than the tolerance to zero

373:   Input Parameters:
374: + A   - The matrix
375: - tol - The zero tolerance

377:   Output Parameters:
378: . A - The chopped matrix

380:   Level: intermediate

382: .seealso: MatCreate(), MatZeroEntries()
383:  @*/
384: PetscErrorCode MatChop(Mat A, PetscReal tol)
385: {
386:   PetscScalar    *newVals;
387:   PetscInt       *newCols;
388:   PetscInt       rStart, rEnd, numRows, maxRows, r, colMax = 0;

392:   MatGetOwnershipRange(A, &rStart, &rEnd);
393:   for (r = rStart; r < rEnd; ++r) {
394:     PetscInt ncols;

396:     MatGetRow(A, r, &ncols, NULL, NULL);
397:     colMax = PetscMax(colMax, ncols);
398:     MatRestoreRow(A, r, &ncols, NULL, NULL);
399:   }
400:   numRows = rEnd - rStart;
401:   MPI_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));
402:   PetscMalloc2(colMax,&newCols,colMax,&newVals);
403:   for (r = rStart; r < rStart+maxRows; ++r) {
404:     const PetscScalar *vals;
405:     const PetscInt    *cols;
406:     PetscInt           ncols, newcols, c;

408:     if (r < rEnd) {
409:       MatGetRow(A, r, &ncols, &cols, &vals);
410:       for (c = 0; c < ncols; ++c) {
411:         newCols[c] = cols[c];
412:         newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c];
413:       }
414:       newcols = ncols;
415:       MatRestoreRow(A, r, &ncols, &cols, &vals);
416:       MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);
417:     }
418:     MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
419:     MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
420:   }
421:   PetscFree2(newCols,newVals);
422:   return(0);
423: }