1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 4674ae819SStefano Zampini #undef __FUNCT__ 5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC" 6674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 7674ae819SStefano Zampini { 8674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 9674ae819SStefano Zampini PetscErrorCode ierr; 10674ae819SStefano Zampini 11674ae819SStefano Zampini PetscFunctionBegin; 12674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 13674ae819SStefano Zampini /* E */ 14674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 15674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 16674ae819SStefano Zampini /* P^-1 */ 17674ae819SStefano Zampini ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 18674ae819SStefano Zampini /* E^T */ 19674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 20674ae819SStefano Zampini ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 21674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 22674ae819SStefano Zampini /* Sum contributions */ 23674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 24*c7017625SStefano Zampini if (pc_ctx->apply_scaling) { 25*c7017625SStefano Zampini ierr = VecScale(y,pc_ctx->scale);CHKERRQ(ierr); 26*c7017625SStefano Zampini } 27674ae819SStefano Zampini PetscFunctionReturn(0); 28674ae819SStefano Zampini } 29674ae819SStefano Zampini 30674ae819SStefano Zampini #undef __FUNCT__ 31674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC" 32674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 33674ae819SStefano Zampini { 34674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 35674ae819SStefano Zampini PetscErrorCode ierr; 36674ae819SStefano Zampini 37674ae819SStefano Zampini PetscFunctionBegin; 38674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 39674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 40674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 41674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 42674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 43674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 44674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 45674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 46674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 48674ae819SStefano Zampini PetscFunctionReturn(0); 49674ae819SStefano Zampini } 50674ae819SStefano Zampini 51674ae819SStefano Zampini #undef __FUNCT__ 52674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection" 53*c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 54674ae819SStefano Zampini { 55674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 56674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 57674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 58*c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 59185763b3SStefano Zampini KSP local_ksp; 60674ae819SStefano Zampini PC newpc; 61674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 62674ae819SStefano Zampini Mat local_mat,local_pmat,small_mat,inv_small_mat; 63674ae819SStefano Zampini const Vec *nullvecs; 64674ae819SStefano Zampini VecScatter scatter_ctx; 65*c7017625SStefano Zampini IS is_aux,local_dofs; 66674ae819SStefano Zampini MatFactorInfo matinfo; 67674ae819SStefano Zampini PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 68674ae819SStefano Zampini PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 69185763b3SStefano Zampini PetscInt basis_dofs,basis_size,nnsp_size,i,k; 70674ae819SStefano Zampini PetscBool nnsp_has_cnst; 71*c7017625SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 72*c7017625SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 73674ae819SStefano Zampini PetscErrorCode ierr; 74674ae819SStefano Zampini 75674ae819SStefano Zampini PetscFunctionBegin; 76*c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 77*c7017625SStefano Zampini if (!NullSpace) PetscFunctionReturn(0); 78674ae819SStefano Zampini /* Infer the local solver */ 79185763b3SStefano Zampini if (isdir) { 80674ae819SStefano Zampini /* Dirichlet solver */ 81185763b3SStefano Zampini local_ksp = pcbddc->ksp_D; 82*c7017625SStefano Zampini local_dofs = pcis->is_I_local; 83674ae819SStefano Zampini } else { 84185763b3SStefano Zampini /* Neumann solver */ 85185763b3SStefano Zampini local_ksp = pcbddc->ksp_R; 86*c7017625SStefano Zampini local_dofs = pcbddc->is_R_local; 87674ae819SStefano Zampini } 88*c7017625SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 89185763b3SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 90674ae819SStefano Zampini 91674ae819SStefano Zampini /* Get null space vecs */ 92*c7017625SStefano Zampini ierr = MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 93674ae819SStefano Zampini basis_size = nnsp_size; 94*c7017625SStefano Zampini if (nnsp_has_cnst) basis_size++; 95674ae819SStefano Zampini 96674ae819SStefano Zampini /* Create shell ctx */ 97854ce69bSBarry Smith ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 98*c7017625SStefano Zampini shell_ctx->apply_scaling = needscaling; 99674ae819SStefano Zampini 100674ae819SStefano Zampini /* Create work vectors in shell context */ 101674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 102674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 103674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 104674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 105674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 106674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 107674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 108674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 109674ae819SStefano Zampini 110674ae819SStefano Zampini /* Allocate workspace */ 111674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);CHKERRQ(ierr); 112674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 113674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 114674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 115674ae819SStefano Zampini 116*c7017625SStefano Zampini /* Restrict local null space on selected dofs 117674ae819SStefano Zampini and compute matrices N and K*N */ 118*c7017625SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 119674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 120*c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 121*c7017625SStefano Zampini ierr = VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122*c7017625SStefano Zampini ierr = VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123*c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 124*c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 125*c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 126*c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 127674ae819SStefano Zampini } 128674ae819SStefano Zampini if (nnsp_has_cnst) { 129*c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 130*c7017625SStefano Zampini ierr = VecSet(shell_ctx->work_full_1,one);CHKERRQ(ierr); 131*c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 132*c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 133*c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 134*c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 135674ae819SStefano Zampini } 136674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 137674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 138674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 139674ae819SStefano Zampini 140674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 141674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 142674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 143674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 144674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 145674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 146785e854fSJed Brown ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr); 147674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 148674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 149674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 150674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 151674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 152674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 153674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 154674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 155674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 156674ae819SStefano Zampini } 157674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 158674ae819SStefano Zampini } 159674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 160674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 161674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 162674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 163674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 164674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 165674ae819SStefano Zampini 166674ae819SStefano Zampini /* Rebuild local PC */ 167185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 168674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 169674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 17023ee1639SBarry Smith ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr); 171674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 172674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 173674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 174674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 175674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 176185763b3SStefano Zampini ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr); 177185763b3SStefano Zampini ierr = KSPSetUp(local_ksp);CHKERRQ(ierr); 178*c7017625SStefano Zampini 179*c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 180*c7017625SStefano Zampini if (needscaling) { 181674ae819SStefano Zampini KSP check_ksp; 182*c7017625SStefano Zampini Vec *workv; 183*c7017625SStefano Zampini 184*c7017625SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 185*c7017625SStefano Zampini ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 186*c7017625SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 187*c7017625SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-4,1.e-12,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr); 188*c7017625SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 189*c7017625SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 190*c7017625SStefano Zampini if (issym) { 191*c7017625SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 192*c7017625SStefano Zampini } 193*c7017625SStefano Zampini ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr); 194*c7017625SStefano Zampini ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr); 195*c7017625SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 196*c7017625SStefano Zampini ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr); 197*c7017625SStefano Zampini ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr); 198*c7017625SStefano Zampini ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr); 199*c7017625SStefano Zampini ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr); 200*c7017625SStefano Zampini ierr = VecNorm(workv[0],NORM_INFINITY,&test_err);CHKERRQ(ierr); 201*c7017625SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 202*c7017625SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 203*c7017625SStefano Zampini if (pcbddc->dbg_flag) { 204*c7017625SStefano Zampini if (isdir) { 205*c7017625SStefano 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); 206*c7017625SStefano Zampini } else { 207*c7017625SStefano 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); 208*c7017625SStefano Zampini } 209*c7017625SStefano Zampini } 210*c7017625SStefano Zampini shell_ctx->scale = 1./lambda_max; 211*c7017625SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 212*c7017625SStefano Zampini ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr); 213*c7017625SStefano Zampini } 214*c7017625SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 215*c7017625SStefano Zampini PetscFunctionReturn(0); 216*c7017625SStefano Zampini } 217*c7017625SStefano Zampini 218*c7017625SStefano Zampini #undef __FUNCT__ 219*c7017625SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceCheckCorrection" 220*c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir) 221*c7017625SStefano Zampini { 222*c7017625SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 223*c7017625SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 224*c7017625SStefano Zampini KSP check_ksp,local_ksp; 225*c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 226*c7017625SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 227674ae819SStefano Zampini PC check_pc; 228*c7017625SStefano Zampini Mat test_mat,local_mat; 229674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 230674ae819SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 231*c7017625SStefano Zampini Vec work1,work2,work3; 232*c7017625SStefano Zampini PetscInt k; 233*c7017625SStefano Zampini PetscErrorCode ierr; 234674ae819SStefano Zampini 235*c7017625SStefano Zampini PetscFunctionBegin; 236*c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 237*c7017625SStefano Zampini if (!NullSpace) PetscFunctionReturn(0); 238*c7017625SStefano Zampini if (!pcbddc->dbg_flag) PetscFunctionReturn(0); 239*c7017625SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; 240*c7017625SStefano Zampini else local_ksp = pcbddc->ksp_R; 241*c7017625SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr); 242185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr); 243*c7017625SStefano Zampini ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr); 244674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 245674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 246674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 247*c7017625SStefano Zampini 248674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 249674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 250674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 251674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 252674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 253*c7017625SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 254674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 255185763b3SStefano Zampini if (isdir) { 256*c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 257674ae819SStefano Zampini } else { 258*c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 259674ae819SStefano Zampini } 260674ae819SStefano Zampini 261674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 262*c7017625SStefano Zampini ierr = MatShift(test_mat,1.);CHKERRQ(ierr); 263674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 264674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 265*c7017625SStefano Zampini if (isdir) { 266*c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 267*c7017625SStefano Zampini } else { 268*c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 269*c7017625SStefano Zampini } 270674ae819SStefano Zampini 271674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 272674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 273422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 27423ee1639SBarry Smith ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 275*c7017625SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-10,1.e-8,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 276674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 277674ae819SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 278674ae819SStefano Zampini if (issym) { 279674ae819SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 280674ae819SStefano Zampini } 281674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 282674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 283674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 284674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 286*c7017625SStefano Zampini ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr); 287674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 288674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 289674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 290*c7017625SStefano Zampini if (isdir) { 291*c7017625SStefano Zampini ierr = 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);CHKERRQ(ierr); 292*c7017625SStefano Zampini } else { 293*c7017625SStefano Zampini ierr = 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);CHKERRQ(ierr); 294*c7017625SStefano Zampini } 295674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 296674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 297674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 298674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 299674ae819SStefano Zampini PetscFunctionReturn(0); 300674ae819SStefano Zampini } 301