xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision c1a3ebd0a0f982a93efc4dd23617f70dafbfbe90)
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 = PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
15bd5e1169SStefano Zampini 
1615579a77SStefano Zampini     ierr = PetscViewerASCIIPrintf(view,"L:\n");CHKERRQ(ierr);
17bd5e1169SStefano Zampini     ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr);
18bd5e1169SStefano Zampini     ierr = MatView(pc_ctx->Lbasis_mat,view);CHKERRQ(ierr);
19bd5e1169SStefano Zampini     ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr);
20bd5e1169SStefano Zampini 
2115579a77SStefano Zampini     ierr = PetscViewerASCIIPrintf(view,"K:\n");CHKERRQ(ierr);
22bd5e1169SStefano Zampini     ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr);
23bd5e1169SStefano Zampini     ierr = MatView(pc_ctx->Kbasis_mat,view);CHKERRQ(ierr);
24bd5e1169SStefano Zampini     ierr = PetscViewerASCIIPopTab(view);CHKERRQ(ierr);
25bd5e1169SStefano Zampini 
26bd5e1169SStefano Zampini     ierr = PetscViewerPopFormat(view);CHKERRQ(ierr);
2715579a77SStefano Zampini 
2815579a77SStefano Zampini     ierr = PetscViewerASCIIPrintf(view,"inner preconditioner:\n");CHKERRQ(ierr);
2915579a77SStefano Zampini     ierr = PetscViewerASCIIPushTab(view);CHKERRQ(ierr);
3015579a77SStefano Zampini     ierr = PCView(pc_ctx->local_pc,view);CHKERRQ(ierr);
3115579a77SStefano Zampini     ierr = PetscViewerASCIIPopTab(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;
892958b453SStefano 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) {
1032958b453SStefano Zampini     ierr = MatGetNearNullSpace(matis->A,&NullSpace);CHKERRQ(ierr);
1042958b453SStefano Zampini   }
1052958b453SStefano Zampini   if (!NullSpace) {
10635529e7bSStefano Zampini     if (pcbddc->dbg_flag) {
1072958b453SStefano 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   }
1112958b453SStefano 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 */
1262958b453SStefano 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);
138c72cccf8SStefano Zampini   ierr = VecSetType(shell_ctx->work_small_1,((PetscObject)pcis->vec1_B)->type_name);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);
142c72cccf8SStefano Zampini   ierr = VecSetType(shell_ctx->work_full_1,((PetscObject)pcis->vec1_B)->type_name);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 */
15335928de7SBarry Smith   ierr = VecScatterCreateWithData(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);
1632958b453SStefano Zampini     ierr = VecResetArray(shell_ctx->work_full_1);CHKERRQ(ierr);
1642958b453SStefano Zampini   }
1652958b453SStefano Zampini 
1662958b453SStefano Zampini   ierr = PetscMalloc1(basis_size,&nullvecs);CHKERRQ(ierr);
1672958b453SStefano Zampini   for (k=0;k<basis_size;k++) {
1682958b453SStefano Zampini     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,basis_dofs,basis_mat + k*basis_dofs,&nullvecs[k]);CHKERRQ(ierr);
1692958b453SStefano Zampini   }
1702958b453SStefano Zampini   ierr = PCBDDCOrthonormalizeVecs(basis_size,nullvecs);CHKERRQ(ierr);
1712958b453SStefano Zampini   ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_FALSE,basis_size,nullvecs,&NullSpace);CHKERRQ(ierr);
172*c1a3ebd0SStefano Zampini   ierr = MatSetNearNullSpace(local_pmat,NullSpace);CHKERRQ(ierr);
1732958b453SStefano Zampini   ierr = MatNullSpaceDestroy(&NullSpace);CHKERRQ(ierr);
1742958b453SStefano Zampini   for (k=0;k<basis_size;k++) {
1752958b453SStefano Zampini     ierr = VecDestroy(&nullvecs[k]);CHKERRQ(ierr);
1762958b453SStefano Zampini   }
1772958b453SStefano Zampini   ierr = PetscFree(nullvecs);CHKERRQ(ierr);
1782958b453SStefano Zampini 
1792958b453SStefano Zampini   for (k=0;k<basis_size;k++) {
1802958b453SStefano 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   }
1862958b453SStefano 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);
220*c1a3ebd0SStefano Zampini   ierr = PCSetOperators(newpc,local_mat,local_pmat);CHKERRQ(ierr);
221674ae819SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
222674ae819SStefano Zampini   ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr);
22315579a77SStefano Zampini   if (isdir) {
22415579a77SStefano Zampini     ierr = PCShellSetName(newpc,"Nullspace corrected interior solve");CHKERRQ(ierr);
22515579a77SStefano Zampini   } else {
22615579a77SStefano Zampini     ierr = PCShellSetName(newpc,"Nullspace corrected correction solve");CHKERRQ(ierr);
22715579a77SStefano Zampini   }
228674ae819SStefano Zampini   ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr);
229bd5e1169SStefano Zampini   ierr = PCShellSetView(newpc,PCBDDCViewNullSpaceCorrectionPC);CHKERRQ(ierr);
230674ae819SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr);
231674ae819SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
232185763b3SStefano Zampini   ierr = KSPSetPC(local_ksp,newpc);CHKERRQ(ierr);
23315579a77SStefano Zampini   ierr = PetscObjectIncrementTabLevel((PetscObject)newpc,(PetscObject)local_ksp,0);CHKERRQ(ierr);
234185763b3SStefano Zampini   ierr = KSPSetUp(local_ksp);CHKERRQ(ierr);
235c7017625SStefano Zampini 
236c7017625SStefano Zampini   /* Create ksp object suitable for extreme eigenvalues' estimation */
237c7017625SStefano Zampini   if (needscaling) {
238674ae819SStefano Zampini     KSP         check_ksp;
239c7017625SStefano Zampini     Vec         *workv;
240bd5e1169SStefano Zampini     const char* prefix;
241c7017625SStefano Zampini 
242c7017625SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
24315579a77SStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,0);CHKERRQ(ierr);
244bd5e1169SStefano Zampini     ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr);
245bd5e1169SStefano Zampini     ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr);
246bd5e1169SStefano Zampini     ierr = KSPAppendOptionsPrefix(check_ksp,"approxscale_");CHKERRQ(ierr);
247bd5e1169SStefano Zampini     ierr = KSPSetErrorIfNotConverged(check_ksp,PETSC_FALSE);CHKERRQ(ierr);
248*c1a3ebd0SStefano Zampini     ierr = KSPSetOperators(check_ksp,local_mat,local_pmat);CHKERRQ(ierr);
249bd5e1169SStefano Zampini     ierr = KSPSetTolerances(check_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
250c7017625SStefano Zampini     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
251c7017625SStefano Zampini     ierr = VecDuplicateVecs(shell_ctx->work_full_1,2,&workv);CHKERRQ(ierr);
252bd5e1169SStefano Zampini     ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr);
253c7017625SStefano Zampini     ierr = KSPSetPC(check_ksp,newpc);CHKERRQ(ierr);
254c7017625SStefano Zampini     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
255c7017625SStefano Zampini     ierr = VecSetRandom(workv[1],NULL);CHKERRQ(ierr);
256c7017625SStefano Zampini     ierr = MatMult(local_mat,workv[1],workv[0]);CHKERRQ(ierr);
257c7017625SStefano Zampini     ierr = KSPSolve(check_ksp,workv[0],workv[0]);CHKERRQ(ierr);
258c0decd05SBarry Smith     ierr = KSPCheckSolve(check_ksp,pc,workv[0]);CHKERRQ(ierr);
259c7017625SStefano Zampini     ierr = VecAXPY(workv[0],-1.,workv[1]);CHKERRQ(ierr);
260c7017625SStefano Zampini     ierr = VecNorm(workv[0],NORM_INFINITY,&test_err);CHKERRQ(ierr);
261c7017625SStefano Zampini     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
262c7017625SStefano Zampini     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
263c7017625SStefano Zampini     if (pcbddc->dbg_flag) {
264c7017625SStefano Zampini       if (isdir) {
265c7017625SStefano 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);
266c7017625SStefano Zampini       } else {
267c7017625SStefano 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);
268c7017625SStefano Zampini       }
269c7017625SStefano Zampini     }
270c7017625SStefano Zampini     shell_ctx->scale = 1./lambda_max;
271c7017625SStefano Zampini     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
272c7017625SStefano Zampini     ierr = VecDestroyVecs(2,&workv);CHKERRQ(ierr);
273c7017625SStefano Zampini   }
274c7017625SStefano Zampini   ierr = PCDestroy(&newpc);CHKERRQ(ierr);
275c7017625SStefano Zampini   PetscFunctionReturn(0);
276c7017625SStefano Zampini }
277c7017625SStefano Zampini 
278c7017625SStefano Zampini PetscErrorCode PCBDDCNullSpaceCheckCorrection(PC pc, PetscBool isdir)
279c7017625SStefano Zampini {
280c7017625SStefano Zampini   PC_BDDC                  *pcbddc = (PC_BDDC*)pc->data;
281c7017625SStefano Zampini   Mat_IS                   *matis = (Mat_IS*)pc->pmat->data;
282c7017625SStefano Zampini   KSP                      check_ksp,local_ksp;
283c7017625SStefano Zampini   MatNullSpace             NullSpace = NULL;
284c7017625SStefano Zampini   NullSpaceCorrection_ctx  shell_ctx;
285674ae819SStefano Zampini   PC                       check_pc;
286c7017625SStefano Zampini   Mat                      test_mat,local_mat;
287674ae819SStefano Zampini   PetscReal                test_err,lambda_min,lambda_max;
288c7017625SStefano Zampini   Vec                      work1,work2,work3;
289c7017625SStefano Zampini   PetscInt                 k;
290bd5e1169SStefano Zampini   const char               *prefix;
291c7017625SStefano Zampini   PetscErrorCode           ierr;
292674ae819SStefano Zampini 
293c7017625SStefano Zampini   PetscFunctionBegin;
294c7017625SStefano Zampini   ierr = MatGetNullSpace(matis->A,&NullSpace);CHKERRQ(ierr);
2952958b453SStefano Zampini   if (!NullSpace) {
2962958b453SStefano Zampini     ierr = MatGetNearNullSpace(matis->A,&NullSpace);CHKERRQ(ierr);
2972958b453SStefano Zampini   }
2982958b453SStefano Zampini   if (!NullSpace) {
2992958b453SStefano Zampini     PetscFunctionReturn(0);
3002958b453SStefano Zampini   }
301c7017625SStefano Zampini   if (!pcbddc->dbg_flag) PetscFunctionReturn(0);
302c7017625SStefano Zampini   if (isdir) local_ksp = pcbddc->ksp_D;
303c7017625SStefano Zampini   else local_ksp = pcbddc->ksp_R;
304c7017625SStefano Zampini   ierr = KSPGetOperators(local_ksp,&local_mat,NULL);CHKERRQ(ierr);
305185763b3SStefano Zampini   ierr = KSPGetPC(local_ksp,&check_pc);CHKERRQ(ierr);
306c7017625SStefano Zampini   ierr = PCShellGetContext(check_pc,(void**)&shell_ctx);CHKERRQ(ierr);
307674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr);
308674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr);
309674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr);
310c7017625SStefano Zampini 
311674ae819SStefano Zampini   ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr);
312674ae819SStefano Zampini   ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr);
313674ae819SStefano Zampini   ierr = VecCopy(work1,work2);CHKERRQ(ierr);
314674ae819SStefano Zampini   ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
315*c1a3ebd0SStefano Zampini   ierr = VecScale(work3,1.0/shell_ctx->scale);CHKERRQ(ierr);
316674ae819SStefano Zampini   ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr);
317c7017625SStefano Zampini   ierr = VecAXPY(work1,-1.,work2);CHKERRQ(ierr);
318674ae819SStefano Zampini   ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
319185763b3SStefano Zampini   if (isdir) {
320c7017625SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
321674ae819SStefano Zampini   } else {
322c7017625SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace correction solver: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
323674ae819SStefano Zampini   }
324674ae819SStefano Zampini 
325674ae819SStefano Zampini   ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr);
326c7017625SStefano Zampini   ierr = MatShift(test_mat,1.);CHKERRQ(ierr);
327674ae819SStefano Zampini   ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr);
328674ae819SStefano Zampini   ierr = MatDestroy(&test_mat);CHKERRQ(ierr);
329c7017625SStefano Zampini   if (isdir) {
330c7017625SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
331c7017625SStefano Zampini   } else {
332c7017625SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann nullspace matrices: %1.14e\n",PetscGlobalRank,test_err);CHKERRQ(ierr);
333c7017625SStefano Zampini   }
334674ae819SStefano Zampini 
335674ae819SStefano Zampini   /* Create ksp object suitable for extreme eigenvalues' estimation */
336674ae819SStefano Zampini   ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
33715579a77SStefano Zampini   ierr = PetscObjectIncrementTabLevel((PetscObject)check_ksp,(PetscObject)local_ksp,1);CHKERRQ(ierr);
338bd5e1169SStefano Zampini   ierr = KSPGetOptionsPrefix(local_ksp,&prefix);CHKERRQ(ierr);
339bd5e1169SStefano Zampini   ierr = KSPSetOptionsPrefix(check_ksp,prefix);CHKERRQ(ierr);
340bd5e1169SStefano Zampini   ierr = KSPAppendOptionsPrefix(check_ksp,"approxcheck_");CHKERRQ(ierr);
341422a814eSBarry Smith   ierr = KSPSetErrorIfNotConverged(check_ksp,pc->erroriffailure);CHKERRQ(ierr);
34223ee1639SBarry Smith   ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr);
343bd5e1169SStefano Zampini   ierr = KSPSetTolerances(check_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
344674ae819SStefano Zampini   ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
345bd5e1169SStefano Zampini   ierr = KSPSetFromOptions(check_ksp);CHKERRQ(ierr);
346674ae819SStefano Zampini   ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr);
347674ae819SStefano Zampini   ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
348674ae819SStefano Zampini   ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr);
349674ae819SStefano Zampini   ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
350674ae819SStefano Zampini   ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr);
351c0decd05SBarry Smith   ierr = KSPCheckSolve(check_ksp,pc,work2);CHKERRQ(ierr);
352c7017625SStefano Zampini   ierr = VecAXPY(work2,-1.,work1);CHKERRQ(ierr);
353674ae819SStefano Zampini   ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr);
354674ae819SStefano Zampini   ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
355674ae819SStefano Zampini   ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
356c7017625SStefano Zampini   if (isdir) {
357c7017625SStefano 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);
358c7017625SStefano Zampini   } else {
359c7017625SStefano 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);
360c7017625SStefano Zampini   }
361674ae819SStefano Zampini   ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
362674ae819SStefano Zampini   ierr = VecDestroy(&work1);CHKERRQ(ierr);
363674ae819SStefano Zampini   ierr = VecDestroy(&work2);CHKERRQ(ierr);
364674ae819SStefano Zampini   ierr = VecDestroy(&work3);CHKERRQ(ierr);
365674ae819SStefano Zampini   PetscFunctionReturn(0);
366674ae819SStefano Zampini }
367