xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 674ae81933ad1188dbe7c87a2be01dd8bb4076e0)
1*674ae819SStefano Zampini #include "bddc.h"
2*674ae819SStefano Zampini #include "bddcprivate.h"
3*674ae819SStefano Zampini 
4*674ae819SStefano Zampini #undef __FUNCT__
5*674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCoarse"
6*674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCoarse(PC pc, MatNullSpace* CoarseNullSpace)
7*674ae819SStefano Zampini {
8*674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
9*674ae819SStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
10*674ae819SStefano Zampini   MatNullSpace   tempCoarseNullSpace;
11*674ae819SStefano Zampini   const Vec      *nsp_vecs;
12*674ae819SStefano Zampini   Vec            *coarse_nsp_vecs,local_vec,local_primal_vec;
13*674ae819SStefano Zampini   PetscInt       nsp_size,coarse_nsp_size,i;
14*674ae819SStefano Zampini   PetscBool      nsp_has_cnst;
15*674ae819SStefano Zampini   PetscReal      test_null;
16*674ae819SStefano Zampini   PetscErrorCode ierr;
17*674ae819SStefano Zampini 
18*674ae819SStefano Zampini   PetscFunctionBegin;
19*674ae819SStefano Zampini   tempCoarseNullSpace = 0;
20*674ae819SStefano Zampini   coarse_nsp_size = 0;
21*674ae819SStefano Zampini   coarse_nsp_vecs = 0;
22*674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr);
23*674ae819SStefano Zampini   if (pcbddc->coarse_mat) {
24*674ae819SStefano Zampini     ierr = PetscMalloc((nsp_size+1)*sizeof(Vec),&coarse_nsp_vecs);CHKERRQ(ierr);
25*674ae819SStefano Zampini     for (i=0;i<nsp_size+1;i++) {
26*674ae819SStefano Zampini       ierr = VecDuplicate(pcbddc->coarse_vec,&coarse_nsp_vecs[i]);CHKERRQ(ierr);
27*674ae819SStefano Zampini     }
28*674ae819SStefano Zampini   }
29*674ae819SStefano Zampini   ierr = MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);CHKERRQ(ierr);
30*674ae819SStefano Zampini   if (nsp_has_cnst) {
31*674ae819SStefano Zampini     ierr = VecSet(local_vec,1.0);CHKERRQ(ierr);
32*674ae819SStefano Zampini     ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr);
33*674ae819SStefano Zampini     ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34*674ae819SStefano Zampini     ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
35*674ae819SStefano Zampini     if (pcbddc->coarse_mat) {
36*674ae819SStefano Zampini       if (pcbddc->dbg_flag ) {
37*674ae819SStefano Zampini         ierr = MatMult(pcbddc->coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr);
38*674ae819SStefano Zampini         ierr = VecNorm(pcbddc->coarse_rhs,NORM_INFINITY,&test_null);CHKERRQ(ierr);
39*674ae819SStefano Zampini         ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);CHKERRQ(ierr);
40*674ae819SStefano Zampini       }
41*674ae819SStefano Zampini       ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr);
42*674ae819SStefano Zampini       coarse_nsp_size++;
43*674ae819SStefano Zampini     }
44*674ae819SStefano Zampini   }
45*674ae819SStefano Zampini   for (i=0;i<nsp_size;i++)  {
46*674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
47*674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48*674ae819SStefano Zampini     ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr);
49*674ae819SStefano Zampini     ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
50*674ae819SStefano Zampini     ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
51*674ae819SStefano Zampini     if (pcbddc->coarse_mat) {
52*674ae819SStefano Zampini       if (pcbddc->dbg_flag ) {
53*674ae819SStefano Zampini         ierr = MatMult(pcbddc->coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr);
54*674ae819SStefano Zampini         ierr = VecNorm(pcbddc->coarse_rhs,NORM_2,&test_null);CHKERRQ(ierr);
55*674ae819SStefano Zampini         ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);CHKERRQ(ierr);
56*674ae819SStefano Zampini       }
57*674ae819SStefano Zampini       ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr);
58*674ae819SStefano Zampini       coarse_nsp_size++;
59*674ae819SStefano Zampini     }
60*674ae819SStefano Zampini   }
61*674ae819SStefano Zampini   if (coarse_nsp_size > 0) {
62*674ae819SStefano Zampini     /* TODO orthonormalize vecs when new_nsp_size > 0! */
63*674ae819SStefano Zampini     MPI_Comm coarse_comm;
64*674ae819SStefano Zampini     ierr = PetscObjectGetComm((PetscObject)(pcbddc->coarse_mat),&coarse_comm);CHKERRQ(ierr);
65*674ae819SStefano Zampini     ierr = VecNormalize(coarse_nsp_vecs[0],NULL);CHKERRQ(ierr);
66*674ae819SStefano Zampini     ierr = MatNullSpaceCreate(coarse_comm,PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);CHKERRQ(ierr);
67*674ae819SStefano Zampini     for (i=0;i<nsp_size+1;i++) {
68*674ae819SStefano Zampini       ierr = VecDestroy(&coarse_nsp_vecs[i]);CHKERRQ(ierr);
69*674ae819SStefano Zampini     }
70*674ae819SStefano Zampini   }
71*674ae819SStefano Zampini   ierr = PetscFree(coarse_nsp_vecs);CHKERRQ(ierr);
72*674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
73*674ae819SStefano Zampini   ierr = VecDestroy(&local_primal_vec);CHKERRQ(ierr);
74*674ae819SStefano Zampini   *CoarseNullSpace = tempCoarseNullSpace;
75*674ae819SStefano Zampini   PetscFunctionReturn(0);
76*674ae819SStefano Zampini }
77*674ae819SStefano Zampini 
78*674ae819SStefano Zampini 
79*674ae819SStefano Zampini #undef __FUNCT__
80*674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC"
81*674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y)
82*674ae819SStefano Zampini {
83*674ae819SStefano Zampini   NullSpaceCorrection_ctx pc_ctx;
84*674ae819SStefano Zampini   PetscErrorCode          ierr;
85*674ae819SStefano Zampini 
86*674ae819SStefano Zampini   PetscFunctionBegin;
87*674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
88*674ae819SStefano Zampini   /* E */
89*674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr);
90*674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr);
91*674ae819SStefano Zampini   /* P^-1 */
92*674ae819SStefano Zampini   ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr);
93*674ae819SStefano Zampini   /* E^T */
94*674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr);
95*674ae819SStefano Zampini   ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr);
96*674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr);
97*674ae819SStefano Zampini   /* Sum contributions */
98*674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr);
99*674ae819SStefano Zampini   PetscFunctionReturn(0);
100*674ae819SStefano Zampini }
101*674ae819SStefano Zampini 
102*674ae819SStefano Zampini #undef __FUNCT__
103*674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC"
104*674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc)
105*674ae819SStefano Zampini {
106*674ae819SStefano Zampini   NullSpaceCorrection_ctx pc_ctx;
107*674ae819SStefano Zampini   PetscErrorCode          ierr;
108*674ae819SStefano Zampini 
109*674ae819SStefano Zampini   PetscFunctionBegin;
110*674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
111*674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr);
112*674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr);
113*674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr);
114*674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr);
115*674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr);
116*674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr);
117*674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr);
118*674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr);
119*674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
120*674ae819SStefano Zampini   PetscFunctionReturn(0);
121*674ae819SStefano Zampini }
122*674ae819SStefano Zampini 
123*674ae819SStefano Zampini /*PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec);
124*674ae819SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC);*/
125*674ae819SStefano Zampini 
126*674ae819SStefano Zampini #undef __FUNCT__
127*674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection"
128*674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs)
129*674ae819SStefano Zampini {
130*674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
131*674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
132*674ae819SStefano Zampini   Mat_IS*        matis = (Mat_IS*)pc->pmat->data;
133*674ae819SStefano Zampini   KSP            *local_ksp;
134*674ae819SStefano Zampini   PC             newpc;
135*674ae819SStefano Zampini   NullSpaceCorrection_ctx  shell_ctx;
136*674ae819SStefano Zampini   Mat            local_mat,local_pmat,small_mat,inv_small_mat;
137*674ae819SStefano Zampini   MatStructure   local_mat_struct;
138*674ae819SStefano Zampini   Vec            work1,work2;
139*674ae819SStefano Zampini   const Vec      *nullvecs;
140*674ae819SStefano Zampini   VecScatter     scatter_ctx;
141*674ae819SStefano Zampini   IS             is_aux;
142*674ae819SStefano Zampini   MatFactorInfo  matinfo;
143*674ae819SStefano Zampini   PetscScalar    *basis_mat,*Kbasis_mat,*array,*array_mat;
144*674ae819SStefano Zampini   PetscScalar    one = 1.0,zero = 0.0, m_one = -1.0;
145*674ae819SStefano Zampini   PetscInt       basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R;
146*674ae819SStefano Zampini   PetscBool      nnsp_has_cnst;
147*674ae819SStefano Zampini   PetscErrorCode ierr;
148*674ae819SStefano Zampini 
149*674ae819SStefano Zampini   PetscFunctionBegin;
150*674ae819SStefano Zampini   /* Infer the local solver */
151*674ae819SStefano Zampini   ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr);
152*674ae819SStefano Zampini   ierr = VecGetSize(pcis->vec1_D,&n_I);CHKERRQ(ierr);
153*674ae819SStefano Zampini   ierr = VecGetSize(pcbddc->vec1_R,&n_R);CHKERRQ(ierr);
154*674ae819SStefano Zampini   if (basis_dofs == n_I) {
155*674ae819SStefano Zampini     /* Dirichlet solver */
156*674ae819SStefano Zampini     local_ksp = &pcbddc->ksp_D;
157*674ae819SStefano Zampini   } else if (basis_dofs == n_R) {
158*674ae819SStefano Zampini     /* Neumann solver */
159*674ae819SStefano Zampini     local_ksp = &pcbddc->ksp_R;
160*674ae819SStefano Zampini   } else {
161*674ae819SStefano Zampini     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in %s: unknown local IS size %d. n_I=%d, n_R=%d)\n",__FUNCT__,basis_dofs,n_I,n_R);
162*674ae819SStefano Zampini   }
163*674ae819SStefano Zampini   ierr = KSPGetOperators(*local_ksp,&local_mat,&local_pmat,&local_mat_struct);CHKERRQ(ierr);
164*674ae819SStefano Zampini 
165*674ae819SStefano Zampini   /* Get null space vecs */
166*674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr);
167*674ae819SStefano Zampini   basis_size = nnsp_size;
168*674ae819SStefano Zampini   if (nnsp_has_cnst) {
169*674ae819SStefano Zampini     basis_size++;
170*674ae819SStefano Zampini   }
171*674ae819SStefano Zampini 
172*674ae819SStefano Zampini   /* Create shell ctx */
173*674ae819SStefano Zampini   ierr = PetscMalloc(sizeof(*shell_ctx),&shell_ctx);CHKERRQ(ierr);
174*674ae819SStefano Zampini 
175*674ae819SStefano Zampini   /* Create work vectors in shell context */
176*674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr);
177*674ae819SStefano Zampini   ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr);
178*674ae819SStefano Zampini   ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr);
179*674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr);
180*674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr);
181*674ae819SStefano Zampini   ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr);
182*674ae819SStefano Zampini   ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr);
183*674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr);
184*674ae819SStefano Zampini 
185*674ae819SStefano Zampini   /* Allocate workspace */
186*674ae819SStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );CHKERRQ(ierr);
187*674ae819SStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr);
188*674ae819SStefano Zampini   ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
189*674ae819SStefano Zampini   ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
190*674ae819SStefano Zampini 
191*674ae819SStefano Zampini   /* Restrict local null space on selected dofs (Dirichlet or Neumann)
192*674ae819SStefano Zampini      and compute matrices N and K*N */
193*674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr);
194*674ae819SStefano Zampini   ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr);
195*674ae819SStefano Zampini   ierr = VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);CHKERRQ(ierr);
196*674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
197*674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
198*674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
199*674ae819SStefano Zampini     ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
200*674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
201*674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
202*674ae819SStefano Zampini     ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
203*674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
204*674ae819SStefano Zampini     ierr = VecResetArray(work1);CHKERRQ(ierr);
205*674ae819SStefano Zampini     ierr = VecResetArray(work2);CHKERRQ(ierr);
206*674ae819SStefano Zampini   }
207*674ae819SStefano Zampini   if (nnsp_has_cnst) {
208*674ae819SStefano Zampini     ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
209*674ae819SStefano Zampini     ierr = VecSet(work1,one);CHKERRQ(ierr);
210*674ae819SStefano Zampini     ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
211*674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
212*674ae819SStefano Zampini     ierr = VecResetArray(work1);CHKERRQ(ierr);
213*674ae819SStefano Zampini     ierr = VecResetArray(work2);CHKERRQ(ierr);
214*674ae819SStefano Zampini   }
215*674ae819SStefano Zampini   ierr = VecDestroy(&work1);CHKERRQ(ierr);
216*674ae819SStefano Zampini   ierr = VecDestroy(&work2);CHKERRQ(ierr);
217*674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
218*674ae819SStefano Zampini   ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
219*674ae819SStefano Zampini   ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
220*674ae819SStefano Zampini 
221*674ae819SStefano Zampini   /* Assemble another Mat object in shell context */
222*674ae819SStefano Zampini   ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr);
223*674ae819SStefano Zampini   ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr);
224*674ae819SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr);
225*674ae819SStefano Zampini   ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr);
226*674ae819SStefano Zampini   ierr = ISDestroy(&is_aux);CHKERRQ(ierr);
227*674ae819SStefano Zampini   ierr = PetscMalloc(basis_size*basis_size*sizeof(PetscScalar),&array_mat);CHKERRQ(ierr);
228*674ae819SStefano Zampini   for (k=0;k<basis_size;k++) {
229*674ae819SStefano Zampini     ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr);
230*674ae819SStefano Zampini     ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr);
231*674ae819SStefano Zampini     ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr);
232*674ae819SStefano Zampini     ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr);
233*674ae819SStefano Zampini     ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr);
234*674ae819SStefano Zampini     ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
235*674ae819SStefano Zampini     for (i=0;i<basis_size;i++) {
236*674ae819SStefano Zampini       array_mat[i*basis_size+k]=array[i];
237*674ae819SStefano Zampini     }
238*674ae819SStefano Zampini     ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
239*674ae819SStefano Zampini   }
240*674ae819SStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr);
241*674ae819SStefano Zampini   ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr);
242*674ae819SStefano Zampini   ierr = PetscFree(array_mat);CHKERRQ(ierr);
243*674ae819SStefano Zampini   ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr);
244*674ae819SStefano Zampini   ierr = MatDestroy(&small_mat);CHKERRQ(ierr);
245*674ae819SStefano Zampini   ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr);
246*674ae819SStefano Zampini 
247*674ae819SStefano Zampini   /* Rebuild local PC */
248*674ae819SStefano Zampini   ierr = KSPGetPC(*local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr);
249*674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr);
250*674ae819SStefano Zampini   ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr);
251*674ae819SStefano Zampini   ierr = PCSetOperators(newpc,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
252*674ae819SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
253*674ae819SStefano Zampini   ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr);
254*674ae819SStefano Zampini   ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr);
255*674ae819SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr);
256*674ae819SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
257*674ae819SStefano Zampini   ierr = KSPSetPC(*local_ksp,newpc);CHKERRQ(ierr);
258*674ae819SStefano Zampini   ierr = PCDestroy(&newpc);CHKERRQ(ierr);
259*674ae819SStefano Zampini   ierr = KSPSetUp(*local_ksp);CHKERRQ(ierr);
260*674ae819SStefano Zampini   /* test */
261*674ae819SStefano Zampini   /* TODO: this cause a deadlock when doing multilevel */
262*674ae819SStefano Zampini #if 0
263*674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
264*674ae819SStefano Zampini     KSP         check_ksp;
265*674ae819SStefano Zampini     PC          check_pc;
266*674ae819SStefano Zampini     Mat         test_mat;
267*674ae819SStefano Zampini     Vec         work3;
268*674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
269*674ae819SStefano Zampini     PetscReal   test_err,lambda_min,lambda_max;
270*674ae819SStefano Zampini     PetscBool   setsym,issym=PETSC_FALSE;
271*674ae819SStefano Zampini 
272*674ae819SStefano Zampini     ierr = KSPGetPC(*local_ksp,&check_pc);CHKERRQ(ierr);
273*674ae819SStefano Zampini     ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr);
274*674ae819SStefano Zampini     ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr);
275*674ae819SStefano Zampini     ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr);
276*674ae819SStefano Zampini     ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr);
277*674ae819SStefano Zampini     ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr);
278*674ae819SStefano Zampini     ierr = VecCopy(work1,work2);CHKERRQ(ierr);
279*674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
280*674ae819SStefano Zampini     ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr);
281*674ae819SStefano Zampini     ierr = VecAXPY(work1,m_one,work2);CHKERRQ(ierr);
282*674ae819SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
283*674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank);
284*674ae819SStefano Zampini     if (basis_dofs == n_I) {
285*674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Dirichlet ");
286*674ae819SStefano Zampini     } else {
287*674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Neumann ");
288*674ae819SStefano Zampini     }
289*674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"solver is :%1.14e\n",test_err);
290*674ae819SStefano Zampini 
291*674ae819SStefano Zampini     ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr);
292*674ae819SStefano Zampini     ierr = MatShift(test_mat,one);CHKERRQ(ierr);
293*674ae819SStefano Zampini     ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr);
294*674ae819SStefano Zampini     ierr = MatDestroy(&test_mat);CHKERRQ(ierr);
295*674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err);
296*674ae819SStefano Zampini 
297*674ae819SStefano Zampini     /* Create ksp object suitable for extreme eigenvalues' estimation */
298*674ae819SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
299*674ae819SStefano Zampini     ierr = KSPSetOperators(check_ksp,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
300*674ae819SStefano Zampini     ierr = KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr);
301*674ae819SStefano Zampini     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
302*674ae819SStefano Zampini     ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr);
303*674ae819SStefano Zampini     if (issym) {
304*674ae819SStefano Zampini       ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
305*674ae819SStefano Zampini     }
306*674ae819SStefano Zampini     ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr);
307*674ae819SStefano Zampini     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
308*674ae819SStefano Zampini     ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr);
309*674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
310*674ae819SStefano Zampini     ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr);
311*674ae819SStefano Zampini     ierr = VecAXPY(work2,m_one,work1);CHKERRQ(ierr);
312*674ae819SStefano Zampini     ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr);
313*674ae819SStefano Zampini     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
314*674ae819SStefano Zampini     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
315*674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d error for adapted KSP %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max);
316*674ae819SStefano Zampini     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
317*674ae819SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
318*674ae819SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
319*674ae819SStefano Zampini     ierr = VecDestroy(&work3);CHKERRQ(ierr);
320*674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
321*674ae819SStefano Zampini   }
322*674ae819SStefano Zampini #endif
323*674ae819SStefano Zampini   PetscFunctionReturn(0);
324*674ae819SStefano Zampini }
325*674ae819SStefano Zampini 
326*674ae819SStefano Zampini #undef __FUNCT__
327*674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAdaptGlobal"
328*674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc)
329*674ae819SStefano Zampini {
330*674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)  (pc->data);
331*674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
332*674ae819SStefano Zampini   KSP            inv_change;
333*674ae819SStefano Zampini   PC             pc_change;
334*674ae819SStefano Zampini   const Vec      *nsp_vecs;
335*674ae819SStefano Zampini   Vec            *new_nsp_vecs;
336*674ae819SStefano Zampini   PetscInt       i,nsp_size,new_nsp_size,start_new;
337*674ae819SStefano Zampini   PetscBool      nsp_has_cnst;
338*674ae819SStefano Zampini   MatNullSpace   new_nsp;
339*674ae819SStefano Zampini   PetscErrorCode ierr;
340*674ae819SStefano Zampini   MPI_Comm       comm;
341*674ae819SStefano Zampini 
342*674ae819SStefano Zampini   PetscFunctionBegin;
343*674ae819SStefano Zampini   /* create KSP for change of basis */
344*674ae819SStefano Zampini   ierr = KSPCreate(PETSC_COMM_SELF,&inv_change);CHKERRQ(ierr);
345*674ae819SStefano Zampini   ierr = KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix,SAME_PRECONDITIONER);CHKERRQ(ierr);
346*674ae819SStefano Zampini   ierr = KSPSetType(inv_change,KSPPREONLY);CHKERRQ(ierr);
347*674ae819SStefano Zampini   ierr = KSPGetPC(inv_change,&pc_change);CHKERRQ(ierr);
348*674ae819SStefano Zampini   ierr = PCSetType(pc_change,PCLU);CHKERRQ(ierr);
349*674ae819SStefano Zampini   ierr = KSPSetUp(inv_change);CHKERRQ(ierr);
350*674ae819SStefano Zampini   /* get nullspace and transform it */
351*674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr);
352*674ae819SStefano Zampini   new_nsp_size = nsp_size;
353*674ae819SStefano Zampini   if (nsp_has_cnst) {
354*674ae819SStefano Zampini     new_nsp_size++;
355*674ae819SStefano Zampini   }
356*674ae819SStefano Zampini   ierr = PetscMalloc(new_nsp_size*sizeof(Vec),&new_nsp_vecs);CHKERRQ(ierr);
357*674ae819SStefano Zampini   for (i=0;i<new_nsp_size;i++) {
358*674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&new_nsp_vecs[i]);CHKERRQ(ierr);
359*674ae819SStefano Zampini   }
360*674ae819SStefano Zampini   start_new = 0;
361*674ae819SStefano Zampini   if (nsp_has_cnst) {
362*674ae819SStefano Zampini     start_new = 1;
363*674ae819SStefano Zampini     ierr = VecSet(new_nsp_vecs[0],1.0);CHKERRQ(ierr);
364*674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr);
365*674ae819SStefano Zampini     ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B);
366*674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
367*674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
368*674ae819SStefano Zampini   }
369*674ae819SStefano Zampini   for (i=0;i<nsp_size;i++) {
370*674ae819SStefano Zampini     ierr = VecCopy(nsp_vecs[i],new_nsp_vecs[i+start_new]);CHKERRQ(ierr);
371*674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
372*674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
373*674ae819SStefano Zampini     ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B);
374*674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
375*674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
376*674ae819SStefano Zampini   }
377*674ae819SStefano Zampini   ierr = VecNormalize(new_nsp_vecs[0],NULL);CHKERRQ(ierr);
378*674ae819SStefano Zampini   /* TODO : Orthonormalize vecs when new_nsp_size > 0! */
379*674ae819SStefano Zampini #if 0
380*674ae819SStefano Zampini   PetscBool nsp_t=PETSC_FALSE;
381*674ae819SStefano Zampini   ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr);
382*674ae819SStefano Zampini   printf("Original Null Space test: %d\n",nsp_t);
383*674ae819SStefano Zampini   Mat temp_mat;
384*674ae819SStefano Zampini   Mat_IS* matis = (Mat_IS*)pc->pmat->data;
385*674ae819SStefano Zampini     temp_mat = matis->A;
386*674ae819SStefano Zampini     matis->A = pcbddc->local_mat;
387*674ae819SStefano Zampini     pcbddc->local_mat = temp_mat;
388*674ae819SStefano Zampini   ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr);
389*674ae819SStefano Zampini   printf("Original Null Space, mat changed test: %d\n",nsp_t);
390*674ae819SStefano Zampini   {
391*674ae819SStefano Zampini     PetscReal test_norm;
392*674ae819SStefano Zampini     for (i=0;i<new_nsp_size;i++) {
393*674ae819SStefano Zampini       ierr = MatMult(pc->pmat,new_nsp_vecs[i],pcis->vec1_global);CHKERRQ(ierr);
394*674ae819SStefano Zampini       ierr = VecNorm(pcis->vec1_global,NORM_2,&test_norm);CHKERRQ(ierr);
395*674ae819SStefano Zampini       if (test_norm > 1.e-12) {
396*674ae819SStefano Zampini         printf("------------ERROR VEC %d------------------\n",i);
397*674ae819SStefano Zampini         ierr = VecView(pcis->vec1_global,PETSC_VIEWER_STDOUT_WORLD);
398*674ae819SStefano Zampini         printf("------------------------------------------\n");
399*674ae819SStefano Zampini       }
400*674ae819SStefano Zampini     }
401*674ae819SStefano Zampini   }
402*674ae819SStefano Zampini #endif
403*674ae819SStefano Zampini 
404*674ae819SStefano Zampini   ierr = KSPDestroy(&inv_change);CHKERRQ(ierr);
405*674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
406*674ae819SStefano Zampini   ierr = MatNullSpaceCreate(comm,PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);CHKERRQ(ierr);
407*674ae819SStefano Zampini   ierr = PCBDDCSetNullSpace(pc,new_nsp);CHKERRQ(ierr);
408*674ae819SStefano Zampini   ierr = MatNullSpaceDestroy(&new_nsp);CHKERRQ(ierr);
409*674ae819SStefano Zampini #if 0
410*674ae819SStefano Zampini   ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr);
411*674ae819SStefano Zampini   printf("New Null Space, mat changed: %d\n",nsp_t);
412*674ae819SStefano Zampini     temp_mat = matis->A;
413*674ae819SStefano Zampini     matis->A = pcbddc->local_mat;
414*674ae819SStefano Zampini     pcbddc->local_mat = temp_mat;
415*674ae819SStefano Zampini   ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr);
416*674ae819SStefano Zampini   printf("New Null Space, mat original: %d\n",nsp_t);
417*674ae819SStefano Zampini #endif
418*674ae819SStefano Zampini 
419*674ae819SStefano Zampini   for (i=0;i<new_nsp_size;i++) {
420*674ae819SStefano Zampini     ierr = VecDestroy(&new_nsp_vecs[i]);CHKERRQ(ierr);
421*674ae819SStefano Zampini   }
422*674ae819SStefano Zampini   ierr = PetscFree(new_nsp_vecs);CHKERRQ(ierr);
423*674ae819SStefano Zampini   PetscFunctionReturn(0);
424*674ae819SStefano Zampini }
425