Actual source code: matnull.c


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

  6: #include <petsc/private/matimpl.h>

  8: PetscClassId MAT_NULLSPACE_CLASSID;

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

 14:    Logically Collective on MatNullSpace

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

 21:    Level: advanced

 23: .seealso: MatNullSpaceDestroy(), MatNullSpaceRemove(), MatSetNullSpace(), MatNullSpace, MatNullSpaceCreate()
 24: @*/
 25: PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace sp, PetscErrorCode (*rem)(MatNullSpace,Vec,void*),void *ctx)
 26: {
 29:   sp->remove = rem;
 30:   sp->rmctx  = ctx;
 31:   return(0);
 32: }

 34: /*@C
 35:    MatNullSpaceGetVecs - get vectors defining the null space

 37:    Not Collective

 39:    Input Arguments:
 40: .  sp - null space object

 42:    Output Arguments:
 43: +  has_cnst - PETSC_TRUE if the null space contains the constant vector, otherwise PETSC_FALSE
 44: .  n - number of vectors (excluding constant vector) in null space
 45: -  vecs - orthonormal vectors that span the null space (excluding the constant vector)

 47:    Level: developer

 49:    Notes:
 50:       These vectors and the array are owned by the MatNullSpace and should not be destroyed or freeded by the caller

 52: .seealso: MatNullSpaceCreate(), MatGetNullSpace(), MatGetNearNullSpace()
 53: @*/
 54: PetscErrorCode MatNullSpaceGetVecs(MatNullSpace sp,PetscBool *has_const,PetscInt *n,const Vec **vecs)
 55: {

 59:   if (has_const) *has_const = sp->has_cnst;
 60:   if (n) *n = sp->n;
 61:   if (vecs) *vecs = sp->vecs;
 62:   return(0);
 63: }

 65: /*@
 66:    MatNullSpaceCreateRigidBody - create rigid body modes from coordinates

 68:    Collective on Vec

 70:    Input Argument:
 71: .  coords - block of coordinates of each node, must have block size set

 73:    Output Argument:
 74: .  sp - the null space

 76:    Level: advanced

 78:    Notes:
 79:      If you are solving an elasticity problem you should likely use this, in conjunction with MatSetNearNullspace(), to provide information that
 80:      the PCGAMG preconditioner can use to construct a much more efficient preconditioner.

 82:      If you are solving an elasticity problem with pure Neumann boundary conditions you can use this in conjunction with MatSetNullspace() to
 83:      provide this information to the linear solver so it can handle the null space appropriately in the linear solution.


 86: .seealso: MatNullSpaceCreate(), MatSetNearNullspace(), MatSetNullspace()
 87: @*/
 88: PetscErrorCode MatNullSpaceCreateRigidBody(Vec coords,MatNullSpace *sp)
 89: {
 90:   PetscErrorCode    ierr;
 91:   const PetscScalar *x;
 92:   PetscScalar       *v[6],dots[5];
 93:   Vec               vec[6];
 94:   PetscInt          n,N,dim,nmodes,i,j;
 95:   PetscReal         sN;

 98:   VecGetBlockSize(coords,&dim);
 99:   VecGetLocalSize(coords,&n);
100:   VecGetSize(coords,&N);
101:   n   /= dim;
102:   N   /= dim;
103:   sN = 1./PetscSqrtReal((PetscReal)N);
104:   switch (dim) {
105:   case 1:
106:     MatNullSpaceCreate(PetscObjectComm((PetscObject)coords),PETSC_TRUE,0,NULL,sp);
107:     break;
108:   case 2:
109:   case 3:
110:     nmodes = (dim == 2) ? 3 : 6;
111:     VecCreate(PetscObjectComm((PetscObject)coords),&vec[0]);
112:     VecSetSizes(vec[0],dim*n,dim*N);
113:     VecSetBlockSize(vec[0],dim);
114:     VecSetUp(vec[0]);
115:     for (i=1; i<nmodes; i++) {VecDuplicate(vec[0],&vec[i]);}
116:     for (i=0; i<nmodes; i++) {VecGetArray(vec[i],&v[i]);}
117:     VecGetArrayRead(coords,&x);
118:     for (i=0; i<n; i++) {
119:       if (dim == 2) {
120:         v[0][i*2+0] = sN;
121:         v[0][i*2+1] = 0.;
122:         v[1][i*2+0] = 0.;
123:         v[1][i*2+1] = sN;
124:         /* Rotations */
125:         v[2][i*2+0] = -x[i*2+1];
126:         v[2][i*2+1] = x[i*2+0];
127:       } else {
128:         v[0][i*3+0] = sN;
129:         v[0][i*3+1] = 0.;
130:         v[0][i*3+2] = 0.;
131:         v[1][i*3+0] = 0.;
132:         v[1][i*3+1] = sN;
133:         v[1][i*3+2] = 0.;
134:         v[2][i*3+0] = 0.;
135:         v[2][i*3+1] = 0.;
136:         v[2][i*3+2] = sN;

138:         v[3][i*3+0] = x[i*3+1];
139:         v[3][i*3+1] = -x[i*3+0];
140:         v[3][i*3+2] = 0.;
141:         v[4][i*3+0] = 0.;
142:         v[4][i*3+1] = -x[i*3+2];
143:         v[4][i*3+2] = x[i*3+1];
144:         v[5][i*3+0] = x[i*3+2];
145:         v[5][i*3+1] = 0.;
146:         v[5][i*3+2] = -x[i*3+0];
147:       }
148:     }
149:     for (i=0; i<nmodes; i++) {VecRestoreArray(vec[i],&v[i]);}
150:     VecRestoreArrayRead(coords,&x);
151:     for (i=dim; i<nmodes; i++) {
152:       /* Orthonormalize vec[i] against vec[0:i-1] */
153:       VecMDot(vec[i],i,vec,dots);
154:       for (j=0; j<i; j++) dots[j] *= -1.;
155:       VecMAXPY(vec[i],i,dots,vec);
156:       VecNormalize(vec[i],NULL);
157:     }
158:     MatNullSpaceCreate(PetscObjectComm((PetscObject)coords),PETSC_FALSE,nmodes,vec,sp);
159:     for (i=0; i<nmodes; i++) {VecDestroy(&vec[i]);}
160:   }
161:   return(0);
162: }

164: /*@C
165:    MatNullSpaceView - Visualizes a null space object.

167:    Collective on MatNullSpace

169:    Input Parameters:
170: +  matnull - the null space
171: -  viewer - visualization context

173:    Level: advanced

175:    Fortran Note:
176:    This routine is not supported in Fortran.

178: .seealso: MatNullSpaceCreate(), PetscViewerASCIIOpen()
179: @*/
180: PetscErrorCode MatNullSpaceView(MatNullSpace sp,PetscViewer viewer)
181: {
183:   PetscBool      iascii;

187:   if (!viewer) {
188:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sp),&viewer);
189:   }

193:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
194:   if (iascii) {
195:     PetscViewerFormat format;
196:     PetscInt          i;
197:     PetscViewerGetFormat(viewer,&format);
198:     PetscObjectPrintClassNamePrefixType((PetscObject)sp,viewer);
199:     PetscViewerASCIIPushTab(viewer);
200:     PetscViewerASCIIPrintf(viewer,"Contains %D vector%s%s\n",sp->n,sp->n==1 ? "" : "s",sp->has_cnst ? " and the constant" : "");
201:     if (sp->remove) {PetscViewerASCIIPrintf(viewer,"Has user-provided removal function\n");}
202:     if (!(format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) {
203:       for (i=0; i<sp->n; i++) {
204:         VecView(sp->vecs[i],viewer);
205:       }
206:     }
207:     PetscViewerASCIIPopTab(viewer);
208:   }
209:   return(0);
210: }

212: /*@C
213:    MatNullSpaceCreate - Creates a data structure used to project vectors
214:    out of null spaces.

216:    Collective

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

227:    Output Parameter:
228: .  SP - the null space context

230:    Level: advanced

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

235:       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
236:        need to pass in a function that eliminates the constant function into MatNullSpaceSetFunction().

238:   Users manual sections:
239: .   sec_singular

241: .seealso: MatNullSpaceDestroy(), MatNullSpaceRemove(), MatSetNullSpace(), MatNullSpace, MatNullSpaceSetFunction()
242: @*/
243: PetscErrorCode  MatNullSpaceCreate(MPI_Comm comm,PetscBool has_cnst,PetscInt n,const Vec vecs[],MatNullSpace *SP)
244: {
245:   MatNullSpace   sp;
247:   PetscInt       i;

250:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",n);
254:   if (n) {
255:     for (i=0; i<n; i++) {
256:       /* prevent the user from changes values in the vector */
257:       VecLockReadPush(vecs[i]);
258:     }
259:   }
260:   if (PetscUnlikelyDebug(n)) {
261:     PetscScalar *dots;
262:     for (i=0; i<n; i++) {
263:       PetscReal norm;
264:       VecNorm(vecs[i],NORM_2,&norm);
265:       if (PetscAbsReal(norm - 1) > PETSC_SQRT_MACHINE_EPSILON) SETERRQ2(PetscObjectComm((PetscObject)vecs[i]),PETSC_ERR_ARG_WRONG,"Vector %D must have 2-norm of 1.0, it is %g",i,(double)norm);
266:     }
267:     if (has_cnst) {
268:       for (i=0; i<n; i++) {
269:         PetscScalar sum;
270:         VecSum(vecs[i],&sum);
271:         if (PetscAbsScalar(sum) > PETSC_SQRT_MACHINE_EPSILON) SETERRQ2(PetscObjectComm((PetscObject)vecs[i]),PETSC_ERR_ARG_WRONG,"Vector %D must be orthogonal to constant vector, inner product is %g",i,(double)PetscAbsScalar(sum));
272:       }
273:     }
274:     PetscMalloc1(n-1,&dots);
275:     for (i=0; i<n-1; i++) {
276:       PetscInt j;
277:       VecMDot(vecs[i],n-i-1,vecs+i+1,dots);
278:       for (j=0;j<n-i-1;j++) {
279:         if (PetscAbsScalar(dots[j]) > PETSC_SQRT_MACHINE_EPSILON) SETERRQ3(PetscObjectComm((PetscObject)vecs[i]),PETSC_ERR_ARG_WRONG,"Vector %D must be orthogonal to vector %D, inner product is %g",i,i+j+1,(double)PetscAbsScalar(dots[j]));
280:       }
281:     }
282:     PetscFree(dots);
283:   }

285:   *SP = NULL;
286:   MatInitializePackage();

288:   PetscHeaderCreate(sp,MAT_NULLSPACE_CLASSID,"MatNullSpace","Null space","Mat",comm,MatNullSpaceDestroy,MatNullSpaceView);

290:   sp->has_cnst = has_cnst;
291:   sp->n        = n;
292:   sp->vecs     = NULL;
293:   sp->alpha    = NULL;
294:   sp->remove   = NULL;
295:   sp->rmctx    = NULL;

297:   if (n) {
298:     PetscMalloc1(n,&sp->vecs);
299:     PetscMalloc1(n,&sp->alpha);
300:     PetscLogObjectMemory((PetscObject)sp,n*(sizeof(Vec)+sizeof(PetscScalar)));
301:     for (i=0; i<n; i++) {
302:       PetscObjectReference((PetscObject)vecs[i]);
303:       sp->vecs[i] = vecs[i];
304:     }
305:   }

307:   *SP = sp;
308:   return(0);
309: }

311: /*@
312:    MatNullSpaceDestroy - Destroys a data structure used to project vectors
313:    out of null spaces.

315:    Collective on MatNullSpace

317:    Input Parameter:
318: .  sp - the null space context to be destroyed

320:    Level: advanced

322: .seealso: MatNullSpaceCreate(), MatNullSpaceRemove(), MatNullSpaceSetFunction()
323: @*/
324: PetscErrorCode  MatNullSpaceDestroy(MatNullSpace *sp)
325: {
327:   PetscInt       i;

330:   if (!*sp) return(0);
332:   if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; return(0);}

334:   for (i=0; i < (*sp)->n; i++) {
335:     VecLockReadPop((*sp)->vecs[i]);
336:   }

338:   VecDestroyVecs((*sp)->n,&(*sp)->vecs);
339:   PetscFree((*sp)->alpha);
340:   PetscHeaderDestroy(sp);
341:   return(0);
342: }

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

347:    Collective on MatNullSpace

349:    Input Parameters:
350: +  sp - the null space context (if this is NULL then no null space is removed)
351: -  vec - the vector from which the null space is to be removed

353:    Level: advanced

355: .seealso: MatNullSpaceCreate(), MatNullSpaceDestroy(), MatNullSpaceSetFunction()
356: @*/
357: PetscErrorCode  MatNullSpaceRemove(MatNullSpace sp,Vec vec)
358: {
359:   PetscScalar    sum;
360:   PetscInt       i,N;

364:   if (!sp) return(0);

368:   if (sp->has_cnst) {
369:     VecGetSize(vec,&N);
370:     if (N > 0) {
371:       VecSum(vec,&sum);
372:       sum  = sum/((PetscScalar)(-1.0*N));
373:       VecShift(vec,sum);
374:     }
375:   }

377:   if (sp->n) {
378:     VecMDot(vec,sp->n,sp->vecs,sp->alpha);
379:     for (i=0; i<sp->n; i++) sp->alpha[i] = -sp->alpha[i];
380:     VecMAXPY(vec,sp->n,sp->alpha,sp->vecs);
381:   }

383:   if (sp->remove) {
384:     (*sp->remove)(sp,vec,sp->rmctx);
385:   }
386:   return(0);
387: }

389: /*@
390:    MatNullSpaceTest  - Tests if the claimed null space is really a
391:      null space of a matrix

393:    Collective on MatNullSpace

395:    Input Parameters:
396: +  sp - the null space context
397: -  mat - the matrix

399:    Output Parameters:
400: .  isNull - PETSC_TRUE if the nullspace is valid for this matrix

402:    Level: advanced

404: .seealso: MatNullSpaceCreate(), MatNullSpaceDestroy(), MatNullSpaceSetFunction()
405: @*/
406: PetscErrorCode  MatNullSpaceTest(MatNullSpace sp,Mat mat,PetscBool  *isNull)
407: {
408:   PetscScalar    sum;
409:   PetscReal      nrm,tol = 10. * PETSC_SQRT_MACHINE_EPSILON;
410:   PetscInt       j,n,N;
412:   Vec            l,r;
413:   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,consistent = PETSC_TRUE;
414:   PetscViewer    viewer;

419:   n    = sp->n;
420:   PetscOptionsGetBool(((PetscObject)sp)->options,((PetscObject)mat)->prefix,"-mat_null_space_test_view",&flg1,NULL);
421:   PetscOptionsGetBool(((PetscObject)sp)->options,((PetscObject)mat)->prefix,"-mat_null_space_test_view_draw",&flg2,NULL);

423:   if (n) {
424:     VecDuplicate(sp->vecs[0],&l);
425:   } else {
426:     MatCreateVecs(mat,&l,NULL);
427:   }

429:   PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sp),&viewer);
430:   if (sp->has_cnst) {
431:     VecDuplicate(l,&r);
432:     VecGetSize(l,&N);
433:     sum  = 1.0/N;
434:     VecSet(l,sum);
435:     MatMult(mat,l,r);
436:     VecNorm(r,NORM_2,&nrm);
437:     if (nrm >= tol) consistent = PETSC_FALSE;
438:     if (flg1) {
439:       if (consistent) {
440:         PetscPrintf(PetscObjectComm((PetscObject)sp),"Constants are likely null vector");
441:       } else {
442:         PetscPrintf(PetscObjectComm((PetscObject)sp),"Constants are unlikely null vector ");
443:       }
444:       PetscPrintf(PetscObjectComm((PetscObject)sp),"|| A * 1/N || = %g\n",(double)nrm);
445:     }
446:     if (!consistent && flg1) {VecView(r,viewer);}
447:     if (!consistent && flg2) {VecView(r,viewer);}
448:     VecDestroy(&r);
449:   }

451:   for (j=0; j<n; j++) {
452:     (*mat->ops->mult)(mat,sp->vecs[j],l);
453:     VecNorm(l,NORM_2,&nrm);
454:     if (nrm >= tol) consistent = PETSC_FALSE;
455:     if (flg1) {
456:       if (consistent) {
457:         PetscPrintf(PetscObjectComm((PetscObject)sp),"Null vector %D is likely null vector",j);
458:       } else {
459:         PetscPrintf(PetscObjectComm((PetscObject)sp),"Null vector %D unlikely null vector ",j);
460:         consistent = PETSC_FALSE;
461:       }
462:       PetscPrintf(PetscObjectComm((PetscObject)sp),"|| A * v[%D] || = %g\n",j,(double)nrm);
463:     }
464:     if (!consistent && flg1) {VecView(l,viewer);}
465:     if (!consistent && flg2) {VecView(l,viewer);}
466:   }

468:   if (sp->remove) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot test a null space provided as a function with MatNullSpaceSetFunction()");
469:   VecDestroy(&l);
470:   if (isNull) *isNull = consistent;
471:   return(0);
472: }