1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <../src/mat/impls/dense/seq/dense.h> 4 5 /* E + small_solve */ 6 static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp,Vec y,Vec x, void* ctx) 7 { 8 NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 9 Mat K; 10 PetscErrorCode ierr; 11 12 PetscFunctionBegin; 13 ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 14 ierr = MatMultTranspose(corr_ctx->basis_mat,y,corr_ctx->sw[0]);CHKERRQ(ierr); 15 ierr = MatMultTranspose(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr); 16 ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr); 17 ierr = MatMult(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0]);CHKERRQ(ierr); 18 ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr); 19 ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr); 20 ierr = MatMultAdd(K,corr_ctx->fw[0],y,y);CHKERRQ(ierr); 21 ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 22 PetscFunctionReturn(0); 23 } 24 25 /* E^t + small */ 26 static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp,Vec y,Vec x, void* ctx) 27 { 28 NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 29 PetscErrorCode ierr; 30 Mat K; 31 32 PetscFunctionBegin; 33 ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 34 ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr); 35 ierr = MatMultTranspose(K,x,corr_ctx->fw[0]);CHKERRQ(ierr); 36 ierr = MatMultTranspose(corr_ctx->basis_mat,corr_ctx->fw[0],corr_ctx->sw[0]);CHKERRQ(ierr); 37 ierr = VecScale(corr_ctx->sw[0],-1.0);CHKERRQ(ierr); 38 ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[2]);CHKERRQ(ierr); 39 ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[2],x,corr_ctx->fw[0]);CHKERRQ(ierr); 40 ierr = VecScale(corr_ctx->fw[0],corr_ctx->scale);CHKERRQ(ierr); 41 /* Sum contributions from approximate solver and projected system */ 42 ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0],x);CHKERRQ(ierr); 43 ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 44 PetscFunctionReturn(0); 45 } 46 47 static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void * ctx) 48 { 49 NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 50 PetscErrorCode ierr; 51 52 PetscFunctionBegin; 53 ierr = VecDestroyVecs(3,&corr_ctx->sw);CHKERRQ(ierr); 54 ierr = VecDestroyVecs(1,&corr_ctx->fw);CHKERRQ(ierr); 55 ierr = MatDestroy(&corr_ctx->basis_mat);CHKERRQ(ierr); 56 ierr = MatDestroy(&corr_ctx->inv_smat);CHKERRQ(ierr); 57 ierr = PetscFree(corr_ctx);CHKERRQ(ierr); 58 PetscFunctionReturn(0); 59 } 60 61 PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 62 { 63 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 64 MatNullSpace NullSpace = NULL; 65 KSP local_ksp; 66 NullSpaceCorrection_ctx shell_ctx; 67 Mat local_mat,local_pmat,dmat,Kbasis_mat; 68 Vec v; 69 PetscContainer c; 70 PetscInt basis_size; 71 IS zerorows; 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 if (isdir) local_ksp = pcbddc->ksp_D; /* Dirichlet solver */ 76 else local_ksp = pcbddc->ksp_R; /* Neumann solver */ 77 ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 78 ierr = MatGetNearNullSpace(local_pmat,&NullSpace);CHKERRQ(ierr); 79 if (!NullSpace) { 80 if (pcbddc->dbg_flag) { 81 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d doesn't have local (near) nullspace: no need for correction in %s solver \n",PetscGlobalRank,isdir ? "Dirichlet" : "Neumann");CHKERRQ(ierr); 82 } 83 PetscFunctionReturn(0); 84 } 85 ierr = PetscObjectQuery((PetscObject)NullSpace,"_PBDDC_Null_dmat",(PetscObject*)&dmat);CHKERRQ(ierr); 86 if (!dmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing dense matrix");CHKERRQ(ierr); 87 ierr = PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr); 88 89 ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 90 shell_ctx->scale = 1.0; 91 ierr = PetscObjectReference((PetscObject)dmat);CHKERRQ(ierr); 92 shell_ctx->basis_mat = dmat; 93 ierr = MatGetSize(dmat,NULL,&basis_size);CHKERRQ(ierr); 94 shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level]; 95 96 /* explicit construct (Phi^T K Phi)^-1 */ 97 ierr = MatMatMult(local_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Kbasis_mat);CHKERRQ(ierr); 98 ierr = MatTransposeMatMult(Kbasis_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->inv_smat);CHKERRQ(ierr); 99 ierr = MatDestroy(&Kbasis_mat);CHKERRQ(ierr); 100 ierr = MatFindZeroRows(shell_ctx->inv_smat,&zerorows);CHKERRQ(ierr); 101 if (zerorows) { /* linearly dependent basis */ 102 const PetscInt *idxs; 103 PetscInt i,nz; 104 105 ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr); 106 ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr); 107 for (i=0;i<nz;i++) { 108 ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],1.0,INSERT_VALUES);CHKERRQ(ierr); 109 } 110 ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr); 111 ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 112 ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 113 } 114 115 ierr = MatLUFactor(shell_ctx->inv_smat,NULL,NULL,NULL);CHKERRQ(ierr); 116 ierr = MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat);CHKERRQ(ierr); 117 if (zerorows) { /* linearly dependent basis */ 118 const PetscInt *idxs; 119 PetscInt i,nz; 120 121 ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr); 122 ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr); 123 for (i=0;i<nz;i++) { 124 ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],0.0,INSERT_VALUES);CHKERRQ(ierr); 125 } 126 ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr); 127 ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128 ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129 } 130 ierr = ISDestroy(&zerorows);CHKERRQ(ierr); 131 132 /* Create work vectors in shell context */ 133 ierr = MatCreateVecs(shell_ctx->inv_smat,&v,NULL);CHKERRQ(ierr); 134 ierr = KSPCreateVecs(local_ksp,1,&shell_ctx->fw,0,NULL);CHKERRQ(ierr); 135 ierr = VecDuplicateVecs(v,3,&shell_ctx->sw);CHKERRQ(ierr); 136 ierr = VecDestroy(&v);CHKERRQ(ierr); 137 138 /* add special pre/post solve to KSP (see [1], eq. 48) */ 139 ierr = KSPSetPreSolve(local_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr); 140 ierr = KSPSetPostSolve(local_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr); 141 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp),&c);CHKERRQ(ierr); 142 ierr = PetscContainerSetPointer(c,shell_ctx);CHKERRQ(ierr); 143 ierr = PetscContainerSetUserDestroy(c,PCBDDCNullSpaceCorrDestroy);CHKERRQ(ierr); 144 ierr = PetscObjectCompose((PetscObject)local_ksp,"_PCBDDC_Null_PrePost_ctx",(PetscObject)c);CHKERRQ(ierr); 145 ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 146 147 /* Create ksp object suitable for extreme eigenvalues' estimation */ 148 if (needscaling || pcbddc->dbg_flag) { 149 KSP check_ksp; 150 PC local_pc; 151 Vec work1,work2; 152 const char* prefix; 153 PetscReal test_err,lambda_min,lambda_max; 154 PetscInt k,maxit; 155 156 ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr); 157 ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr); 158 ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 159 if (local_mat->spd) { 160 ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 161 } 162 ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0);CHKERRQ(ierr); 163 ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 164 ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 165 ierr = KSPAppendOptionsPrefix(check_ksp,"approximate_scale_");CHKERRQ(ierr); 166 ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr); 167 ierr = KSPSetOperators(check_ksp,local_mat,local_pmat);CHKERRQ(ierr); 168 ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 169 ierr = KSPSetPreSolve(check_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr); 170 ierr = KSPSetPostSolve(check_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr); 171 ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 172 ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 173 /* setup with default maxit, then set maxit to min(10,any_set_from_command_line) (bug in computing eigenvalues when chaning the number of iterations */ 174 ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 175 ierr = KSPGetPC(local_ksp,&local_pc);CHKERRQ(ierr); 176 ierr = KSPSetPC(check_ksp,local_pc);CHKERRQ(ierr); 177 ierr = KSPGetTolerances(check_ksp,NULL,NULL,NULL,&maxit);CHKERRQ(ierr); 178 ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PetscMin(10,maxit));CHKERRQ(ierr); 179 ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr); 180 ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr); 181 ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr); 182 ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr); 183 ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 184 ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 185 ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 186 ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 187 if (pcbddc->dbg_flag) { 188 if (isdir) { 189 ierr = 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);CHKERRQ(ierr); 190 } else { 191 ierr = 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);CHKERRQ(ierr); 192 } 193 } 194 if (needscaling) shell_ctx->scale = 1.0/lambda_max; 195 196 if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */ 197 PC new_pc; 198 199 ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr); 200 ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr); 201 ierr = PCCreate(PetscObjectComm((PetscObject)check_ksp),&new_pc);CHKERRQ(ierr); 202 ierr = PCSetType(new_pc,PCKSP);CHKERRQ(ierr); 203 ierr = PCSetOperators(new_pc,local_mat,local_pmat);CHKERRQ(ierr); 204 ierr = PCKSPSetKSP(new_pc,local_ksp);CHKERRQ(ierr); 205 ierr = KSPSetPC(check_ksp,new_pc);CHKERRQ(ierr); 206 ierr = PCDestroy(&new_pc);CHKERRQ(ierr); 207 ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);CHKERRQ(ierr); 208 ierr = KSPSetPreSolve(check_ksp,NULL,NULL);CHKERRQ(ierr); 209 ierr = KSPSetPostSolve(check_ksp,NULL,NULL);CHKERRQ(ierr); 210 ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr); 211 ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr); 212 ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 213 ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 214 ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 215 ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 216 if (pcbddc->dbg_flag) { 217 if (isdir) { 218 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet adapted solver (scale %g) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,(double)PetscRealPart(shell_ctx->scale),test_err,k,lambda_min,lambda_max);CHKERRQ(ierr); 219 } else { 220 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann adapted solver (scale %g) %1.14e (it %D, eigs %1.6e %1.6e)\n",PetscGlobalRank,(double)PetscRealPart(shell_ctx->scale),test_err,k,lambda_min,lambda_max);CHKERRQ(ierr); 221 } 222 } 223 } 224 ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 225 ierr = VecDestroy(&work1);CHKERRQ(ierr); 226 ierr = VecDestroy(&work2);CHKERRQ(ierr); 227 } 228 ierr = PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr); 229 230 if (pcbddc->dbg_flag) { 231 Vec work1,work2,work3; 232 PetscReal test_err; 233 234 /* check nullspace basis is solved exactly */ 235 ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr); 236 ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr); 237 ierr = VecDuplicate(shell_ctx->fw[0],&work3);CHKERRQ(ierr); 238 ierr = VecSetRandom(shell_ctx->sw[0],NULL);CHKERRQ(ierr); 239 ierr = MatMult(shell_ctx->basis_mat,shell_ctx->sw[0],work1);CHKERRQ(ierr); 240 ierr = VecCopy(work1,work2);CHKERRQ(ierr); 241 ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 242 ierr = KSPSolve(local_ksp,work3,work1);CHKERRQ(ierr); 243 ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 244 ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 245 if (isdir) { 246 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 247 } else { 248 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 249 } 250 ierr = VecDestroy(&work1);CHKERRQ(ierr); 251 ierr = VecDestroy(&work2);CHKERRQ(ierr); 252 ierr = VecDestroy(&work3);CHKERRQ(ierr); 253 } 254 PetscFunctionReturn(0); 255 } 256