1 #include "bddc.h" 2 #include "bddcprivate.h" 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "PCBDDCNullSpaceAssembleCoarse" 6 PetscErrorCode PCBDDCNullSpaceAssembleCoarse(PC pc, 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; 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 tempCoarseNullSpace = 0; 20 coarse_nsp_size = 0; 21 coarse_nsp_vecs = 0; 22 ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr); 23 if (pcbddc->coarse_mat) { 24 ierr = PetscMalloc((nsp_size+1)*sizeof(Vec),&coarse_nsp_vecs);CHKERRQ(ierr); 25 for (i=0;i<nsp_size+1;i++) { 26 ierr = VecDuplicate(pcbddc->coarse_vec,&coarse_nsp_vecs[i]);CHKERRQ(ierr); 27 } 28 } 29 ierr = MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);CHKERRQ(ierr); 30 if (nsp_has_cnst) { 31 ierr = VecSet(local_vec,1.0);CHKERRQ(ierr); 32 ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 33 ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34 ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35 if (pcbddc->coarse_mat) { 36 if (pcbddc->dbg_flag) { 37 ierr = MatMult(pcbddc->coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr); 38 ierr = VecNorm(pcbddc->coarse_rhs,NORM_INFINITY,&test_null);CHKERRQ(ierr); 39 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);CHKERRQ(ierr); 40 } 41 ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr); 42 coarse_nsp_size++; 43 } 44 } 45 for (i=0;i<nsp_size;i++) { 46 ierr = VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47 ierr = VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48 ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 49 ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 50 ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 51 if (pcbddc->coarse_mat) { 52 if (pcbddc->dbg_flag) { 53 ierr = MatMult(pcbddc->coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr); 54 ierr = VecNorm(pcbddc->coarse_rhs,NORM_2,&test_null);CHKERRQ(ierr); 55 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);CHKERRQ(ierr); 56 } 57 ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr); 58 coarse_nsp_size++; 59 } 60 } 61 if (coarse_nsp_size > 0) { 62 ierr = PCBDDCOrthonormalizeVecs(coarse_nsp_size,coarse_nsp_vecs);CHKERRQ(ierr); 63 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)(pcbddc->coarse_mat)),PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);CHKERRQ(ierr); 64 for (i=0;i<nsp_size+1;i++) { 65 ierr = VecDestroy(&coarse_nsp_vecs[i]);CHKERRQ(ierr); 66 } 67 } 68 ierr = PetscFree(coarse_nsp_vecs);CHKERRQ(ierr); 69 ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 70 ierr = VecDestroy(&local_primal_vec);CHKERRQ(ierr); 71 *CoarseNullSpace = tempCoarseNullSpace; 72 PetscFunctionReturn(0); 73 } 74 75 76 #undef __FUNCT__ 77 #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC" 78 static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 79 { 80 NullSpaceCorrection_ctx pc_ctx; 81 PetscErrorCode ierr; 82 83 PetscFunctionBegin; 84 ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 85 /* E */ 86 ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 87 ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 88 /* P^-1 */ 89 ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 90 /* E^T */ 91 ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 92 ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 93 ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 94 /* Sum contributions */ 95 ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC" 101 static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 102 { 103 NullSpaceCorrection_ctx pc_ctx; 104 PetscErrorCode ierr; 105 106 PetscFunctionBegin; 107 ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 108 ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 109 ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 110 ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 111 ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 112 ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 113 ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 114 ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 115 ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 116 ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 /*PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec); 121 PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC);*/ 122 123 #undef __FUNCT__ 124 #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection" 125 PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs) 126 { 127 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 128 PC_IS *pcis = (PC_IS*)pc->data; 129 Mat_IS* matis = (Mat_IS*)pc->pmat->data; 130 KSP *local_ksp; 131 PC newpc; 132 NullSpaceCorrection_ctx shell_ctx; 133 Mat local_mat,local_pmat,small_mat,inv_small_mat; 134 MatStructure local_mat_struct; 135 Vec work1,work2; 136 const Vec *nullvecs; 137 VecScatter scatter_ctx; 138 IS is_aux; 139 MatFactorInfo matinfo; 140 PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 141 PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 142 PetscInt basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R; 143 PetscBool nnsp_has_cnst; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 /* Infer the local solver */ 148 ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 149 ierr = VecGetSize(pcis->vec1_D,&n_I);CHKERRQ(ierr); 150 ierr = VecGetSize(pcbddc->vec1_R,&n_R);CHKERRQ(ierr); 151 if (basis_dofs == n_I) { 152 /* Dirichlet solver */ 153 local_ksp = &pcbddc->ksp_D; 154 } else if (basis_dofs == n_R) { 155 /* Neumann solver */ 156 local_ksp = &pcbddc->ksp_R; 157 } else { 158 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); 159 } 160 ierr = KSPGetOperators(*local_ksp,&local_mat,&local_pmat,&local_mat_struct);CHKERRQ(ierr); 161 162 /* Get null space vecs */ 163 ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 164 basis_size = nnsp_size; 165 if (nnsp_has_cnst) { 166 basis_size++; 167 } 168 169 /* Create shell ctx */ 170 ierr = PetscMalloc(sizeof(*shell_ctx),&shell_ctx);CHKERRQ(ierr); 171 172 /* Create work vectors in shell context */ 173 ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 174 ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 175 ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 176 ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 177 ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 178 ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 179 ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 180 ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 181 182 /* Allocate workspace */ 183 ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );CHKERRQ(ierr); 184 ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 185 ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 186 ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 187 188 /* Restrict local null space on selected dofs (Dirichlet or Neumann) 189 and compute matrices N and K*N */ 190 ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 191 ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 192 ierr = VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 193 for (k=0;k<nnsp_size;k++) { 194 ierr = VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 195 ierr = VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196 ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 197 ierr = VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 198 ierr = VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 199 ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 200 ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 201 ierr = VecResetArray(work1);CHKERRQ(ierr); 202 ierr = VecResetArray(work2);CHKERRQ(ierr); 203 } 204 if (nnsp_has_cnst) { 205 ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 206 ierr = VecSet(work1,one);CHKERRQ(ierr); 207 ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 208 ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 209 ierr = VecResetArray(work1);CHKERRQ(ierr); 210 ierr = VecResetArray(work2);CHKERRQ(ierr); 211 } 212 ierr = VecDestroy(&work1);CHKERRQ(ierr); 213 ierr = VecDestroy(&work2);CHKERRQ(ierr); 214 ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 215 ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 216 ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 217 218 /* Assemble another Mat object in shell context */ 219 ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 220 ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 221 ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 222 ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 223 ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 224 ierr = PetscMalloc(basis_size*basis_size*sizeof(PetscScalar),&array_mat);CHKERRQ(ierr); 225 for (k=0;k<basis_size;k++) { 226 ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 227 ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 228 ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 229 ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 230 ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 231 ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 232 for (i=0;i<basis_size;i++) { 233 array_mat[i*basis_size+k]=array[i]; 234 } 235 ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 236 } 237 ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 238 ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 239 ierr = PetscFree(array_mat);CHKERRQ(ierr); 240 ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 241 ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 242 ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 243 244 /* Rebuild local PC */ 245 ierr = KSPGetPC(*local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 246 ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 247 ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 248 ierr = PCSetOperators(newpc,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 249 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 250 ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 251 ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 252 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 253 ierr = PCSetUp(newpc);CHKERRQ(ierr); 254 ierr = KSPSetPC(*local_ksp,newpc);CHKERRQ(ierr); 255 ierr = PCDestroy(&newpc);CHKERRQ(ierr); 256 ierr = KSPSetUp(*local_ksp);CHKERRQ(ierr); 257 /* test */ 258 /* TODO: this cause a deadlock when doing multilevel */ 259 #if 0 260 if (pcbddc->dbg_flag) { 261 KSP check_ksp; 262 PC check_pc; 263 Mat test_mat; 264 Vec work3; 265 PetscViewer viewer=pcbddc->dbg_viewer; 266 PetscReal test_err,lambda_min,lambda_max; 267 PetscBool setsym,issym=PETSC_FALSE; 268 269 ierr = KSPGetPC(*local_ksp,&check_pc);CHKERRQ(ierr); 270 ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 271 ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 272 ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 273 ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 274 ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 275 ierr = VecCopy(work1,work2);CHKERRQ(ierr); 276 ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 277 ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 278 ierr = VecAXPY(work1,m_one,work2);CHKERRQ(ierr); 279 ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 280 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank); 281 if (basis_dofs == n_I) { 282 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Dirichlet "); 283 } else { 284 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Neumann "); 285 } 286 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"solver is :%1.14e\n",test_err); 287 288 ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 289 ierr = MatShift(test_mat,one);CHKERRQ(ierr); 290 ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 291 ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 292 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err); 293 294 /* Create ksp object suitable for extreme eigenvalues' estimation */ 295 ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 296 ierr = KSPSetOperators(check_ksp,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 297 ierr = KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr); 298 ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 299 ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 300 if (issym) { 301 ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 302 } 303 ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 304 ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 305 ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 306 ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 307 ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 308 ierr = VecAXPY(work2,m_one,work1);CHKERRQ(ierr); 309 ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 310 ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 311 ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 312 ierr = PetscViewerASCIISynchronizedPrintf(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); 313 ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 314 ierr = VecDestroy(&work1);CHKERRQ(ierr); 315 ierr = VecDestroy(&work2);CHKERRQ(ierr); 316 ierr = VecDestroy(&work3);CHKERRQ(ierr); 317 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 318 } 319 #endif 320 PetscFunctionReturn(0); 321 } 322 323 #undef __FUNCT__ 324 #define __FUNCT__ "PCBDDCNullSpaceAdaptGlobal" 325 PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc) 326 { 327 PC_IS* pcis = (PC_IS*) (pc->data); 328 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 329 KSP inv_change; 330 PC pc_change; 331 const Vec *nsp_vecs; 332 Vec *new_nsp_vecs; 333 PetscInt i,nsp_size,new_nsp_size,start_new; 334 PetscBool nsp_has_cnst; 335 MatNullSpace new_nsp; 336 PetscErrorCode ierr; 337 MPI_Comm comm; 338 339 PetscFunctionBegin; 340 /* create KSP for change of basis */ 341 ierr = KSPCreate(PETSC_COMM_SELF,&inv_change);CHKERRQ(ierr); 342 ierr = KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix,SAME_PRECONDITIONER);CHKERRQ(ierr); 343 ierr = KSPSetType(inv_change,KSPPREONLY);CHKERRQ(ierr); 344 ierr = KSPGetPC(inv_change,&pc_change);CHKERRQ(ierr); 345 ierr = PCSetType(pc_change,PCLU);CHKERRQ(ierr); 346 ierr = KSPSetUp(inv_change);CHKERRQ(ierr); 347 /* get nullspace and transform it */ 348 ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr); 349 new_nsp_size = nsp_size; 350 if (nsp_has_cnst) { 351 new_nsp_size++; 352 } 353 ierr = PetscMalloc(new_nsp_size*sizeof(Vec),&new_nsp_vecs);CHKERRQ(ierr); 354 for (i=0;i<new_nsp_size;i++) { 355 ierr = VecDuplicate(pcis->vec1_global,&new_nsp_vecs[i]);CHKERRQ(ierr); 356 } 357 start_new = 0; 358 if (nsp_has_cnst) { 359 start_new = 1; 360 ierr = VecSet(new_nsp_vecs[0],1.0);CHKERRQ(ierr); 361 ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr); 362 ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B); 363 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 364 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 365 } 366 for (i=0;i<nsp_size;i++) { 367 ierr = VecCopy(nsp_vecs[i],new_nsp_vecs[i+start_new]);CHKERRQ(ierr); 368 ierr = VecScatterBegin(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 369 ierr = VecScatterEnd (pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 370 ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B); 371 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 372 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 373 } 374 ierr = PCBDDCOrthonormalizeVecs(new_nsp_size,new_nsp_vecs);CHKERRQ(ierr); 375 #if 0 376 PetscBool nsp_t=PETSC_FALSE; 377 ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 378 printf("Original Null Space test: %d\n",nsp_t); 379 Mat temp_mat; 380 Mat_IS* matis = (Mat_IS*)pc->pmat->data; 381 temp_mat = matis->A; 382 matis->A = pcbddc->local_mat; 383 pcbddc->local_mat = temp_mat; 384 ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 385 printf("Original Null Space, mat changed test: %d\n",nsp_t); 386 { 387 PetscReal test_norm; 388 for (i=0;i<new_nsp_size;i++) { 389 ierr = MatMult(pc->pmat,new_nsp_vecs[i],pcis->vec1_global);CHKERRQ(ierr); 390 ierr = VecNorm(pcis->vec1_global,NORM_2,&test_norm);CHKERRQ(ierr); 391 if (test_norm > 1.e-12) { 392 printf("------------ERROR VEC %d------------------\n",i); 393 ierr = VecView(pcis->vec1_global,PETSC_VIEWER_STDOUT_WORLD); 394 printf("------------------------------------------\n"); 395 } 396 } 397 } 398 #endif 399 400 ierr = KSPDestroy(&inv_change);CHKERRQ(ierr); 401 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 402 ierr = MatNullSpaceCreate(comm,PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);CHKERRQ(ierr); 403 ierr = PCBDDCSetNullSpace(pc,new_nsp);CHKERRQ(ierr); 404 ierr = MatNullSpaceDestroy(&new_nsp);CHKERRQ(ierr); 405 #if 0 406 ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 407 printf("New Null Space, mat changed: %d\n",nsp_t); 408 temp_mat = matis->A; 409 matis->A = pcbddc->local_mat; 410 pcbddc->local_mat = temp_mat; 411 ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 412 printf("New Null Space, mat original: %d\n",nsp_t); 413 #endif 414 415 for (i=0;i<new_nsp_size;i++) { 416 ierr = VecDestroy(&new_nsp_vecs[i]);CHKERRQ(ierr); 417 } 418 ierr = PetscFree(new_nsp_vecs);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421