1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 4674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 5674ae819SStefano Zampini { 6674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 7674ae819SStefano Zampini PetscErrorCode ierr; 8674ae819SStefano Zampini 9674ae819SStefano Zampini PetscFunctionBegin; 10674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 11674ae819SStefano Zampini /* E */ 12674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 13674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 14674ae819SStefano Zampini /* P^-1 */ 15674ae819SStefano Zampini ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 16674ae819SStefano Zampini /* E^T */ 17674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 18674ae819SStefano Zampini ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 19674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 20674ae819SStefano Zampini /* Sum contributions */ 21674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 22c7017625SStefano Zampini if (pc_ctx->apply_scaling) { 23c7017625SStefano Zampini ierr = VecScale(y,pc_ctx->scale);CHKERRQ(ierr); 24c7017625SStefano Zampini } 25674ae819SStefano Zampini PetscFunctionReturn(0); 26674ae819SStefano Zampini } 27674ae819SStefano Zampini 28674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 29674ae819SStefano Zampini { 30674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 31674ae819SStefano Zampini PetscErrorCode ierr; 32674ae819SStefano Zampini 33674ae819SStefano Zampini PetscFunctionBegin; 34674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 35674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 36674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 37674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 38674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 39674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 40674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 41674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 42674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 43674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 44674ae819SStefano Zampini PetscFunctionReturn(0); 45674ae819SStefano Zampini } 46674ae819SStefano Zampini 47c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 48674ae819SStefano Zampini { 49674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 50674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 51674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 52c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 53185763b3SStefano Zampini KSP local_ksp; 54674ae819SStefano Zampini PC newpc; 55674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 56674ae819SStefano Zampini Mat local_mat,local_pmat,small_mat,inv_small_mat; 57674ae819SStefano Zampini const Vec *nullvecs; 58674ae819SStefano Zampini VecScatter scatter_ctx; 59c7017625SStefano Zampini IS is_aux,local_dofs; 60674ae819SStefano Zampini MatFactorInfo matinfo; 61674ae819SStefano Zampini PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 62674ae819SStefano Zampini PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 63185763b3SStefano Zampini PetscInt basis_dofs,basis_size,nnsp_size,i,k; 64674ae819SStefano Zampini PetscBool nnsp_has_cnst; 65c7017625SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 66c7017625SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 67674ae819SStefano Zampini PetscErrorCode ierr; 68674ae819SStefano Zampini 69674ae819SStefano Zampini PetscFunctionBegin; 70c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 71*35529e7bSStefano Zampini if (!NullSpace) { 72*35529e7bSStefano Zampini if (pcbddc->dbg_flag) { 73*35529e7bSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d doesn't have local nullspace: no need for correction in %s solver \n",PetscGlobalRank,isdir ? "Dirichlet" : "Neumann");CHKERRQ(ierr); 74*35529e7bSStefano Zampini } 75*35529e7bSStefano Zampini PetscFunctionReturn(0); 76*35529e7bSStefano Zampini } 77674ae819SStefano Zampini /* Infer the local solver */ 78185763b3SStefano Zampini if (isdir) { 79674ae819SStefano Zampini /* Dirichlet solver */ 80185763b3SStefano Zampini local_ksp = pcbddc->ksp_D; 81c7017625SStefano Zampini local_dofs = pcis->is_I_local; 82674ae819SStefano Zampini } else { 83185763b3SStefano Zampini /* Neumann solver */ 84185763b3SStefano Zampini local_ksp = pcbddc->ksp_R; 85c7017625SStefano Zampini local_dofs = pcbddc->is_R_local; 86674ae819SStefano Zampini } 87c7017625SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 88185763b3SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 89674ae819SStefano Zampini 90674ae819SStefano Zampini /* Get null space vecs */ 91c7017625SStefano Zampini ierr = MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 92674ae819SStefano Zampini basis_size = nnsp_size; 93c7017625SStefano Zampini if (nnsp_has_cnst) basis_size++; 94674ae819SStefano Zampini 95674ae819SStefano Zampini /* Create shell ctx */ 96854ce69bSBarry Smith ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 97c7017625SStefano Zampini shell_ctx->apply_scaling = needscaling; 98674ae819SStefano Zampini 99674ae819SStefano Zampini /* Create work vectors in shell context */ 100674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 101674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 102674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 103674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 104674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 105674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 106674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 107674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 108674ae819SStefano Zampini 109674ae819SStefano Zampini /* Allocate workspace */ 110674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);CHKERRQ(ierr); 111674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 112674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 113674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 114674ae819SStefano Zampini 115c7017625SStefano Zampini /* Restrict local null space on selected dofs 116674ae819SStefano Zampini and compute matrices N and K*N */ 117c7017625SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 118674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 119c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 120c7017625SStefano Zampini ierr = VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121c7017625SStefano Zampini ierr = VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 123c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 124c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 125c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 126674ae819SStefano Zampini } 127674ae819SStefano Zampini if (nnsp_has_cnst) { 128c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 129c7017625SStefano Zampini ierr = VecSet(shell_ctx->work_full_1,one);CHKERRQ(ierr); 130c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 131c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 132c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 133c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 134674ae819SStefano Zampini } 135674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 136674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 137674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 138674ae819SStefano Zampini 139674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 140674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 141674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 142674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 143674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 144674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 145785e854fSJed Brown ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr); 146674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 147674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 148674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 149674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 150674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 151674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 152674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 153674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 154674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 155674ae819SStefano Zampini } 156674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 157674ae819SStefano Zampini } 158674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 159674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 160674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 161674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 162674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 163674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 164674ae819SStefano Zampini 165674ae819SStefano Zampini /* Rebuild local PC */ 166185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 167674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 168674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 16923ee1639SBarry Smith ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr); 170674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 171674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 172674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 173674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 174674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 175185763b3SStefano Zampini ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr); 176185763b3SStefano Zampini ierr = KSPSetUp(local_ksp);CHKERRQ(ierr); 177c7017625SStefano Zampini 178c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 179c7017625SStefano Zampini if (needscaling) { 180674ae819SStefano Zampini KSP check_ksp; 181c7017625SStefano Zampini Vec *workv; 182c7017625SStefano Zampini 183c7017625SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 184c7017625SStefano Zampini ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 185c7017625SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 186c7017625SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-4,1.e-12,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr); 187c7017625SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 188c7017625SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 189c7017625SStefano Zampini if (issym) { 190c7017625SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 191c7017625SStefano Zampini } 192c7017625SStefano Zampini ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr); 193c7017625SStefano Zampini ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr); 194c7017625SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 195c7017625SStefano Zampini ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr); 196c7017625SStefano Zampini ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr); 197c7017625SStefano Zampini ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr); 198c7017625SStefano Zampini ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr); 199c7017625SStefano Zampini ierr = VecNorm(workv[0],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) { 204c7017625SStefano 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 { 206c7017625SStefano 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 } 209c7017625SStefano Zampini shell_ctx->scale = 1./lambda_max; 210c7017625SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 211c7017625SStefano Zampini ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr); 212c7017625SStefano Zampini } 213c7017625SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 214c7017625SStefano Zampini PetscFunctionReturn(0); 215c7017625SStefano Zampini } 216c7017625SStefano Zampini 217c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir) 218c7017625SStefano Zampini { 219c7017625SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 220c7017625SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 221c7017625SStefano Zampini KSP check_ksp,local_ksp; 222c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 223c7017625SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 224674ae819SStefano Zampini PC check_pc; 225c7017625SStefano Zampini Mat test_mat,local_mat; 226674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 227674ae819SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 228c7017625SStefano Zampini Vec work1,work2,work3; 229c7017625SStefano Zampini PetscInt k; 230c7017625SStefano Zampini PetscErrorCode ierr; 231674ae819SStefano Zampini 232c7017625SStefano Zampini PetscFunctionBegin; 233c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 234c7017625SStefano Zampini if (!NullSpace) PetscFunctionReturn(0); 235c7017625SStefano Zampini if (!pcbddc->dbg_flag) PetscFunctionReturn(0); 236c7017625SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; 237c7017625SStefano Zampini else local_ksp = pcbddc->ksp_R; 238c7017625SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr); 239185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr); 240c7017625SStefano Zampini ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr); 241674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 242674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 243674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 244c7017625SStefano Zampini 245674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 246674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 247674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 248674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 249674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 250c7017625SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 251674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 252185763b3SStefano Zampini if (isdir) { 253c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 254674ae819SStefano Zampini } else { 255c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 256674ae819SStefano Zampini } 257674ae819SStefano Zampini 258674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 259c7017625SStefano Zampini ierr = MatShift(test_mat,1.);CHKERRQ(ierr); 260674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 261674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 262c7017625SStefano Zampini if (isdir) { 263c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 264c7017625SStefano Zampini } else { 265c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 266c7017625SStefano Zampini } 267674ae819SStefano Zampini 268674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 269674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 270422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 27123ee1639SBarry Smith ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 272c7017625SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-10,1.e-8,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 273674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 274674ae819SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 275674ae819SStefano Zampini if (issym) { 276674ae819SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 277674ae819SStefano Zampini } 278674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 279674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 280674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 281674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 282674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 283c7017625SStefano Zampini ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr); 284674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 286674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 287c7017625SStefano Zampini if (isdir) { 288c7017625SStefano 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); 289c7017625SStefano Zampini } else { 290c7017625SStefano 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); 291c7017625SStefano Zampini } 292674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 293674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 294674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 295674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 296674ae819SStefano Zampini PetscFunctionReturn(0); 297674ae819SStefano Zampini } 298