Actual source code: bddcnullspace.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  1:  #include <../src/ksp/pc/impls/bddc/bddc.h>
  2:  #include <../src/ksp/pc/impls/bddc/bddcprivate.h>

  4: static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y)
  5: {
  6:   NullSpaceCorrection_ctx pc_ctx;
  7:   PetscErrorCode          ierr;

 10:   PCShellGetContext(pc,(void**)&pc_ctx);
 11:   /* E */
 12:   MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);
 13:   MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);
 14:   /* P^-1 */
 15:   PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);
 16:   /* E^T */
 17:   MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);
 18:   VecScale(pc_ctx->work_small_1,-1.0);
 19:   MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);
 20:   /* Sum contributions */
 21:   MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);
 22:   if (pc_ctx->apply_scaling) {
 23:     VecScale(y,pc_ctx->scale);
 24:   }
 25:   return(0);
 26: }

 28: static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc)
 29: {
 30:   NullSpaceCorrection_ctx pc_ctx;
 31:   PetscErrorCode          ierr;

 34:   PCShellGetContext(pc,(void**)&pc_ctx);
 35:   VecDestroy(&pc_ctx->work_small_1);
 36:   VecDestroy(&pc_ctx->work_small_2);
 37:   VecDestroy(&pc_ctx->work_full_1);
 38:   VecDestroy(&pc_ctx->work_full_2);
 39:   MatDestroy(&pc_ctx->basis_mat);
 40:   MatDestroy(&pc_ctx->Lbasis_mat);
 41:   MatDestroy(&pc_ctx->Kbasis_mat);
 42:   PCDestroy(&pc_ctx->local_pc);
 43:   PetscFree(pc_ctx);
 44:   return(0);
 45: }

 47: PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling)
 48: {
 49:   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
 50:   PC_IS                    *pcis = (PC_IS*)pc->data;
 51:   Mat_IS                   *matis = (Mat_IS*)pc->pmat->data;
 52:   MatNullSpace             NullSpace = NULL;
 53:   KSP                      local_ksp;
 54:   PC                       newpc;
 55:   NullSpaceCorrection_ctx  shell_ctx;
 56:   Mat                      local_mat,local_pmat,small_mat,inv_small_mat;
 57:   const Vec                *nullvecs;
 58:   VecScatter               scatter_ctx;
 59:   IS                       is_aux,local_dofs;
 60:   MatFactorInfo            matinfo;
 61:   PetscScalar              *basis_mat,*Kbasis_mat,*array,*array_mat;
 62:   PetscScalar              one = 1.0,zero = 0.0, m_one = -1.0;
 63:   PetscInt                 basis_dofs,basis_size,nnsp_size,i,k;
 64:   PetscBool                nnsp_has_cnst;
 65:   PetscReal                test_err,lambda_min,lambda_max;
 66:   PetscBool                setsym,issym=PETSC_FALSE;
 67:   PetscErrorCode           ierr;

 70:   MatGetNullSpace(matis->A,&NullSpace);
 71:   if (!NullSpace) return(0);
 72:   /* Infer the local solver */
 73:   if (isdir) {
 74:     /* Dirichlet solver */
 75:     local_ksp = pcbddc->ksp_D;
 76:     local_dofs = pcis->is_I_local;
 77:   } else {
 78:     /* Neumann solver */
 79:     local_ksp = pcbddc->ksp_R;
 80:     local_dofs = pcbddc->is_R_local;
 81:   }
 82:   ISGetSize(local_dofs,&basis_dofs);
 83:   KSPGetOperators(local_ksp,&local_mat,&local_pmat);

 85:   /* Get null space vecs */
 86:   MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);
 87:   basis_size = nnsp_size;
 88:   if (nnsp_has_cnst) basis_size++;

 90:    /* Create shell ctx */
 91:   PetscNew(&shell_ctx);
 92:   shell_ctx->apply_scaling = needscaling;

 94:   /* Create work vectors in shell context */
 95:   VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);
 96:   VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);
 97:   VecSetType(shell_ctx->work_small_1,VECSEQ);
 98:   VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);
 99:   VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);
100:   VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);
101:   VecSetType(shell_ctx->work_full_1,VECSEQ);
102:   VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);

104:   /* Allocate workspace */
105:   MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);
106:   MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);
107:   MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);
108:   MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);

110:   /* Restrict local null space on selected dofs
111:      and compute matrices N and K*N */
112:   VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);
113:   for (k=0;k<nnsp_size;k++) {
114:     VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);
115:     VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);
116:     VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);
117:     VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);
118:     MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);
119:     VecResetArray(shell_ctx->work_full_1);
120:     VecResetArray(shell_ctx->work_full_2);
121:   }
122:   if (nnsp_has_cnst) {
123:     VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);
124:     VecSet(shell_ctx->work_full_1,one);
125:     VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);
126:     MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);
127:     VecResetArray(shell_ctx->work_full_1);
128:     VecResetArray(shell_ctx->work_full_2);
129:   }
130:   VecScatterDestroy(&scatter_ctx);
131:   MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);
132:   MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);

134:   /* Assemble another Mat object in shell context */
135:   MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);
136:   MatFactorInfoInitialize(&matinfo);
137:   ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);
138:   MatLUFactor(small_mat,is_aux,is_aux,&matinfo);
139:   ISDestroy(&is_aux);
140:   PetscMalloc1(basis_size*basis_size,&array_mat);
141:   for (k=0;k<basis_size;k++) {
142:     VecSet(shell_ctx->work_small_1,zero);
143:     VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);
144:     VecAssemblyBegin(shell_ctx->work_small_1);
145:     VecAssemblyEnd(shell_ctx->work_small_1);
146:     MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);
147:     VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);
148:     for (i=0;i<basis_size;i++) {
149:       array_mat[i*basis_size+k]=array[i];
150:     }
151:     VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);
152:   }
153:   MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);
154:   MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);
155:   PetscFree(array_mat);
156:   MatDestroy(&inv_small_mat);
157:   MatDestroy(&small_mat);
158:   MatScale(shell_ctx->Kbasis_mat,m_one);

160:   /* Rebuild local PC */
161:   KSPGetPC(local_ksp,&shell_ctx->local_pc);
162:   PetscObjectReference((PetscObject)shell_ctx->local_pc);
163:   PCCreate(PETSC_COMM_SELF,&newpc);
164:   PCSetOperators(newpc,local_mat,local_mat);
165:   PCSetType(newpc,PCSHELL);
166:   PCShellSetContext(newpc,shell_ctx);
167:   PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);
168:   PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);
169:   PCSetUp(newpc);
170:   KSPSetPC(local_ksp,newpc);
171:   KSPSetUp(local_ksp);

173:   /* Create ksp object suitable for extreme eigenvalues' estimation */
174:   if (needscaling) {
175:     KSP check_ksp;
176:     Vec *workv;

178:     KSPCreate(PETSC_COMM_SELF,&check_ksp);
179:     KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);
180:     KSPSetOperators(check_ksp,local_mat,local_mat);
181:     KSPSetTolerances(check_ksp,1.e-4,1.e-12,PETSC_DEFAULT,basis_dofs);
182:     KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);
183:     MatIsSymmetricKnown(pc->pmat,&setsym,&issym);
184:     if (issym) {
185:       KSPSetType(check_ksp,KSPCG);
186:     }
187:     VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);
188:     KSPSetPC(check_ksp,newpc);
189:     KSPSetUp(check_ksp);
190:     VecSetRandom(workv[1],NULL);
191:     MatMult(local_mat,workv[1],workv[0]);
192:     KSPSolve(check_ksp,workv[0],workv[0]);
193:     VecAXPY(workv[0],-1.,workv[1]);
194:     VecNorm(workv[0],NORM_INFINITY,&test_err);
195:     KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);
196:     KSPGetIterationNumber(check_ksp,&k);
197:     if (pcbddc->dbg_flag) {
198:       if (isdir) {
199:         PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted solver (no scale) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max);
200:       } else {
201:         PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted solver (no scale) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max);
202:       }
203:     }
204:     shell_ctx->scale = 1./lambda_max;
205:     KSPDestroy(&check_ksp);
206:     VecDestroyVecs(2,&workv);
207:   }
208:   PCDestroy(&newpc);
209:   return(0);
210: }

212: PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir)
213: {
214:   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
215:   Mat_IS                   *matis = (Mat_IS*)pc->pmat->data;
216:   KSP                      check_ksp,local_ksp;
217:   MatNullSpace             NullSpace = NULL;
218:   NullSpaceCorrection_ctx  shell_ctx;
219:   PC                       check_pc;
220:   Mat                      test_mat,local_mat;
221:   PetscReal                test_err,lambda_min,lambda_max;
222:   PetscBool                setsym,issym=PETSC_FALSE;
223:   Vec                      work1,work2,work3;
224:   PetscInt                 k;
225:   PetscErrorCode           ierr;

228:   MatGetNullSpace(matis->A,&NullSpace);
229:   if (!NullSpace) return(0);
230:   if (!pcbddc->dbg_flag) return(0);
231:   if (isdir) local_ksp = pcbddc->ksp_D;
232:   else local_ksp = pcbddc->ksp_R;
233:   KSPGetOperators(local_ksp,&local_mat,NULL);
234:   KSPGetPC(local_ksp,&check_pc);
235:   PCShellGetContext(check_pc,(void**)&shell_ctx);
236:   VecDuplicate(shell_ctx->work_full_1,&work1);
237:   VecDuplicate(shell_ctx->work_full_1,&work2);
238:   VecDuplicate(shell_ctx->work_full_1,&work3);

240:   VecSetRandom(shell_ctx->work_small_1,NULL);
241:   MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);
242:   VecCopy(work1,work2);
243:   MatMult(local_mat,work1,work3);
244:   PCApply(check_pc,work3,work1);
245:   VecAXPY(work1,-1.,work2);
246:   VecNorm(work1,NORM_INFINITY,&test_err);
247:   if (isdir) {
248:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);
249:   } else {
250:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);
251:   }

253:   MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);
254:   MatShift(test_mat,1.);
255:   MatNorm(test_mat,NORM_INFINITY,&test_err);
256:   MatDestroy(&test_mat);
257:   if (isdir) {
258:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);
259:   } else {
260:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);
261:   }

263:   /* Create ksp object suitable for extreme eigenvalues' estimation */
264:   KSPCreate(PETSC_COMM_SELF,&check_ksp);
265:   KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);
266:   KSPSetOperators(check_ksp,local_mat,local_mat);
267:   KSPSetTolerances(check_ksp,1.e-10,1.e-8,PETSC_DEFAULT,PETSC_DEFAULT);
268:   KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);
269:   MatIsSymmetricKnown(pc->pmat,&setsym,&issym);
270:   if (issym) {
271:     KSPSetType(check_ksp,KSPCG);
272:   }
273:   KSPSetPC(check_ksp,check_pc);
274:   KSPSetUp(check_ksp);
275:   VecSetRandom(work1,NULL);
276:   MatMult(local_mat,work1,work2);
277:   KSPSolve(check_ksp,work2,work2);
278:   VecAXPY(work2,-1.,work1);
279:   VecNorm(work2,NORM_INFINITY,&test_err);
280:   KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);
281:   KSPGetIterationNumber(check_ksp,&k);
282:   if (isdir) {
283:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted KSP (scale %d) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,shell_ctx->apply_scaling,test_err,k,lambda_min,lambda_max);
284:   } else {
285:     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted KSP (scale %d) %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,shell_ctx->apply_scaling,test_err,k,lambda_min,lambda_max);
286:   }
287:   KSPDestroy(&check_ksp);
288:   VecDestroy(&work1);
289:   VecDestroy(&work2);
290:   VecDestroy(&work3);
291:   return(0);
292: }