1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 392cccca0SStefano Zampini #include <../src/mat/impls/dense/seq/dense.h> 4674ae819SStefano Zampini 592cccca0SStefano Zampini /* E + small_solve */ 692cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPreSolve(KSP ksp,Vec y,Vec x, void* ctx) 7bd5e1169SStefano Zampini { 892cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 992cccca0SStefano Zampini Mat K; 10bd5e1169SStefano Zampini PetscErrorCode ierr; 11bd5e1169SStefano Zampini 12bd5e1169SStefano Zampini PetscFunctionBegin; 138ead10e4SStefano Zampini ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 1492cccca0SStefano Zampini ierr = MatMultTranspose(corr_ctx->basis_mat,y,corr_ctx->sw[0]);CHKERRQ(ierr); 157ebab0bbSStefano Zampini if (corr_ctx->symm) { 167ebab0bbSStefano Zampini ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr); 177ebab0bbSStefano Zampini } else { 1892cccca0SStefano Zampini ierr = MatMultTranspose(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[1]);CHKERRQ(ierr); 197ebab0bbSStefano Zampini } 2092cccca0SStefano Zampini ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr); 2192cccca0SStefano Zampini ierr = MatMult(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0]);CHKERRQ(ierr); 2292cccca0SStefano Zampini ierr = VecScale(corr_ctx->sw[1],-1.0);CHKERRQ(ierr); 2392cccca0SStefano Zampini ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr); 2492cccca0SStefano Zampini ierr = MatMultAdd(K,corr_ctx->fw[0],y,y);CHKERRQ(ierr); 258ead10e4SStefano Zampini ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 26bd5e1169SStefano Zampini PetscFunctionReturn(0); 27bd5e1169SStefano Zampini } 28bd5e1169SStefano Zampini 2992cccca0SStefano Zampini /* E^t + small */ 3092cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrPostSolve(KSP ksp,Vec y,Vec x, void* ctx) 31674ae819SStefano Zampini { 3292cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 33674ae819SStefano Zampini PetscErrorCode ierr; 3492cccca0SStefano Zampini Mat K; 35674ae819SStefano Zampini 36674ae819SStefano Zampini PetscFunctionBegin; 378ead10e4SStefano Zampini ierr = PetscLogEventBegin(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 3892cccca0SStefano Zampini ierr = KSPGetOperators(ksp,&K,NULL);CHKERRQ(ierr); 397ebab0bbSStefano Zampini if (corr_ctx->symm) { 407ebab0bbSStefano Zampini ierr = MatMult(K,x,corr_ctx->fw[0]);CHKERRQ(ierr); 417ebab0bbSStefano Zampini } else { 4292cccca0SStefano Zampini ierr = MatMultTranspose(K,x,corr_ctx->fw[0]);CHKERRQ(ierr); 437ebab0bbSStefano Zampini } 4492cccca0SStefano Zampini ierr = MatMultTranspose(corr_ctx->basis_mat,corr_ctx->fw[0],corr_ctx->sw[0]);CHKERRQ(ierr); 4592cccca0SStefano Zampini ierr = VecScale(corr_ctx->sw[0],-1.0);CHKERRQ(ierr); 4692cccca0SStefano Zampini ierr = MatMult(corr_ctx->inv_smat,corr_ctx->sw[0],corr_ctx->sw[2]);CHKERRQ(ierr); 4792cccca0SStefano Zampini ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[2],x,corr_ctx->fw[0]);CHKERRQ(ierr); 4892cccca0SStefano Zampini ierr = VecScale(corr_ctx->fw[0],corr_ctx->scale);CHKERRQ(ierr); 4992cccca0SStefano Zampini /* Sum contributions from approximate solver and projected system */ 5092cccca0SStefano Zampini ierr = MatMultAdd(corr_ctx->basis_mat,corr_ctx->sw[1],corr_ctx->fw[0],x);CHKERRQ(ierr); 518ead10e4SStefano Zampini ierr = PetscLogEventEnd(corr_ctx->evapply,ksp,0,0,0);CHKERRQ(ierr); 52674ae819SStefano Zampini PetscFunctionReturn(0); 53674ae819SStefano Zampini } 54674ae819SStefano Zampini 5592cccca0SStefano Zampini static PetscErrorCode PCBDDCNullSpaceCorrDestroy(void * ctx) 56674ae819SStefano Zampini { 5792cccca0SStefano Zampini NullSpaceCorrection_ctx corr_ctx = (NullSpaceCorrection_ctx)ctx; 58674ae819SStefano Zampini PetscErrorCode ierr; 59674ae819SStefano Zampini 60674ae819SStefano Zampini PetscFunctionBegin; 6192cccca0SStefano Zampini ierr = VecDestroyVecs(3,&corr_ctx->sw);CHKERRQ(ierr); 6292cccca0SStefano Zampini ierr = VecDestroyVecs(1,&corr_ctx->fw);CHKERRQ(ierr); 6392cccca0SStefano Zampini ierr = MatDestroy(&corr_ctx->basis_mat);CHKERRQ(ierr); 6492cccca0SStefano Zampini ierr = MatDestroy(&corr_ctx->inv_smat);CHKERRQ(ierr); 6592cccca0SStefano Zampini ierr = PetscFree(corr_ctx);CHKERRQ(ierr); 66674ae819SStefano Zampini PetscFunctionReturn(0); 67674ae819SStefano Zampini } 68674ae819SStefano Zampini 69c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 70674ae819SStefano Zampini { 71674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 72c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 73185763b3SStefano Zampini KSP local_ksp; 74674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 7592cccca0SStefano Zampini Mat local_mat,local_pmat,dmat,Kbasis_mat; 7692cccca0SStefano Zampini Vec v; 7792cccca0SStefano Zampini PetscContainer c; 7892cccca0SStefano Zampini PetscInt basis_size; 7992cccca0SStefano Zampini IS zerorows; 80*80fdaca0SStefano Zampini PetscBool iscusp; 81674ae819SStefano Zampini PetscErrorCode ierr; 82674ae819SStefano Zampini 83674ae819SStefano Zampini PetscFunctionBegin; 848ead10e4SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; /* Dirichlet solver */ 858ead10e4SStefano Zampini else local_ksp = pcbddc->ksp_R; /* Neumann solver */ 8692cccca0SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 8792cccca0SStefano Zampini ierr = MatGetNearNullSpace(local_pmat,&NullSpace);CHKERRQ(ierr); 882958b453SStefano Zampini if (!NullSpace) { 8935529e7bSStefano Zampini if (pcbddc->dbg_flag) { 902958b453SStefano Zampini 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); 9135529e7bSStefano Zampini } 9235529e7bSStefano Zampini PetscFunctionReturn(0); 9335529e7bSStefano Zampini } 9492cccca0SStefano Zampini ierr = PetscObjectQuery((PetscObject)NullSpace,"_PBDDC_Null_dmat",(PetscObject*)&dmat);CHKERRQ(ierr); 9592cccca0SStefano Zampini if (!dmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing dense matrix");CHKERRQ(ierr); 968ead10e4SStefano Zampini ierr = PetscLogEventBegin(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr); 972958b453SStefano Zampini 98854ce69bSBarry Smith ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 99bd5e1169SStefano Zampini shell_ctx->scale = 1.0; 10092cccca0SStefano Zampini ierr = PetscObjectReference((PetscObject)dmat);CHKERRQ(ierr); 10192cccca0SStefano Zampini shell_ctx->basis_mat = dmat; 10292cccca0SStefano Zampini ierr = MatGetSize(dmat,NULL,&basis_size);CHKERRQ(ierr); 1038ead10e4SStefano Zampini shell_ctx->evapply = PC_BDDC_ApproxApply[pcbddc->current_level]; 10492cccca0SStefano Zampini 1057ebab0bbSStefano Zampini ierr = MatGetOption(local_mat,MAT_SYMMETRIC,&shell_ctx->symm);CHKERRQ(ierr); 1067ebab0bbSStefano Zampini 10792cccca0SStefano Zampini /* explicit construct (Phi^T K Phi)^-1 */ 108*80fdaca0SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)local_mat,MATSEQAIJCUSPARSE,&iscusp);CHKERRQ(ierr); 109*80fdaca0SStefano Zampini if (iscusp) { 110*80fdaca0SStefano Zampini ierr = MatConvert(shell_ctx->basis_mat,MATSEQDENSECUDA,MAT_INPLACE_MATRIX,&shell_ctx->basis_mat);CHKERRQ(ierr); 111*80fdaca0SStefano Zampini } 11292cccca0SStefano Zampini ierr = MatMatMult(local_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Kbasis_mat);CHKERRQ(ierr); 11392cccca0SStefano Zampini ierr = MatTransposeMatMult(Kbasis_mat,shell_ctx->basis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->inv_smat);CHKERRQ(ierr); 11492cccca0SStefano Zampini ierr = MatDestroy(&Kbasis_mat);CHKERRQ(ierr); 115*80fdaca0SStefano Zampini ierr = MatBindToCPU(shell_ctx->inv_smat,PETSC_TRUE);CHKERRQ(ierr); 11692cccca0SStefano Zampini ierr = MatFindZeroRows(shell_ctx->inv_smat,&zerorows);CHKERRQ(ierr); 11792cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 11892cccca0SStefano Zampini const PetscInt *idxs; 11992cccca0SStefano Zampini PetscInt i,nz; 12092cccca0SStefano Zampini 12192cccca0SStefano Zampini ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr); 12292cccca0SStefano Zampini ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr); 12392cccca0SStefano Zampini for (i=0;i<nz;i++) { 12492cccca0SStefano Zampini ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],1.0,INSERT_VALUES);CHKERRQ(ierr); 12592cccca0SStefano Zampini } 12692cccca0SStefano Zampini ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr); 12792cccca0SStefano Zampini ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12892cccca0SStefano Zampini ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12992cccca0SStefano Zampini } 13092cccca0SStefano Zampini ierr = MatLUFactor(shell_ctx->inv_smat,NULL,NULL,NULL);CHKERRQ(ierr); 13192cccca0SStefano Zampini ierr = MatSeqDenseInvertFactors_Private(shell_ctx->inv_smat);CHKERRQ(ierr); 13292cccca0SStefano Zampini if (zerorows) { /* linearly dependent basis */ 13392cccca0SStefano Zampini const PetscInt *idxs; 13492cccca0SStefano Zampini PetscInt i,nz; 13592cccca0SStefano Zampini 13692cccca0SStefano Zampini ierr = ISGetLocalSize(zerorows,&nz);CHKERRQ(ierr); 13792cccca0SStefano Zampini ierr = ISGetIndices(zerorows,&idxs);CHKERRQ(ierr); 13892cccca0SStefano Zampini for (i=0;i<nz;i++) { 13992cccca0SStefano Zampini ierr = MatSetValue(shell_ctx->inv_smat,idxs[i],idxs[i],0.0,INSERT_VALUES);CHKERRQ(ierr); 14092cccca0SStefano Zampini } 14192cccca0SStefano Zampini ierr = ISRestoreIndices(zerorows,&idxs);CHKERRQ(ierr); 14292cccca0SStefano Zampini ierr = MatAssemblyBegin(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14392cccca0SStefano Zampini ierr = MatAssemblyEnd(shell_ctx->inv_smat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14492cccca0SStefano Zampini } 14592cccca0SStefano Zampini ierr = ISDestroy(&zerorows);CHKERRQ(ierr); 146674ae819SStefano Zampini 147674ae819SStefano Zampini /* Create work vectors in shell context */ 14892cccca0SStefano Zampini ierr = MatCreateVecs(shell_ctx->inv_smat,&v,NULL);CHKERRQ(ierr); 14992cccca0SStefano Zampini ierr = KSPCreateVecs(local_ksp,1,&shell_ctx->fw,0,NULL);CHKERRQ(ierr); 15092cccca0SStefano Zampini ierr = VecDuplicateVecs(v,3,&shell_ctx->sw);CHKERRQ(ierr); 15192cccca0SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 152674ae819SStefano Zampini 15392cccca0SStefano Zampini /* add special pre/post solve to KSP (see [1], eq. 48) */ 15492cccca0SStefano Zampini ierr = KSPSetPreSolve(local_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr); 15592cccca0SStefano Zampini ierr = KSPSetPostSolve(local_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr); 15692cccca0SStefano Zampini ierr = PetscContainerCreate(PetscObjectComm((PetscObject)local_ksp),&c);CHKERRQ(ierr); 15792cccca0SStefano Zampini ierr = PetscContainerSetPointer(c,shell_ctx);CHKERRQ(ierr); 15892cccca0SStefano Zampini ierr = PetscContainerSetUserDestroy(c,PCBDDCNullSpaceCorrDestroy);CHKERRQ(ierr); 15992cccca0SStefano Zampini ierr = PetscObjectCompose((PetscObject)local_ksp,"_PCBDDC_Null_PrePost_ctx",(PetscObject)c);CHKERRQ(ierr); 16092cccca0SStefano Zampini ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 161c7017625SStefano Zampini 162c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 16392cccca0SStefano Zampini if (needscaling || pcbddc->dbg_flag) { 164674ae819SStefano Zampini KSP check_ksp; 16592cccca0SStefano Zampini PC local_pc; 16692cccca0SStefano Zampini Vec work1,work2; 167bd5e1169SStefano Zampini const char* prefix; 16892cccca0SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 16992cccca0SStefano Zampini PetscInt k,maxit; 170c7017625SStefano Zampini 17192cccca0SStefano Zampini ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr); 17292cccca0SStefano Zampini ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr); 173c7017625SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 17492cccca0SStefano Zampini if (local_mat->spd) { 17592cccca0SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 17692cccca0SStefano Zampini } 17715579a77SStefano Zampini ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0);CHKERRQ(ierr); 178bd5e1169SStefano Zampini ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 179bd5e1169SStefano Zampini ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 18092cccca0SStefano Zampini ierr = KSPAppendOptionsPrefix(check_ksp,"approximate_scale_");CHKERRQ(ierr); 181bd5e1169SStefano Zampini ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr); 182c1a3ebd0SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_pmat);CHKERRQ(ierr); 183c7017625SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 18492cccca0SStefano Zampini ierr = KSPSetPreSolve(check_ksp,PCBDDCNullSpaceCorrPreSolve,shell_ctx);CHKERRQ(ierr); 18592cccca0SStefano Zampini ierr = KSPSetPostSolve(check_ksp,PCBDDCNullSpaceCorrPostSolve,shell_ctx);CHKERRQ(ierr); 18692cccca0SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 187bd5e1169SStefano Zampini ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 18892cccca0SStefano Zampini /* 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 */ 189c7017625SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 19092cccca0SStefano Zampini ierr = KSPGetPC(local_ksp,&local_pc);CHKERRQ(ierr); 19192cccca0SStefano Zampini ierr = KSPSetPC(check_ksp,local_pc);CHKERRQ(ierr); 19292cccca0SStefano Zampini ierr = KSPGetTolerances(check_ksp,NULL,NULL,NULL,&maxit);CHKERRQ(ierr); 19392cccca0SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PetscMin(10,maxit));CHKERRQ(ierr); 19492cccca0SStefano Zampini ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr); 19592cccca0SStefano Zampini ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr); 19692cccca0SStefano Zampini ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr); 19792cccca0SStefano Zampini ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr); 19892cccca0SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 19992cccca0SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 200c7017625SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 201c7017625SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 202c7017625SStefano Zampini if (pcbddc->dbg_flag) { 203c7017625SStefano Zampini if (isdir) { 20492cccca0SStefano Zampini 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); 205c7017625SStefano Zampini } else { 20692cccca0SStefano Zampini 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); 207c7017625SStefano Zampini } 208c7017625SStefano Zampini } 20992cccca0SStefano Zampini if (needscaling) shell_ctx->scale = 1.0/lambda_max; 21092cccca0SStefano Zampini 21192cccca0SStefano Zampini if (needscaling && pcbddc->dbg_flag) { /* test for scaling factor */ 21292cccca0SStefano Zampini PC new_pc; 21392cccca0SStefano Zampini 21492cccca0SStefano Zampini ierr = VecSetRandom(work2,NULL);CHKERRQ(ierr); 21592cccca0SStefano Zampini ierr = MatMult(local_mat,work2,work1);CHKERRQ(ierr); 21692cccca0SStefano Zampini ierr = PCCreate(PetscObjectComm((PetscObject)check_ksp),&new_pc);CHKERRQ(ierr); 21792cccca0SStefano Zampini ierr = PCSetType(new_pc,PCKSP);CHKERRQ(ierr); 21892cccca0SStefano Zampini ierr = PCSetOperators(new_pc,local_mat,local_pmat);CHKERRQ(ierr); 21992cccca0SStefano Zampini ierr = PCKSPSetKSP(new_pc,local_ksp);CHKERRQ(ierr); 22092cccca0SStefano Zampini ierr = KSPSetPC(check_ksp,new_pc);CHKERRQ(ierr); 22192cccca0SStefano Zampini ierr = PCDestroy(&new_pc);CHKERRQ(ierr); 22292cccca0SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);CHKERRQ(ierr); 22392cccca0SStefano Zampini ierr = KSPSetPreSolve(check_ksp,NULL,NULL);CHKERRQ(ierr); 22492cccca0SStefano Zampini ierr = KSPSetPostSolve(check_ksp,NULL,NULL);CHKERRQ(ierr); 22592cccca0SStefano Zampini ierr = KSPSolve(check_ksp,work1,work1);CHKERRQ(ierr); 22692cccca0SStefano Zampini ierr = KSPCheckSolve(check_ksp,pc,work1);CHKERRQ(ierr); 22792cccca0SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 22892cccca0SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 22992cccca0SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 23092cccca0SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 23192cccca0SStefano Zampini if (pcbddc->dbg_flag) { 23292cccca0SStefano Zampini if (isdir) { 23392cccca0SStefano Zampini 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); 23492cccca0SStefano Zampini } else { 23592cccca0SStefano Zampini 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); 23692cccca0SStefano Zampini } 23792cccca0SStefano Zampini } 23892cccca0SStefano Zampini } 239c7017625SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 24092cccca0SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 24192cccca0SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 242c7017625SStefano Zampini } 2438ead10e4SStefano Zampini ierr = PetscLogEventEnd(PC_BDDC_ApproxSetUp[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr); 244c7017625SStefano Zampini 24592cccca0SStefano Zampini if (pcbddc->dbg_flag) { 246c7017625SStefano Zampini Vec work1,work2,work3; 24792cccca0SStefano Zampini PetscReal test_err; 248674ae819SStefano Zampini 24992cccca0SStefano Zampini /* check nullspace basis is solved exactly */ 25092cccca0SStefano Zampini ierr = VecDuplicate(shell_ctx->fw[0],&work1);CHKERRQ(ierr); 25192cccca0SStefano Zampini ierr = VecDuplicate(shell_ctx->fw[0],&work2);CHKERRQ(ierr); 25292cccca0SStefano Zampini ierr = VecDuplicate(shell_ctx->fw[0],&work3);CHKERRQ(ierr); 25392cccca0SStefano Zampini ierr = VecSetRandom(shell_ctx->sw[0],NULL);CHKERRQ(ierr); 25492cccca0SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->sw[0],work1);CHKERRQ(ierr); 255674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 256674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 25792cccca0SStefano Zampini ierr = KSPSolve(local_ksp,work3,work1);CHKERRQ(ierr); 258c7017625SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 259674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 260185763b3SStefano Zampini if (isdir) { 261c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 262674ae819SStefano Zampini } else { 263c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 264674ae819SStefano Zampini } 265674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 266674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 267674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 26892cccca0SStefano Zampini } 269674ae819SStefano Zampini PetscFunctionReturn(0); 270674ae819SStefano Zampini } 271