Actual source code: matnull.c

  2: /*
  3:     Routines to project vectors out of null spaces.
  4: */

  6: #include <private/matimpl.h>      /*I "petscmat.h" I*/

  8: PetscClassId  MAT_NULLSPACE_CLASSID;

 12: /*@C
 13:    MatNullSpaceSetFunction - set a function that removes a null space from a vector
 14:    out of null spaces.

 16:    Logically Collective on MatNullSpace

 18:    Input Parameters:
 19: +  sp - the null space object
 20: .  rem - the function that removes the null space
 21: -  ctx - context for the remove function

 23:    Level: advanced

 25: .keywords: PC, null space, create

 27: .seealso: MatNullSpaceDestroy(), MatNullSpaceRemove(), KSPSetNullSpace(), MatNullSpace, MatNullSpaceCreate()
 28: @*/
 29: PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace sp, PetscErrorCode (*rem)(MatNullSpace,Vec,void*),void *ctx)
 30: {
 33:   sp->remove = rem;
 34:   sp->rmctx  = ctx;
 35:   return(0);
 36: }

 40: /*@C
 41:    MatNullSpaceView - Visualizes a null space object.

 43:    Collective on MatNullSpace

 45:    Input Parameters:
 46: +  matnull - the null space
 47: -  viewer - visualization context

 49:    Level: advanced

 51:    Fortran Note:
 52:    This routine is not supported in Fortran.

 54: .seealso: MatNullSpaceCreate(), PetscViewerASCIIOpen()
 55: @*/
 56: PetscErrorCode MatNullSpaceView(MatNullSpace sp,PetscViewer viewer)
 57: {
 59:   PetscBool      iascii;

 63:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)sp)->comm);

 67:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 68:   if (iascii) {
 69:     PetscViewerFormat format;
 70:     PetscInt i;
 71:     PetscViewerGetFormat(viewer,&format);
 72:     PetscObjectPrintClassNamePrefixType((PetscObject)sp,viewer,"MatNullSpace Object");
 73:     PetscViewerASCIIPushTab(viewer);
 74:     PetscViewerASCIIPrintf(viewer,"Contains %D vector%s%s\n",sp->n,sp->n==1?"":"s",sp->has_cnst?" and the constant":"");
 75:     if (sp->remove) {PetscViewerASCIIPrintf(viewer,"Has user-provided removal function\n");}
 76:     if (!(format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) {
 77:       for (i=0; i<sp->n; i++) {
 78:         VecView(sp->vecs[i],viewer);
 79:       }
 80:     }
 81:     PetscViewerASCIIPopTab(viewer);
 82:   }
 83:   return(0);
 84: }

 88: /*@
 89:    MatNullSpaceCreate - Creates a data structure used to project vectors 
 90:    out of null spaces.

 92:    Collective on MPI_Comm

 94:    Input Parameters:
 95: +  comm - the MPI communicator associated with the object
 96: .  has_cnst - PETSC_TRUE if the null space contains the constant vector; otherwise PETSC_FALSE
 97: .  n - number of vectors (excluding constant vector) in null space
 98: -  vecs - the vectors that span the null space (excluding the constant vector);
 99:           these vectors must be orthonormal. These vectors are NOT copied, so do not change them
100:           after this call. You should free the array that you pass in and destroy the vectors (this will reduce the reference count 
101:           for them by one).

103:    Output Parameter:
104: .  SP - the null space context

106:    Level: advanced

108:    Notes: See MatNullSpaceSetFunction() as an alternative way of providing the null space information instead of setting vecs.

110:       If has_cnst is PETSC_TRUE you do not need to pass a constant vector in as a fourth argument to this routine, nor do you 
111:        need to pass in a function that eliminates the constant function into MatNullSpaceSetFunction().

113:   Users manual sections:
114: .   Section 4.16 Solving Singular Systems

116: .keywords: PC, null space, create

118: .seealso: MatNullSpaceDestroy(), MatNullSpaceRemove(), KSPSetNullSpace(), MatNullSpace, MatNullSpaceSetFunction()
119: @*/
120: PetscErrorCode  MatNullSpaceCreate(MPI_Comm comm,PetscBool  has_cnst,PetscInt n,const Vec vecs[],MatNullSpace *SP)
121: {
122:   MatNullSpace   sp;
124:   PetscInt       i;

127:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",n);
131: 
132:   *SP = PETSC_NULL;
133: #ifndef PETSC_USE_DYNAMIC_LIBRARIES 
134:   MatInitializePackage(PETSC_NULL);
135: #endif 

137:   PetscHeaderCreate(sp,_p_MatNullSpace,int,MAT_NULLSPACE_CLASSID,0,"MatNullSpace","Null space","Mat",comm,MatNullSpaceDestroy,MatNullSpaceView);

139:   sp->has_cnst = has_cnst;
140:   sp->n        = n;
141:   sp->vecs     = 0;
142:   sp->alpha    = 0;
143:   sp->vec      = 0;
144:   sp->remove   = 0;
145:   sp->rmctx    = 0;

147:   if (n) {
148:     PetscMalloc(n*sizeof(Vec),&sp->vecs);
149:     PetscMalloc(n*sizeof(PetscScalar),&sp->alpha);
150:     PetscLogObjectMemory(sp,n*(sizeof(Vec)+sizeof(PetscScalar)));
151:     for (i=0; i<n; i++) {
152:       PetscObjectReference((PetscObject)vecs[i]);
153:       sp->vecs[i] = vecs[i];
154:     }
155:   }

157:   *SP          = sp;
158:   return(0);
159: }

163: /*@
164:    MatNullSpaceDestroy - Destroys a data structure used to project vectors 
165:    out of null spaces.

167:    Collective on MatNullSpace

169:    Input Parameter:
170: .  sp - the null space context to be destroyed

172:    Level: advanced

174: .keywords: PC, null space, destroy

176: .seealso: MatNullSpaceCreate(), MatNullSpaceRemove(), MatNullSpaceSetFunction()
177: @*/
178: PetscErrorCode  MatNullSpaceDestroy(MatNullSpace *sp)
179: {

183:   if (!*sp) return(0);
185:   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; return(0);}

187:   VecDestroy(&(*sp)->vec);
188:   VecDestroyVecs((*sp)->n,&(*sp)->vecs);
189:   PetscFree((*sp)->alpha);
190:   PetscHeaderDestroy(sp);
191:   return(0);
192: }

196: /*@C
197:    MatNullSpaceRemove - Removes all the components of a null space from a vector.

199:    Collective on MatNullSpace

201:    Input Parameters:
202: +  sp - the null space context
203: .  vec - the vector from which the null space is to be removed 
204: -  out - if this is requested (not PETSC_NULL) then this is a vector with the null space removed otherwise
205:          the removal is done in-place (in vec)

207:    Note: The user is not responsible for the vector returned and should not destroy it.

209:    Level: advanced

211: .keywords: PC, null space, remove

213: .seealso: MatNullSpaceCreate(), MatNullSpaceDestroy(), MatNullSpaceSetFunction()
214: @*/
215: PetscErrorCode  MatNullSpaceRemove(MatNullSpace sp,Vec vec,Vec *out)
216: {
217:   PetscScalar    sum;
218:   PetscInt       i,N;


225:   if (out) {
227:     if (!sp->vec) {
228:       VecDuplicate(vec,&sp->vec);
229:       PetscLogObjectParent(sp,sp->vec);
230:     }
231:     VecCopy(vec,sp->vec);
232:     vec = *out = sp->vec;
233:   }
234: 
235:   if (sp->has_cnst) {
236:     VecGetSize(vec,&N);
237:     if (N > 0) {
238:       VecSum(vec,&sum);
239:       sum  = sum/((PetscScalar)(-1.0*N));
240:       VecShift(vec,sum);
241:     }
242:   }
243: 
244:   if (sp->n) {
245:     VecMDot(vec,sp->n,sp->vecs,sp->alpha);
246:     for (i=0; i<sp->n; i++) sp->alpha[i] = -sp->alpha[i];
247:     VecMAXPY(vec,sp->n,sp->alpha,sp->vecs);
248:   }

250:   if (sp->remove){
251:     (*sp->remove)(sp,vec,sp->rmctx);
252:   }
253:   return(0);
254: }

258: /*@
259:    MatNullSpaceTest  - Tests if the claimed null space is really a
260:      null space of a matrix

262:    Collective on MatNullSpace

264:    Input Parameters:
265: +  sp - the null space context
266: -  mat - the matrix

268:    Output Parameters:
269: .  isNull - PETSC_TRUE if the nullspace is valid for this matrix

271:    Level: advanced

273: .keywords: PC, null space, remove

275: .seealso: MatNullSpaceCreate(), MatNullSpaceDestroy(), MatNullSpaceSetFunction()
276: @*/
277: PetscErrorCode  MatNullSpaceTest(MatNullSpace sp,Mat mat,PetscBool  *isNull)
278: {
279:   PetscScalar    sum;
280:   PetscReal      nrm;
281:   PetscInt       j,n,N;
283:   Vec            l,r;
284:   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,consistent = PETSC_TRUE;
285:   PetscViewer    viewer;

290:   n = sp->n;
291:   PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test_view",&flg1,PETSC_NULL);
292:   PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test_view_draw",&flg2,PETSC_NULL);

294:   if (!sp->vec) {
295:     if (n) {
296:       VecDuplicate(sp->vecs[0],&sp->vec);
297:     } else {
298:       MatGetVecs(mat,&sp->vec,PETSC_NULL);
299:     }
300:   }
301:   l    = sp->vec;

303:   PetscViewerASCIIGetStdout(((PetscObject)sp)->comm,&viewer);
304:   if (sp->has_cnst) {
305:     VecDuplicate(l,&r);
306:     VecGetSize(l,&N);
307:     sum  = 1.0/N;
308:     VecSet(l,sum);
309:     MatMult(mat,l,r);
310:     VecNorm(r,NORM_2,&nrm);
311:     if (nrm < 1.e-7) {
312:       PetscPrintf(((PetscObject)sp)->comm,"Constants are likely null vector");
313:     } else {
314:       PetscPrintf(((PetscObject)sp)->comm,"Constants are unlikely null vector ");
315:       consistent = PETSC_FALSE;
316:     }
317:     PetscPrintf(((PetscObject)sp)->comm,"|| A * 1/N || = %G\n",nrm);
318:     if (nrm > 1.e-7 && flg1) {VecView(r,viewer);}
319:     if (nrm > 1.e-7 && flg2) {VecView(r,viewer);}
320:     VecDestroy(&r);
321:   }

323:   for (j=0; j<n; j++) {
324:     (*mat->ops->mult)(mat,sp->vecs[j],l);
325:     VecNorm(l,NORM_2,&nrm);
326:     if (nrm < 1.e-7) {
327:       PetscPrintf(((PetscObject)sp)->comm,"Null vector %D is likely null vector",j);
328:     } else {
329:       PetscPrintf(((PetscObject)sp)->comm,"Null vector %D unlikely null vector ",j);
330:       consistent = PETSC_FALSE;
331:     }
332:     PetscPrintf(((PetscObject)sp)->comm,"|| A * v[%D] || = %G\n",j,nrm);
333:     if (nrm > 1.e-7 && flg1) {VecView(l,viewer);}
334:     if (nrm > 1.e-7 && flg2) {VecView(l,viewer);}
335:   }

337:   if (sp->remove) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot test a null space provided as a function with MatNullSpaceSetFunction()");
338:   if (isNull) *isNull = consistent;
339:   return(0);
340: }