1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 4*bd5e1169SStefano Zampini static PetscErrorCode PCBDDCViewNullSpaceCorrectionPC(PC pc,PetscViewer view) 5*bd5e1169SStefano Zampini { 6*bd5e1169SStefano Zampini PetscErrorCode ierr; 7*bd5e1169SStefano Zampini PetscBool isascii; 8*bd5e1169SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 9*bd5e1169SStefano Zampini 10*bd5e1169SStefano Zampini PetscFunctionBegin; 11*bd5e1169SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 12*bd5e1169SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)view,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 13*bd5e1169SStefano Zampini if (isascii) { 14*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"inner preconditioner:\n");CHKERRQ(ierr); 15*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 16*bd5e1169SStefano Zampini ierr = PCView(pc_ctx->local_pc,view);CHKERRQ(ierr); 17*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 18*bd5e1169SStefano Zampini 19*bd5e1169SStefano Zampini ierr = PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 20*bd5e1169SStefano Zampini 21*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Lbasis:\n");CHKERRQ(ierr); 22*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 23*bd5e1169SStefano Zampini ierr = MatView(pc_ctx->Lbasis_mat,view);CHKERRQ(ierr); 24*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 25*bd5e1169SStefano Zampini 26*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Kbasis:\n");CHKERRQ(ierr); 27*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 28*bd5e1169SStefano Zampini ierr = MatView(pc_ctx->Kbasis_mat,view);CHKERRQ(ierr); 29*bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 30*bd5e1169SStefano Zampini 31*bd5e1169SStefano Zampini ierr = PetscViewerPopFormat(view);CHKERRQ(ierr); 32*bd5e1169SStefano Zampini } 33*bd5e1169SStefano Zampini PetscFunctionReturn(0); 34*bd5e1169SStefano Zampini } 35*bd5e1169SStefano Zampini 36674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 37674ae819SStefano Zampini { 38674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 39674ae819SStefano Zampini PetscErrorCode ierr; 40674ae819SStefano Zampini 41674ae819SStefano Zampini PetscFunctionBegin; 42674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 43674ae819SStefano Zampini /* E */ 44674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 45674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 46674ae819SStefano Zampini /* P^-1 */ 47674ae819SStefano Zampini ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 48674ae819SStefano Zampini /* E^T */ 49674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 50674ae819SStefano Zampini ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 51674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 52674ae819SStefano Zampini /* Sum contributions */ 53674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 54c7017625SStefano Zampini if (pc_ctx->apply_scaling) { 55c7017625SStefano Zampini ierr = VecScale(y,pc_ctx->scale);CHKERRQ(ierr); 56c7017625SStefano Zampini } 57674ae819SStefano Zampini PetscFunctionReturn(0); 58674ae819SStefano Zampini } 59674ae819SStefano Zampini 60674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 61674ae819SStefano Zampini { 62674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 63674ae819SStefano Zampini PetscErrorCode ierr; 64674ae819SStefano Zampini 65674ae819SStefano Zampini PetscFunctionBegin; 66674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 67674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 68674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 69674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 70674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 71674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 72674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 73674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 74674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 75674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 76674ae819SStefano Zampini PetscFunctionReturn(0); 77674ae819SStefano Zampini } 78674ae819SStefano Zampini 79c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc, PetscBool isdir, PetscBool needscaling) 80674ae819SStefano Zampini { 81674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 82674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 83674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 84c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 85185763b3SStefano Zampini KSP local_ksp; 86674ae819SStefano Zampini PC newpc; 87674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 88674ae819SStefano Zampini Mat local_mat,local_pmat,small_mat,inv_small_mat; 89674ae819SStefano Zampini const Vec *nullvecs; 90674ae819SStefano Zampini VecScatter scatter_ctx; 91c7017625SStefano Zampini IS is_aux,local_dofs; 92674ae819SStefano Zampini MatFactorInfo matinfo; 93674ae819SStefano Zampini PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 94674ae819SStefano Zampini PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 95185763b3SStefano Zampini PetscInt basis_dofs,basis_size,nnsp_size,i,k; 96674ae819SStefano Zampini PetscBool nnsp_has_cnst; 97c7017625SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 98674ae819SStefano Zampini PetscErrorCode ierr; 99674ae819SStefano Zampini 100674ae819SStefano Zampini PetscFunctionBegin; 101c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 10235529e7bSStefano Zampini if (!NullSpace) { 10335529e7bSStefano Zampini if (pcbddc->dbg_flag) { 10435529e7bSStefano 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); 10535529e7bSStefano Zampini } 10635529e7bSStefano Zampini PetscFunctionReturn(0); 10735529e7bSStefano Zampini } 108674ae819SStefano Zampini /* Infer the local solver */ 109185763b3SStefano Zampini if (isdir) { 110674ae819SStefano Zampini /* Dirichlet solver */ 111185763b3SStefano Zampini local_ksp = pcbddc->ksp_D; 112c7017625SStefano Zampini local_dofs = pcis->is_I_local; 113674ae819SStefano Zampini } else { 114185763b3SStefano Zampini /* Neumann solver */ 115185763b3SStefano Zampini local_ksp = pcbddc->ksp_R; 116c7017625SStefano Zampini local_dofs = pcbddc->is_R_local; 117674ae819SStefano Zampini } 118c7017625SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 119185763b3SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 120674ae819SStefano Zampini 121674ae819SStefano Zampini /* Get null space vecs */ 122c7017625SStefano Zampini ierr = MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 123674ae819SStefano Zampini basis_size = nnsp_size; 124c7017625SStefano Zampini if (nnsp_has_cnst) basis_size++; 125674ae819SStefano Zampini 126674ae819SStefano Zampini /* Create shell ctx */ 127854ce69bSBarry Smith ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 128c7017625SStefano Zampini shell_ctx->apply_scaling = needscaling; 129*bd5e1169SStefano Zampini shell_ctx->scale = 1.0; 130674ae819SStefano Zampini 131674ae819SStefano Zampini /* Create work vectors in shell context */ 132674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 133674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 134674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 135674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 136674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 137674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 138674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 139674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 140674ae819SStefano Zampini 141674ae819SStefano Zampini /* Allocate workspace */ 142674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);CHKERRQ(ierr); 143674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 144674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 145674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 146674ae819SStefano Zampini 147c7017625SStefano Zampini /* Restrict local null space on selected dofs 148674ae819SStefano Zampini and compute matrices N and K*N */ 149c7017625SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 150674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 151c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 152c7017625SStefano Zampini ierr = VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 153c7017625SStefano Zampini ierr = VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 154c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 155c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 156c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 157c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 158674ae819SStefano Zampini } 159674ae819SStefano Zampini if (nnsp_has_cnst) { 160c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 161c7017625SStefano Zampini ierr = VecSet(shell_ctx->work_full_1,one);CHKERRQ(ierr); 162c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 163c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 164c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 165c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 166674ae819SStefano Zampini } 167674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 168674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 169674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 170674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 171674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 172674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 173674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 174674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 175674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 176785e854fSJed Brown ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr); 177674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 178674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 179674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 180674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 181674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 182674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 183674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 184674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 185674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 186674ae819SStefano Zampini } 187674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 188674ae819SStefano Zampini } 189674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 190674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 191674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 192674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 193674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 194674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 195674ae819SStefano Zampini 196674ae819SStefano Zampini /* Rebuild local PC */ 197185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 198674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 199674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 20023ee1639SBarry Smith ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr); 201674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 202674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 203674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 204*bd5e1169SStefano Zampini ierr = PCShellSetView(newpc,PCBDDCViewNullSpaceCorrectionPC);CHKERRQ(ierr); 205674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 206*bd5e1169SStefano Zampini ierr = PCSetUp(shell_ctx->local_pc);CHKERRQ(ierr); 207674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 208185763b3SStefano Zampini ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr); 209185763b3SStefano Zampini ierr = KSPSetUp(local_ksp);CHKERRQ(ierr); 210c7017625SStefano Zampini 211c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 212c7017625SStefano Zampini if (needscaling) { 213674ae819SStefano Zampini KSP check_ksp; 214c7017625SStefano Zampini Vec *workv; 215*bd5e1169SStefano Zampini const char* prefix; 216c7017625SStefano Zampini 217c7017625SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 218*bd5e1169SStefano Zampini ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 219*bd5e1169SStefano Zampini ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 220*bd5e1169SStefano Zampini ierr = KSPAppendOptionsPrefix(check_ksp,"approxscale_");CHKERRQ(ierr); 221*bd5e1169SStefano Zampini ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr); 222c7017625SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 223*bd5e1169SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 224c7017625SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 225c7017625SStefano Zampini ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr); 226*bd5e1169SStefano Zampini ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 227c7017625SStefano Zampini ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr); 228c7017625SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 229c7017625SStefano Zampini ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr); 230c7017625SStefano Zampini ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr); 231c7017625SStefano Zampini ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr); 232c7017625SStefano Zampini ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr); 233c7017625SStefano Zampini ierr = VecNorm(workv[0],NORM_INFINITY,&test_err);CHKERRQ(ierr); 234c7017625SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 235c7017625SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 236c7017625SStefano Zampini if (pcbddc->dbg_flag) { 237c7017625SStefano Zampini if (isdir) { 238c7017625SStefano 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); 239c7017625SStefano Zampini } else { 240c7017625SStefano 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); 241c7017625SStefano Zampini } 242c7017625SStefano Zampini } 243c7017625SStefano Zampini shell_ctx->scale = 1./lambda_max; 244c7017625SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 245c7017625SStefano Zampini ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr); 246c7017625SStefano Zampini } 247c7017625SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 248c7017625SStefano Zampini PetscFunctionReturn(0); 249c7017625SStefano Zampini } 250c7017625SStefano Zampini 251c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir) 252c7017625SStefano Zampini { 253c7017625SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 254c7017625SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 255c7017625SStefano Zampini KSP check_ksp,local_ksp; 256c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 257c7017625SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 258674ae819SStefano Zampini PC check_pc; 259c7017625SStefano Zampini Mat test_mat,local_mat; 260674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 261c7017625SStefano Zampini Vec work1,work2,work3; 262c7017625SStefano Zampini PetscInt k; 263*bd5e1169SStefano Zampini const char *prefix; 264c7017625SStefano Zampini PetscErrorCode ierr; 265674ae819SStefano Zampini 266c7017625SStefano Zampini PetscFunctionBegin; 267c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 268c7017625SStefano Zampini if (!NullSpace) PetscFunctionReturn(0); 269c7017625SStefano Zampini if (!pcbddc->dbg_flag) PetscFunctionReturn(0); 270c7017625SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; 271c7017625SStefano Zampini else local_ksp = pcbddc->ksp_R; 272c7017625SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr); 273185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr); 274c7017625SStefano Zampini ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr); 275674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 276674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 277674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 278c7017625SStefano Zampini 279674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 280674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 281674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 282674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 283674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 284c7017625SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 286185763b3SStefano Zampini if (isdir) { 287c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 288674ae819SStefano Zampini } else { 289c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 290674ae819SStefano Zampini } 291674ae819SStefano Zampini 292674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 293c7017625SStefano Zampini ierr = MatShift(test_mat,1.);CHKERRQ(ierr); 294674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 295674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 296c7017625SStefano Zampini if (isdir) { 297c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 298c7017625SStefano Zampini } else { 299c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 300c7017625SStefano Zampini } 301674ae819SStefano Zampini 302674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 303674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 304*bd5e1169SStefano Zampini ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 305*bd5e1169SStefano Zampini ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 306*bd5e1169SStefano Zampini ierr = KSPAppendOptionsPrefix(check_ksp,"approxcheck_");CHKERRQ(ierr); 307422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 30823ee1639SBarry Smith ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 309*bd5e1169SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 310674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 311*bd5e1169SStefano Zampini ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 312674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 313674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 314674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 315674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 316674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 317c7017625SStefano Zampini ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr); 318674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 319674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 320674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 321c7017625SStefano Zampini if (isdir) { 322c7017625SStefano 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); 323c7017625SStefano Zampini } else { 324c7017625SStefano 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); 325c7017625SStefano Zampini } 326674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 327674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 328674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 329674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 330674ae819SStefano Zampini PetscFunctionReturn(0); 331674ae819SStefano Zampini } 332