Actual source code: bddcnullspace.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <bddc.h>
  2: #include <bddcprivate.h>

  6: PetscErrorCode PCBDDCNullSpaceAssembleCoarse(PC pc, Mat coarse_mat, MatNullSpace* CoarseNullSpace)
  7: {
  8:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
  9:   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
 10:   MatNullSpace   tempCoarseNullSpace;
 11:   const Vec      *nsp_vecs;
 12:   Vec            *coarse_nsp_vecs,local_vec,local_primal_vec;
 13:   PetscInt       nsp_size,coarse_nsp_size,i;
 14:   PetscBool      nsp_has_cnst;
 15:   PetscReal      test_null;

 19:   tempCoarseNullSpace = 0;
 20:   coarse_nsp_size = 0;
 21:   coarse_nsp_vecs = 0;
 22:   MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);
 23:   if (coarse_mat) {
 24:     PetscMalloc1((nsp_size+1),&coarse_nsp_vecs);
 25:     for (i=0;i<nsp_size+1;i++) {
 26:       VecDuplicate(pcbddc->coarse_vec,&coarse_nsp_vecs[i]);
 27:     }
 28:   }
 29:   MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);
 30:   if (nsp_has_cnst) {
 31:     VecSet(local_vec,1.0);
 32:     MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);
 33:     PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);
 34:     PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);
 35:     if (coarse_mat) {
 36:       if (pcbddc->dbg_flag) {
 37:         MatMult(coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);
 38:         VecNorm(pcbddc->coarse_rhs,NORM_INFINITY,&test_null);
 39:         PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);
 40:       }
 41:       VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);
 42:       coarse_nsp_size++;
 43:     }
 44:   }
 45:   for (i=0;i<nsp_size;i++)  {
 46:     VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);
 47:     VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);
 48:     MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);
 49:     PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);
 50:     PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);
 51:     if (coarse_mat) {
 52:       if (pcbddc->dbg_flag) {
 53:         MatMult(coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);
 54:         VecNorm(pcbddc->coarse_rhs,NORM_2,&test_null);
 55:         PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);
 56:       }
 57:       VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);
 58:       coarse_nsp_size++;
 59:     }
 60:   }
 61:   if (coarse_nsp_size > 0) {
 62:     PCBDDCOrthonormalizeVecs(coarse_nsp_size,coarse_nsp_vecs);
 63:     MatNullSpaceCreate(PetscObjectComm((PetscObject)coarse_mat),PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);
 64:     for (i=0;i<nsp_size+1;i++) {
 65:       VecDestroy(&coarse_nsp_vecs[i]);
 66:     }
 67:   }
 68:   PetscFree(coarse_nsp_vecs);
 69:   VecDestroy(&local_vec);
 70:   VecDestroy(&local_primal_vec);
 71:   *CoarseNullSpace = tempCoarseNullSpace;
 72:   return(0);
 73: }


 78: static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y)
 79: {
 80:   NullSpaceCorrection_ctx pc_ctx;
 81:   PetscErrorCode          ierr;

 84:   PCShellGetContext(pc,(void**)&pc_ctx);
 85:   /* E */
 86:   MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);
 87:   MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);
 88:   /* P^-1 */
 89:   PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);
 90:   /* E^T */
 91:   MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);
 92:   VecScale(pc_ctx->work_small_1,-1.0);
 93:   MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);
 94:   /* Sum contributions */
 95:   MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);
 96:   return(0);
 97: }

101: static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc)
102: {
103:   NullSpaceCorrection_ctx pc_ctx;
104:   PetscErrorCode          ierr;

107:   PCShellGetContext(pc,(void**)&pc_ctx);
108:   VecDestroy(&pc_ctx->work_small_1);
109:   VecDestroy(&pc_ctx->work_small_2);
110:   VecDestroy(&pc_ctx->work_full_1);
111:   VecDestroy(&pc_ctx->work_full_2);
112:   MatDestroy(&pc_ctx->basis_mat);
113:   MatDestroy(&pc_ctx->Lbasis_mat);
114:   MatDestroy(&pc_ctx->Kbasis_mat);
115:   PCDestroy(&pc_ctx->local_pc);
116:   PetscFree(pc_ctx);
117:   return(0);
118: }

120: /*
121: PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec);
122: PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC);
123: */

127: PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs)
128: {
129:   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
130:   PC_IS          *pcis = (PC_IS*)pc->data;
131:   Mat_IS*        matis = (Mat_IS*)pc->pmat->data;
132:   KSP            *local_ksp;
133:   PC             newpc;
134:   NullSpaceCorrection_ctx  shell_ctx;
135:   Mat            local_mat,local_pmat,small_mat,inv_small_mat;
136:   MatStructure   local_mat_struct;
137:   Vec            work1,work2;
138:   const Vec      *nullvecs;
139:   VecScatter     scatter_ctx;
140:   IS             is_aux;
141:   MatFactorInfo  matinfo;
142:   PetscScalar    *basis_mat,*Kbasis_mat,*array,*array_mat;
143:   PetscScalar    one = 1.0,zero = 0.0, m_one = -1.0;
144:   PetscInt       basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R;
145:   PetscBool      nnsp_has_cnst;

149:   /* Infer the local solver */
150:   ISGetSize(local_dofs,&basis_dofs);
151:   VecGetSize(pcis->vec1_D,&n_I);
152:   VecGetSize(pcbddc->vec1_R,&n_R);
153:   if (basis_dofs == n_I) {
154:     /* Dirichlet solver */
155:     local_ksp = &pcbddc->ksp_D;
156:   } else if (basis_dofs == n_R) {
157:     /* Neumann solver */
158:     local_ksp = &pcbddc->ksp_R;
159:   } else {
160:     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in %s: unknown local IS size %d. n_I=%d, n_R=%d)\n",__FUNCT__,basis_dofs,n_I,n_R);
161:   }
162:   KSPGetOperators(*local_ksp,&local_mat,&local_pmat,&local_mat_struct);

164:   /* Get null space vecs */
165:   MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);
166:   basis_size = nnsp_size;
167:   if (nnsp_has_cnst) {
168:     basis_size++;
169:   }

171:   if (basis_dofs) {
172:      /* Create shell ctx */
173:      PetscMalloc(sizeof(*shell_ctx),&shell_ctx);

175:      /* Create work vectors in shell context */
176:      VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);
177:      VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);
178:      VecSetType(shell_ctx->work_small_1,VECSEQ);
179:      VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);
180:      VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);
181:      VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);
182:      VecSetType(shell_ctx->work_full_1,VECSEQ);
183:      VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);

185:      /* Allocate workspace */
186:      MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );
187:      MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);
188:      MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);
189:      MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);

191:      /* Restrict local null space on selected dofs (Dirichlet or Neumann)
192:         and compute matrices N and K*N */
193:      VecDuplicate(shell_ctx->work_full_1,&work1);
194:      VecDuplicate(shell_ctx->work_full_1,&work2);
195:      VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);
196:   }

198:   for (k=0;k<nnsp_size;k++) {
199:     VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
200:     VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
201:     if (basis_dofs) {
202:       VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);
203:       VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);
204:       VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);
205:       VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);
206:       MatMult(local_mat,work1,work2);
207:       VecResetArray(work1);
208:       VecResetArray(work2);
209:     }
210:   }

212:   if (basis_dofs) {
213:     if (nnsp_has_cnst) {
214:       VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);
215:       VecSet(work1,one);
216:       VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);
217:       MatMult(local_mat,work1,work2);
218:       VecResetArray(work1);
219:       VecResetArray(work2);
220:     }
221:     VecDestroy(&work1);
222:     VecDestroy(&work2);
223:     VecScatterDestroy(&scatter_ctx);
224:     MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);
225:     MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);

227:     /* Assemble another Mat object in shell context */
228:     MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);
229:     MatFactorInfoInitialize(&matinfo);
230:     ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);
231:     MatLUFactor(small_mat,is_aux,is_aux,&matinfo);
232:     ISDestroy(&is_aux);
233:     PetscMalloc1(basis_size*basis_size,&array_mat);
234:     for (k=0;k<basis_size;k++) {
235:       VecSet(shell_ctx->work_small_1,zero);
236:       VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);
237:       VecAssemblyBegin(shell_ctx->work_small_1);
238:       VecAssemblyEnd(shell_ctx->work_small_1);
239:       MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);
240:       VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);
241:       for (i=0;i<basis_size;i++) {
242:         array_mat[i*basis_size+k]=array[i];
243:       }
244:       VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);
245:     }
246:     MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);
247:     MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);
248:     PetscFree(array_mat);
249:     MatDestroy(&inv_small_mat);
250:     MatDestroy(&small_mat);
251:     MatScale(shell_ctx->Kbasis_mat,m_one);

253:     /* Rebuild local PC */
254:     KSPGetPC(*local_ksp,&shell_ctx->local_pc);
255:     PetscObjectReference((PetscObject)shell_ctx->local_pc);
256:     PCCreate(PETSC_COMM_SELF,&newpc);
257:     PCSetOperators(newpc,local_mat,local_mat,SAME_PRECONDITIONER);
258:     PCSetType(newpc,PCSHELL);
259:     PCShellSetContext(newpc,shell_ctx);
260:     PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);
261:     PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);
262:     PCSetUp(newpc);
263:     KSPSetPC(*local_ksp,newpc);
264:     PCDestroy(&newpc);
265:     KSPSetUp(*local_ksp);
266:   }
267:   /* test */
268:   if (pcbddc->dbg_flag && basis_dofs) {
269:     KSP         check_ksp;
270:     PC          check_pc;
271:     Mat         test_mat;
272:     Vec         work3;
273:     PetscReal   test_err,lambda_min,lambda_max;
274:     PetscBool   setsym,issym=PETSC_FALSE;
275:     PetscInt    tabs;

277:     PetscViewerASCIIGetTab(pcbddc->dbg_viewer,&tabs);
278:     KSPGetPC(*local_ksp,&check_pc);
279:     VecDuplicate(shell_ctx->work_full_1,&work1);
280:     VecDuplicate(shell_ctx->work_full_1,&work2);
281:     VecDuplicate(shell_ctx->work_full_1,&work3);
282:     VecSetRandom(shell_ctx->work_small_1,NULL);
283:     MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);
284:     VecCopy(work1,work2);
285:     MatMult(local_mat,work1,work3);
286:     PCApply(check_pc,work3,work1);
287:     VecAXPY(work1,m_one,work2);
288:     VecNorm(work1,NORM_INFINITY,&test_err);
289:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank);
290:     PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_FALSE);
291:     if (basis_dofs == n_I) {
292:       PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Dirichlet ");
293:     } else {
294:       PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Neumann ");
295:     }
296:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"solver is :%1.14e\n",test_err);
297:     PetscViewerASCIISetTab(pcbddc->dbg_viewer,tabs);
298:     PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_TRUE);

300:     MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);
301:     MatShift(test_mat,one);
302:     MatNorm(test_mat,NORM_INFINITY,&test_err);
303:     MatDestroy(&test_mat);
304:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err);

306:     /* Create ksp object suitable for extreme eigenvalues' estimation */
307:     KSPCreate(PETSC_COMM_SELF,&check_ksp);
308:     KSPSetOperators(check_ksp,local_mat,local_mat,SAME_PRECONDITIONER);
309:     KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);
310:     KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);
311:     MatIsSymmetricKnown(pc->pmat,&setsym,&issym);
312:     if (issym) {
313:       KSPSetType(check_ksp,KSPCG);
314:     }
315:     KSPSetPC(check_ksp,check_pc);
316:     KSPSetUp(check_ksp);
317:     VecSetRandom(work1,NULL);
318:     MatMult(local_mat,work1,work2);
319:     KSPSolve(check_ksp,work2,work2);
320:     VecAXPY(work2,m_one,work1);
321:     VecNorm(work2,NORM_INFINITY,&test_err);
322:     KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);
323:     KSPGetIterationNumber(check_ksp,&k);
324:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for adapted KSP %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max);
325:     KSPDestroy(&check_ksp);
326:     VecDestroy(&work1);
327:     VecDestroy(&work2);
328:     VecDestroy(&work3);
329:   }
330:   /* all processes shoud call this, even the void ones */
331:   if (pcbddc->dbg_flag) {
332:     PetscViewerFlush(pcbddc->dbg_viewer);
333:   }
334:   return(0);
335: }

339: PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc)
340: {
341:   PC_IS*         pcis = (PC_IS*)  (pc->data);
342:   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
343:   KSP            inv_change;
344:   PC             pc_change;
345:   const Vec      *nsp_vecs;
346:   Vec            *new_nsp_vecs;
347:   PetscInt       i,nsp_size,new_nsp_size,start_new;
348:   PetscBool      nsp_has_cnst;
349:   MatNullSpace   new_nsp;
351:   MPI_Comm       comm;

354:   /* create KSP for change of basis */
355:   KSPCreate(PETSC_COMM_SELF,&inv_change);
356:   KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix,SAME_PRECONDITIONER);
357:   KSPSetType(inv_change,KSPPREONLY);
358:   KSPGetPC(inv_change,&pc_change);
359:   PCSetType(pc_change,PCLU);
360:   KSPSetUp(inv_change);
361:   /* get nullspace and transform it */
362:   MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);
363:   new_nsp_size = nsp_size;
364:   if (nsp_has_cnst) {
365:     new_nsp_size++;
366:   }
367:   PetscMalloc1(new_nsp_size,&new_nsp_vecs);
368:   for (i=0;i<new_nsp_size;i++) {
369:     VecDuplicate(pcis->vec1_global,&new_nsp_vecs[i]);
370:   }
371:   start_new = 0;
372:   if (nsp_has_cnst) {
373:     start_new = 1;
374:     VecSet(new_nsp_vecs[0],1.0);
375:     VecSet(pcis->vec1_B,1.0);
376:     KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B);
377:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);
378:     VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);
379:   }
380:   for (i=0;i<nsp_size;i++) {
381:     VecCopy(nsp_vecs[i],new_nsp_vecs[i+start_new]);
382:     VecScatterBegin(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
383:     VecScatterEnd  (pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
384:     KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B);
385:     VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);
386:     VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);
387:   }
388:   PCBDDCOrthonormalizeVecs(new_nsp_size,new_nsp_vecs);
389: #if 0
390:   PetscBool nsp_t=PETSC_FALSE;
391:   MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);
392:   printf("Original Null Space test: %d\n",nsp_t);
393:   Mat temp_mat;
394:   Mat_IS* matis = (Mat_IS*)pc->pmat->data;
395:     temp_mat = matis->A;
396:     matis->A = pcbddc->local_mat;
397:     pcbddc->local_mat = temp_mat;
398:   MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);
399:   printf("Original Null Space, mat changed test: %d\n",nsp_t);
400:   {
401:     PetscReal test_norm;
402:     for (i=0;i<new_nsp_size;i++) {
403:       MatMult(pc->pmat,new_nsp_vecs[i],pcis->vec1_global);
404:       VecNorm(pcis->vec1_global,NORM_2,&test_norm);
405:       if (test_norm > 1.e-12) {
406:         printf("------------ERROR VEC %d------------------\n",i);
407:         VecView(pcis->vec1_global,PETSC_VIEWER_STDOUT_WORLD);
408:         printf("------------------------------------------\n");
409:       }
410:     }
411:   }
412: #endif

414:   KSPDestroy(&inv_change);
415:   PetscObjectGetComm((PetscObject)pc,&comm);
416:   MatNullSpaceCreate(comm,PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);
417:   PCBDDCSetNullSpace(pc,new_nsp);
418:   MatNullSpaceDestroy(&new_nsp);
419: #if 0
420:   MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);
421:   printf("New Null Space, mat changed: %d\n",nsp_t);
422:     temp_mat = matis->A;
423:     matis->A = pcbddc->local_mat;
424:     pcbddc->local_mat = temp_mat;
425:   MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);
426:   printf("New Null Space, mat original: %d\n",nsp_t);
427: #endif

429:   for (i=0;i<new_nsp_size;i++) {
430:     VecDestroy(&new_nsp_vecs[i]);
431:   }
432:   PetscFree(new_nsp_vecs);
433:   return(0);
434: }