1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 4bd5e1169SStefano Zampini static PetscErrorCode PCBDDCViewNullSpaceCorrectionPC(PC pc,PetscViewer view) 5bd5e1169SStefano Zampini { 6bd5e1169SStefano Zampini PetscErrorCode ierr; 7bd5e1169SStefano Zampini PetscBool isascii; 8bd5e1169SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 9bd5e1169SStefano Zampini 10bd5e1169SStefano Zampini PetscFunctionBegin; 11bd5e1169SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 12bd5e1169SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)view,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 13bd5e1169SStefano Zampini if (isascii) { 14bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"inner preconditioner:\n");CHKERRQ(ierr); 15bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 16bd5e1169SStefano Zampini ierr = PCView(pc_ctx->local_pc,view);CHKERRQ(ierr); 17bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 18bd5e1169SStefano Zampini 19bd5e1169SStefano Zampini ierr = PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 20bd5e1169SStefano Zampini 21bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Lbasis:\n");CHKERRQ(ierr); 22bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 23bd5e1169SStefano Zampini ierr = MatView(pc_ctx->Lbasis_mat,view);CHKERRQ(ierr); 24bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 25bd5e1169SStefano Zampini 26bd5e1169SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Kbasis:\n");CHKERRQ(ierr); 27bd5e1169SStefano Zampini ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr); 28bd5e1169SStefano Zampini ierr = MatView(pc_ctx->Kbasis_mat,view);CHKERRQ(ierr); 29bd5e1169SStefano Zampini ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr); 30bd5e1169SStefano Zampini 31bd5e1169SStefano Zampini ierr = PetscViewerPopFormat(view);CHKERRQ(ierr); 32bd5e1169SStefano Zampini } 33bd5e1169SStefano Zampini PetscFunctionReturn(0); 34bd5e1169SStefano Zampini } 35bd5e1169SStefano 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; 89*2958b453SStefano Zampini 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) { 103*2958b453SStefano Zampini ierr = MatGetNearNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 104*2958b453SStefano Zampini } 105*2958b453SStefano Zampini if (!NullSpace) { 10635529e7bSStefano Zampini if (pcbddc->dbg_flag) { 107*2958b453SStefano 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); 10835529e7bSStefano Zampini } 10935529e7bSStefano Zampini PetscFunctionReturn(0); 11035529e7bSStefano Zampini } 111*2958b453SStefano Zampini 112674ae819SStefano Zampini /* Infer the local solver */ 113185763b3SStefano Zampini if (isdir) { 114674ae819SStefano Zampini /* Dirichlet solver */ 115185763b3SStefano Zampini local_ksp = pcbddc->ksp_D; 116c7017625SStefano Zampini local_dofs = pcis->is_I_local; 117674ae819SStefano Zampini } else { 118185763b3SStefano Zampini /* Neumann solver */ 119185763b3SStefano Zampini local_ksp = pcbddc->ksp_R; 120c7017625SStefano Zampini local_dofs = pcbddc->is_R_local; 121674ae819SStefano Zampini } 122c7017625SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 123185763b3SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr); 124674ae819SStefano Zampini 125674ae819SStefano Zampini /* Get null space vecs */ 126*2958b453SStefano Zampini ierr = MatNullSpaceGetVecs(NullSpace,&nnsp_has_cnst,&nnsp_size,(const Vec**)&nullvecs);CHKERRQ(ierr); 127674ae819SStefano Zampini basis_size = nnsp_size; 128c7017625SStefano Zampini if (nnsp_has_cnst) basis_size++; 129674ae819SStefano Zampini 130674ae819SStefano Zampini /* Create shell ctx */ 131854ce69bSBarry Smith ierr = PetscNew(&shell_ctx);CHKERRQ(ierr); 132c7017625SStefano Zampini shell_ctx->apply_scaling = needscaling; 133bd5e1169SStefano Zampini shell_ctx->scale = 1.0; 134674ae819SStefano Zampini 135674ae819SStefano Zampini /* Create work vectors in shell context */ 136674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 137674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 138674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 139674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 140674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 141674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 142674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 143674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 144674ae819SStefano Zampini 145674ae819SStefano Zampini /* Allocate workspace */ 146674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat);CHKERRQ(ierr); 147674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 148674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 149674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 150674ae819SStefano Zampini 151c7017625SStefano Zampini /* Restrict local null space on selected dofs 152674ae819SStefano Zampini and compute matrices N and K*N */ 153c7017625SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,shell_ctx->work_full_1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 154674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 155c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 156c7017625SStefano Zampini ierr = VecScatterBegin(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 157c7017625SStefano Zampini ierr = VecScatterEnd(scatter_ctx,nullvecs[k],shell_ctx->work_full_1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 158c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 159674ae819SStefano Zampini } 160674ae819SStefano Zampini if (nnsp_has_cnst) { 161c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 162c7017625SStefano Zampini ierr = VecSet(shell_ctx->work_full_1,one);CHKERRQ(ierr); 163*2958b453SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 164*2958b453SStefano Zampini } 165*2958b453SStefano Zampini 166*2958b453SStefano Zampini ierr = PetscMalloc1(basis_size,&nullvecs);CHKERRQ(ierr); 167*2958b453SStefano Zampini for (k=0;k<basis_size;k++) { 168*2958b453SStefano Zampini ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,basis_dofs,basis_mat + k*basis_dofs,&nullvecs[k]);CHKERRQ(ierr); 169*2958b453SStefano Zampini } 170*2958b453SStefano Zampini ierr = PCBDDCOrthonormalizeVecs(basis_size,nullvecs);CHKERRQ(ierr); 171*2958b453SStefano Zampini ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_FALSE,basis_size,nullvecs,&NullSpace);CHKERRQ(ierr); 172*2958b453SStefano Zampini ierr = MatSetNearNullSpace(local_mat,NullSpace);CHKERRQ(ierr); 173*2958b453SStefano Zampini ierr = MatNullSpaceDestroy(&NullSpace);CHKERRQ(ierr); 174*2958b453SStefano Zampini for (k=0;k<basis_size;k++) { 175*2958b453SStefano Zampini ierr = VecDestroy(&nullvecs[k]);CHKERRQ(ierr); 176*2958b453SStefano Zampini } 177*2958b453SStefano Zampini ierr = PetscFree(nullvecs);CHKERRQ(ierr); 178*2958b453SStefano Zampini 179*2958b453SStefano Zampini for (k=0;k<basis_size;k++) { 180*2958b453SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 181c7017625SStefano Zampini ierr = VecPlaceArray(shell_ctx->work_full_2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 182c7017625SStefano Zampini ierr = MatMult(local_mat,shell_ctx->work_full_1,shell_ctx->work_full_2);CHKERRQ(ierr); 183c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr); 184c7017625SStefano Zampini ierr = VecResetArray(shell_ctx->work_full_2);CHKERRQ(ierr); 185674ae819SStefano Zampini } 186*2958b453SStefano Zampini 187674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 188674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 189674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 190674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 191674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 192674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 193674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 194674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 195674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 196785e854fSJed Brown ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr); 197674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 198674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 199674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 200674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 201674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 202674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 203674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 204674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 205674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 206674ae819SStefano Zampini } 207674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 208674ae819SStefano Zampini } 209674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 210674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 211674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 212674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 213674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 214674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 215674ae819SStefano Zampini 216674ae819SStefano Zampini /* Rebuild local PC */ 217185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 218674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 219674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 22023ee1639SBarry Smith ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr); 221674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 222674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 223674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 224bd5e1169SStefano Zampini ierr = PCShellSetView(newpc,PCBDDCViewNullSpaceCorrectionPC);CHKERRQ(ierr); 225674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 226674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 227185763b3SStefano Zampini ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr); 228185763b3SStefano Zampini ierr = KSPSetUp(local_ksp);CHKERRQ(ierr); 229c7017625SStefano Zampini 230c7017625SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 231c7017625SStefano Zampini if (needscaling) { 232674ae819SStefano Zampini KSP check_ksp; 233c7017625SStefano Zampini Vec *workv; 234bd5e1169SStefano Zampini const char* prefix; 235c7017625SStefano Zampini 236c7017625SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 237bd5e1169SStefano Zampini ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 238bd5e1169SStefano Zampini ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 239bd5e1169SStefano Zampini ierr = KSPAppendOptionsPrefix(check_ksp,"approxscale_");CHKERRQ(ierr); 240bd5e1169SStefano Zampini ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr); 241c7017625SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 242bd5e1169SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 243c7017625SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 244c7017625SStefano Zampini ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr); 245bd5e1169SStefano Zampini ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 246c7017625SStefano Zampini ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr); 247c7017625SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 248c7017625SStefano Zampini ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr); 249c7017625SStefano Zampini ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr); 250c7017625SStefano Zampini ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr); 251c7017625SStefano Zampini ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr); 252c7017625SStefano Zampini ierr = VecNorm(workv[0],NORM_INFINITY,&test_err);CHKERRQ(ierr); 253c7017625SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 254c7017625SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 255c7017625SStefano Zampini if (pcbddc->dbg_flag) { 256c7017625SStefano Zampini if (isdir) { 257c7017625SStefano 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); 258c7017625SStefano Zampini } else { 259c7017625SStefano 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); 260c7017625SStefano Zampini } 261c7017625SStefano Zampini } 262c7017625SStefano Zampini shell_ctx->scale = 1./lambda_max; 263c7017625SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 264c7017625SStefano Zampini ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr); 265c7017625SStefano Zampini } 266c7017625SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 267c7017625SStefano Zampini PetscFunctionReturn(0); 268c7017625SStefano Zampini } 269c7017625SStefano Zampini 270c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir) 271c7017625SStefano Zampini { 272c7017625SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 273c7017625SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 274c7017625SStefano Zampini KSP check_ksp,local_ksp; 275c7017625SStefano Zampini MatNullSpace NullSpace = NULL; 276c7017625SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 277674ae819SStefano Zampini PC check_pc; 278c7017625SStefano Zampini Mat test_mat,local_mat; 279674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 280c7017625SStefano Zampini Vec work1,work2,work3; 281c7017625SStefano Zampini PetscInt k; 282bd5e1169SStefano Zampini const char *prefix; 283c7017625SStefano Zampini PetscErrorCode ierr; 284674ae819SStefano Zampini 285c7017625SStefano Zampini PetscFunctionBegin; 286c7017625SStefano Zampini ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 287*2958b453SStefano Zampini if (!NullSpace) { 288*2958b453SStefano Zampini ierr = MatGetNearNullSpace(matis->A,&NullSpace);CHKERRQ(ierr); 289*2958b453SStefano Zampini } 290*2958b453SStefano Zampini if (!NullSpace) { 291*2958b453SStefano Zampini PetscFunctionReturn(0); 292*2958b453SStefano Zampini } 293c7017625SStefano Zampini if (!pcbddc->dbg_flag) PetscFunctionReturn(0); 294c7017625SStefano Zampini if (isdir) local_ksp = pcbddc->ksp_D; 295c7017625SStefano Zampini else local_ksp = pcbddc->ksp_R; 296c7017625SStefano Zampini ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr); 297185763b3SStefano Zampini ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr); 298c7017625SStefano Zampini ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr); 299674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 300674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 301674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 302c7017625SStefano Zampini 303674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 304674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 305674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 306674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 307674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 308c7017625SStefano Zampini ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr); 309674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 310185763b3SStefano Zampini if (isdir) { 311c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 312674ae819SStefano Zampini } else { 313c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 314674ae819SStefano Zampini } 315674ae819SStefano Zampini 316674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 317c7017625SStefano Zampini ierr = MatShift(test_mat,1.);CHKERRQ(ierr); 318674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 319674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 320c7017625SStefano Zampini if (isdir) { 321c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 322c7017625SStefano Zampini } else { 323c7017625SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr); 324c7017625SStefano Zampini } 325674ae819SStefano Zampini 326674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 327674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 328bd5e1169SStefano Zampini ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr); 329bd5e1169SStefano Zampini ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr); 330bd5e1169SStefano Zampini ierr = KSPAppendOptionsPrefix(check_ksp,"approxcheck_");CHKERRQ(ierr); 331422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr); 33223ee1639SBarry Smith ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr); 333bd5e1169SStefano Zampini ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 334674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 335bd5e1169SStefano Zampini ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr); 336674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 337674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 338674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 339674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 340674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 341c7017625SStefano Zampini ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr); 342674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 343674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 344674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 345c7017625SStefano Zampini if (isdir) { 346c7017625SStefano 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); 347c7017625SStefano Zampini } else { 348c7017625SStefano 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); 349c7017625SStefano Zampini } 350674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 351674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 352674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 353674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 354674ae819SStefano Zampini PetscFunctionReturn(0); 355674ae819SStefano Zampini } 356