Actual source code: axpy.c

  1: #define PETSCMAT_DLL

 3:  #include private/matimpl.h

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

 10:    Collective on Mat

 12:    Input Parameters:
 13: +  a - the scalar multiplier
 14: .  X - the first matrix
 15: .  Y - the second matrix
 16: -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 
 17:          or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's)

 19:    Notes:
 20:      Will only be efficient if one has the SAME_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN

 22:    Level: intermediate

 24: .keywords: matrix, add

 26: .seealso: MatAYPX()
 27:  @*/
 28: PetscErrorCode  MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str)
 29: {
 31:   PetscInt       m1,m2,n1,n2;


 37:   MatGetSize(X,&m1,&n1);
 38:   MatGetSize(Y,&m2,&n2);
 39:   if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %D %D %D %D",m1,m2,n1,n2);

 41:   if (Y->ops->axpy) {
 42:     (*Y->ops->axpy)(Y,a,X,str);
 43:   } else {
 44:     MatAXPY_Basic(Y,a,X,str);
 45:   }
 46:   return(0);
 47: }


 52: PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str)
 53: {
 54:   PetscInt          i,start,end,j,ncols,m,n;
 55:   PetscErrorCode    ierr;
 56:   const PetscInt    *row;
 57:   PetscScalar       *val;
 58:   const PetscScalar *vals;

 61:   MatGetSize(X,&m,&n);
 62:   MatGetOwnershipRange(X,&start,&end);
 63:   if (a == 1.0) {
 64:     for (i = start; i < end; i++) {
 65:       MatGetRow(X,i,&ncols,&row,&vals);
 66:       MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);
 67:       MatRestoreRow(X,i,&ncols,&row,&vals);
 68:     }
 69:   } else {
 70:     PetscMalloc((n+1)*sizeof(PetscScalar),&val);
 71:     for (i=start; i<end; i++) {
 72:       MatGetRow(X,i,&ncols,&row,&vals);
 73:       for (j=0; j<ncols; j++) {
 74:         val[j] = a*vals[j];
 75:       }
 76:       MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);
 77:       MatRestoreRow(X,i,&ncols,&row,&vals);
 78:     }
 79:     PetscFree(val);
 80:   }
 81:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
 83:   return(0);
 84: }

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

 91:    Collective on Mat

 93:    Input Parameters:
 94: +  Y - the matrices
 95: -  a - the PetscScalar 

 97:    Level: intermediate

 99: .keywords: matrix, add, shift

101: .seealso: MatDiagonalSet()
102:  @*/
103: PetscErrorCode  MatShift(Mat Y,PetscScalar a)
104: {
106:   PetscInt       i,start,end;

110:   if (!Y->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
111:   if (Y->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
112:   MatPreallocated(Y);

114:   if (Y->ops->shift) {
115:     (*Y->ops->shift)(Y,a);
116:   } else {
117:     PetscScalar alpha = a;
118:     MatGetOwnershipRange(Y,&start,&end);
119:     for (i=start; i<end; i++) {
120:       MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);
121:     }
122:     MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
123:     MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
124:   }
125:   return(0);
126: }

130: PetscErrorCode  MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is)
131: {
133:   PetscInt       i,start,end,vstart,vend;
134:   PetscScalar    *v;

137:   VecGetOwnershipRange(D,&vstart,&vend);
138:   MatGetOwnershipRange(Y,&start,&end);
139:   if (vstart != start || vend != end) {
140:     SETERRQ4(PETSC_ERR_ARG_SIZ,"Vector ownership range not compatible with matrix: %D %D vec %D %D mat",vstart,vend,start,end);
141:   }
142:   VecGetArray(D,&v);
143:   for (i=start; i<end; i++) {
144:     MatSetValues(Y,1,&i,1,&i,v+i-start,is);
145:   }
146:   VecRestoreArray(D,&v);
147:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
148:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
149:   return(0);
150: }

154: /*@
155:    MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix
156:    that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is
157:    INSERT_VALUES.

159:    Input Parameters:
160: +  Y - the input matrix
161: .  D - the diagonal matrix, represented as a vector
162: -  i - INSERT_VALUES or ADD_VALUES

164:    Collective on Mat and Vec

166:    Level: intermediate

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

170: .seealso: MatShift()
171: @*/
172: PetscErrorCode  MatDiagonalSet(Mat Y,Vec D,InsertMode is)
173: {

179:   if (Y->ops->diagonalset) {
180:     (*Y->ops->diagonalset)(Y,D,is);
181:   } else {
182:     MatDiagonalSet_Default(Y,D,is);
183:   }
184:   return(0);
185: }

189: /*@
190:    MatAYPX - Computes Y = a*Y + X.

192:    Collective on Mat

194:    Input Parameters:
195: +  a - the PetscScalar multiplier
196: .  Y - the first matrix
197: .  X - the second matrix
198: -  str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 

200:    Notes:
201:      Will only be efficient if one has the SAME_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN

203:    Level: intermediate

205: .keywords: matrix, add

207: .seealso: MatAXPY()
208:  @*/
209: PetscErrorCode  MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str)
210: {
211:   PetscScalar    one = 1.0;
213:   PetscInt       mX,mY,nX,nY;


219:   MatGetSize(X,&mX,&nX);
220:   MatGetSize(X,&mY,&nY);
221:   if (mX != mY || nX != nY) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrices: %D %D first %D %D second",mX,mY,nX,nY);

223:   MatScale(Y,a);
224:   MatAXPY(Y,one,X,str);CHKERRQ(ierr)
225:   return(0);
226: }

230: /*@
231:     MatComputeExplicitOperator - Computes the explicit matrix

233:     Collective on Mat

235:     Input Parameter:
236: .   inmat - the matrix

238:     Output Parameter:
239: .   mat - the explict preconditioned operator

241:     Notes:
242:     This computation is done by applying the operators to columns of the 
243:     identity matrix.

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

249:     Level: advanced
250:    
251: .keywords: Mat, compute, explicit, operator

253: @*/
254: PetscErrorCode  MatComputeExplicitOperator(Mat inmat,Mat *mat)
255: {
256:   Vec            in,out;
258:   PetscInt       i,m,n,M,N,*rows,start,end;
259:   MPI_Comm       comm;
260:   PetscScalar    *array,zero = 0.0,one = 1.0;
261:   PetscMPIInt    size;


267:   comm = ((PetscObject)inmat)->comm;
268:   MPI_Comm_size(comm,&size);

270:   MatGetLocalSize(inmat,&m,&n);
271:   MatGetSize(inmat,&M,&N);
272:   MatGetVecs(inmat,&in,&out);
273:   VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
274:   VecGetOwnershipRange(out,&start,&end);
275:   PetscMalloc(m*sizeof(PetscInt),&rows);
276:   for (i=0; i<m; i++) {rows[i] = start + i;}

278:   MatCreate(comm,mat);
279:   MatSetSizes(*mat,m,n,M,N);
280:   if (size == 1) {
281:     MatSetType(*mat,MATSEQDENSE);
282:     MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
283:   } else {
284:     MatSetType(*mat,MATMPIAIJ);
285:     MatMPIAIJSetPreallocation(*mat,n,PETSC_NULL,N-n,PETSC_NULL);
286:   }

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

290:     VecSet(in,zero);
291:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
292:     VecAssemblyBegin(in);
293:     VecAssemblyEnd(in);

295:     MatMult(inmat,in,out);

297:     VecGetArray(out,&array);
298:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
299:     VecRestoreArray(out,&array);

301:   }
302:   PetscFree(rows);
303:   VecDestroy(out);
304:   VecDestroy(in);
305:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
306:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
307:   return(0);
308: }

310: /* Get the map xtoy which is used by MatAXPY() in the case of SUBSET_NONZERO_PATTERN */
313: PetscErrorCode MatAXPYGetxtoy_Private(PetscInt m,PetscInt *xi,PetscInt *xj,PetscInt *xgarray, PetscInt *yi,PetscInt *yj,PetscInt *ygarray, PetscInt **xtoy)
314: {
316:   PetscInt       row,i,nz,xcol,ycol,jx,jy,*x2y;

319:   PetscMalloc(xi[m]*sizeof(PetscInt),&x2y);
320:   i = 0;
321:   for (row=0; row<m; row++){
322:     nz = xi[1] - xi[0];
323:     jy = 0;
324:     for (jx=0; jx<nz; jx++,jy++){
325:       if (xgarray && ygarray){
326:         xcol = xgarray[xj[*xi + jx]];
327:         ycol = ygarray[yj[*yi + jy]];
328:       } else {
329:         xcol = xj[*xi + jx];
330:         ycol = yj[*yi + jy];  /* col index for y */
331:       }
332:       while ( ycol < xcol ) {
333:         jy++;
334:         if (ygarray){
335:           ycol = ygarray[yj[*yi + jy]];
336:         } else {
337:           ycol = yj[*yi + jy];
338:         }
339:       }
340:       if (xcol != ycol) SETERRQ2(PETSC_ERR_ARG_WRONG,"X matrix entry (%D,%D) is not in Y matrix",row,ycol);
341:       x2y[i++] = *yi + jy;
342:     }
343:     xi++; yi++;
344:   }
345:   *xtoy = x2y;
346:   return(0);
347: }