xref: /petsc/src/ksp/pc/impls/bddc/bddcnullspace.c (revision 906d46d482912ff65fe95caf6add84f15a3f7e20)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3674ae819SStefano Zampini 
4674ae819SStefano Zampini #undef __FUNCT__
5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCoarse"
6fdc635d7SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCoarse(PC pc, Mat coarse_mat, MatNullSpace* CoarseNullSpace)
7674ae819SStefano Zampini {
8674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
9674ae819SStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
1098a51de6SStefano Zampini   MatNullSpace   tempCoarseNullSpace=NULL;
11674ae819SStefano Zampini   const Vec      *nsp_vecs;
1298a51de6SStefano Zampini   Vec            *coarse_nsp_vecs,local_vec,local_primal_vec,wcoarse_vec,wcoarse_rhs;
13674ae819SStefano Zampini   PetscInt       nsp_size,coarse_nsp_size,i;
14674ae819SStefano Zampini   PetscBool      nsp_has_cnst;
15674ae819SStefano Zampini   PetscReal      test_null;
16674ae819SStefano Zampini   PetscErrorCode ierr;
17674ae819SStefano Zampini 
18674ae819SStefano Zampini   PetscFunctionBegin;
19674ae819SStefano Zampini   tempCoarseNullSpace = 0;
20674ae819SStefano Zampini   coarse_nsp_size = 0;
21674ae819SStefano Zampini   coarse_nsp_vecs = 0;
22674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr);
23fdc635d7SStefano Zampini   if (coarse_mat) {
24785e854fSJed Brown     ierr = PetscMalloc1((nsp_size+1),&coarse_nsp_vecs);CHKERRQ(ierr);
25674ae819SStefano Zampini     for (i=0;i<nsp_size+1;i++) {
2698a51de6SStefano Zampini       ierr = MatGetVecs(coarse_mat,&coarse_nsp_vecs[i],NULL);CHKERRQ(ierr);
2798a51de6SStefano Zampini     }
2898a51de6SStefano Zampini     if (pcbddc->dbg_flag) {
2998a51de6SStefano Zampini       ierr = MatGetVecs(coarse_mat,&wcoarse_vec,&wcoarse_rhs);CHKERRQ(ierr);
30674ae819SStefano Zampini     }
31674ae819SStefano Zampini   }
32674ae819SStefano Zampini   ierr = MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);CHKERRQ(ierr);
33674ae819SStefano Zampini   if (nsp_has_cnst) {
34674ae819SStefano Zampini     ierr = VecSet(local_vec,1.0);CHKERRQ(ierr);
35674ae819SStefano Zampini     ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr);
3612edc857SStefano Zampini     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3712edc857SStefano Zampini     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
38fdc635d7SStefano Zampini     if (coarse_mat) {
3998a51de6SStefano Zampini       PetscScalar *array_out;
4098a51de6SStefano Zampini       const PetscScalar *array_in;
4198a51de6SStefano Zampini       PetscInt lsize;
42674ae819SStefano Zampini       if (pcbddc->dbg_flag) {
4398a51de6SStefano Zampini         PetscViewer dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)coarse_mat));
4498a51de6SStefano Zampini         ierr = MatMult(coarse_mat,wcoarse_vec,wcoarse_rhs);CHKERRQ(ierr);
4598a51de6SStefano Zampini         ierr = VecNorm(wcoarse_rhs,NORM_INFINITY,&test_null);CHKERRQ(ierr);
4698a51de6SStefano Zampini         ierr = PetscViewerASCIIPrintf(dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);CHKERRQ(ierr);
4798a51de6SStefano Zampini         ierr = PetscViewerFlush(dbg_viewer);CHKERRQ(ierr);
48674ae819SStefano Zampini       }
4998a51de6SStefano Zampini       ierr = VecGetLocalSize(pcbddc->coarse_vec,&lsize);CHKERRQ(ierr);
5098a51de6SStefano Zampini       ierr = VecGetArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr);
5198a51de6SStefano Zampini       ierr = VecGetArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr);
5298a51de6SStefano Zampini       ierr = PetscMemcpy(array_out,array_in,lsize*sizeof(PetscScalar));CHKERRQ(ierr);
5398a51de6SStefano Zampini       ierr = VecRestoreArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr);
5498a51de6SStefano Zampini       ierr = VecRestoreArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr);
55674ae819SStefano Zampini       coarse_nsp_size++;
56674ae819SStefano Zampini     }
57674ae819SStefano Zampini   }
58674ae819SStefano Zampini   for (i=0;i<nsp_size;i++)  {
59674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
60674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
61674ae819SStefano Zampini     ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr);
6212edc857SStefano Zampini     ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6312edc857SStefano Zampini     ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64fdc635d7SStefano Zampini     if (coarse_mat) {
6598a51de6SStefano Zampini       PetscScalar *array_out;
6698a51de6SStefano Zampini       const PetscScalar *array_in;
6798a51de6SStefano Zampini       PetscInt lsize;
68674ae819SStefano Zampini       if (pcbddc->dbg_flag) {
6998a51de6SStefano Zampini         PetscViewer dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)coarse_mat));
7098a51de6SStefano Zampini         ierr = MatMult(coarse_mat,wcoarse_vec,wcoarse_rhs);CHKERRQ(ierr);
7198a51de6SStefano Zampini         ierr = VecNorm(wcoarse_rhs,NORM_2,&test_null);CHKERRQ(ierr);
7298a51de6SStefano Zampini         ierr = PetscViewerASCIIPrintf(dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);CHKERRQ(ierr);
7398a51de6SStefano Zampini         ierr = PetscViewerFlush(dbg_viewer);CHKERRQ(ierr);
74674ae819SStefano Zampini       }
7598a51de6SStefano Zampini       ierr = VecGetLocalSize(pcbddc->coarse_vec,&lsize);CHKERRQ(ierr);
7698a51de6SStefano Zampini       ierr = VecGetArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr);
7798a51de6SStefano Zampini       ierr = VecGetArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr);
7898a51de6SStefano Zampini       ierr = PetscMemcpy(array_out,array_in,lsize*sizeof(PetscScalar));CHKERRQ(ierr);
7998a51de6SStefano Zampini       ierr = VecRestoreArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr);
8098a51de6SStefano Zampini       ierr = VecRestoreArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr);
81674ae819SStefano Zampini       coarse_nsp_size++;
82674ae819SStefano Zampini     }
83674ae819SStefano Zampini   }
84674ae819SStefano Zampini   if (coarse_nsp_size > 0) {
859a7d3425SStefano Zampini     ierr = PCBDDCOrthonormalizeVecs(coarse_nsp_size,coarse_nsp_vecs);CHKERRQ(ierr);
86fdc635d7SStefano Zampini     ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)coarse_mat),PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);CHKERRQ(ierr);
87674ae819SStefano Zampini     for (i=0;i<nsp_size+1;i++) {
88674ae819SStefano Zampini       ierr = VecDestroy(&coarse_nsp_vecs[i]);CHKERRQ(ierr);
89674ae819SStefano Zampini     }
90674ae819SStefano Zampini   }
9198a51de6SStefano Zampini   if (coarse_mat) {
92674ae819SStefano Zampini     ierr = PetscFree(coarse_nsp_vecs);CHKERRQ(ierr);
9398a51de6SStefano Zampini     if (pcbddc->dbg_flag) {
9498a51de6SStefano Zampini       ierr = VecDestroy(&wcoarse_vec);CHKERRQ(ierr);
9598a51de6SStefano Zampini       ierr = VecDestroy(&wcoarse_rhs);CHKERRQ(ierr);
9698a51de6SStefano Zampini     }
9798a51de6SStefano Zampini   }
98674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
99674ae819SStefano Zampini   ierr = VecDestroy(&local_primal_vec);CHKERRQ(ierr);
100674ae819SStefano Zampini   *CoarseNullSpace = tempCoarseNullSpace;
101674ae819SStefano Zampini   PetscFunctionReturn(0);
102674ae819SStefano Zampini }
103674ae819SStefano Zampini 
104674ae819SStefano Zampini 
105674ae819SStefano Zampini #undef __FUNCT__
106674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC"
107674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y)
108674ae819SStefano Zampini {
109674ae819SStefano Zampini   NullSpaceCorrection_ctx pc_ctx;
110674ae819SStefano Zampini   PetscErrorCode          ierr;
111674ae819SStefano Zampini 
112674ae819SStefano Zampini   PetscFunctionBegin;
113674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
114674ae819SStefano Zampini   /* E */
115674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr);
116674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr);
117674ae819SStefano Zampini   /* P^-1 */
118674ae819SStefano Zampini   ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr);
119674ae819SStefano Zampini   /* E^T */
120674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr);
121674ae819SStefano Zampini   ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr);
122674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr);
123674ae819SStefano Zampini   /* Sum contributions */
124674ae819SStefano Zampini   ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr);
125674ae819SStefano Zampini   PetscFunctionReturn(0);
126674ae819SStefano Zampini }
127674ae819SStefano Zampini 
128674ae819SStefano Zampini #undef __FUNCT__
129674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC"
130674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc)
131674ae819SStefano Zampini {
132674ae819SStefano Zampini   NullSpaceCorrection_ctx pc_ctx;
133674ae819SStefano Zampini   PetscErrorCode          ierr;
134674ae819SStefano Zampini 
135674ae819SStefano Zampini   PetscFunctionBegin;
136674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
137674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr);
138674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr);
139674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr);
140674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr);
141674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr);
142674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr);
143674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr);
144674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr);
145674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
146674ae819SStefano Zampini   PetscFunctionReturn(0);
147674ae819SStefano Zampini }
148674ae819SStefano Zampini 
14928509bceSStefano Zampini /*
15028509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec);
15128509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC);
15228509bceSStefano Zampini */
153674ae819SStefano Zampini 
154674ae819SStefano Zampini #undef __FUNCT__
155674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection"
156674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs)
157674ae819SStefano Zampini {
158674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
159674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
160674ae819SStefano Zampini   Mat_IS*        matis = (Mat_IS*)pc->pmat->data;
161674ae819SStefano Zampini   KSP            *local_ksp;
162674ae819SStefano Zampini   PC             newpc;
163674ae819SStefano Zampini   NullSpaceCorrection_ctx  shell_ctx;
164674ae819SStefano Zampini   Mat            local_mat,local_pmat,small_mat,inv_small_mat;
165674ae819SStefano Zampini   Vec            work1,work2;
166674ae819SStefano Zampini   const Vec      *nullvecs;
167674ae819SStefano Zampini   VecScatter     scatter_ctx;
168674ae819SStefano Zampini   IS             is_aux;
169674ae819SStefano Zampini   MatFactorInfo  matinfo;
170674ae819SStefano Zampini   PetscScalar    *basis_mat,*Kbasis_mat,*array,*array_mat;
171674ae819SStefano Zampini   PetscScalar    one = 1.0,zero = 0.0, m_one = -1.0;
172674ae819SStefano Zampini   PetscInt       basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R;
173674ae819SStefano Zampini   PetscBool      nnsp_has_cnst;
174674ae819SStefano Zampini   PetscErrorCode ierr;
175674ae819SStefano Zampini 
176674ae819SStefano Zampini   PetscFunctionBegin;
177674ae819SStefano Zampini   /* Infer the local solver */
178674ae819SStefano Zampini   ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr);
179674ae819SStefano Zampini   ierr = VecGetSize(pcis->vec1_D,&n_I);CHKERRQ(ierr);
180674ae819SStefano Zampini   ierr = VecGetSize(pcbddc->vec1_R,&n_R);CHKERRQ(ierr);
181674ae819SStefano Zampini   if (basis_dofs == n_I) {
182674ae819SStefano Zampini     /* Dirichlet solver */
183674ae819SStefano Zampini     local_ksp = &pcbddc->ksp_D;
184674ae819SStefano Zampini   } else if (basis_dofs == n_R) {
185674ae819SStefano Zampini     /* Neumann solver */
186674ae819SStefano Zampini     local_ksp = &pcbddc->ksp_R;
187674ae819SStefano Zampini   } else {
188674ae819SStefano 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);
189674ae819SStefano Zampini   }
19081d9aea3SBarry Smith   ierr = KSPGetOperators(*local_ksp,&local_mat,&local_pmat);CHKERRQ(ierr);
191674ae819SStefano Zampini 
192674ae819SStefano Zampini   /* Get null space vecs */
193674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr);
194674ae819SStefano Zampini   basis_size = nnsp_size;
195674ae819SStefano Zampini   if (nnsp_has_cnst) {
196674ae819SStefano Zampini     basis_size++;
197674ae819SStefano Zampini   }
198674ae819SStefano Zampini 
199b8ffe317SStefano Zampini   if (basis_dofs) {
200674ae819SStefano Zampini      /* Create shell ctx */
201674ae819SStefano Zampini      ierr = PetscMalloc(sizeof(*shell_ctx),&shell_ctx);CHKERRQ(ierr);
202674ae819SStefano Zampini 
203674ae819SStefano Zampini      /* Create work vectors in shell context */
204674ae819SStefano Zampini      ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr);
205674ae819SStefano Zampini      ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr);
206674ae819SStefano Zampini      ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr);
207674ae819SStefano Zampini      ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr);
208674ae819SStefano Zampini      ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr);
209674ae819SStefano Zampini      ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr);
210674ae819SStefano Zampini      ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr);
211674ae819SStefano Zampini      ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr);
212674ae819SStefano Zampini 
213674ae819SStefano Zampini      /* Allocate workspace */
214674ae819SStefano Zampini      ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );CHKERRQ(ierr);
215674ae819SStefano Zampini      ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr);
216674ae819SStefano Zampini      ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
217674ae819SStefano Zampini      ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
218674ae819SStefano Zampini 
219674ae819SStefano Zampini      /* Restrict local null space on selected dofs (Dirichlet or Neumann)
220674ae819SStefano Zampini         and compute matrices N and K*N */
221674ae819SStefano Zampini      ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr);
222674ae819SStefano Zampini      ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr);
223674ae819SStefano Zampini      ierr = VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);CHKERRQ(ierr);
224b8ffe317SStefano Zampini   }
225b8ffe317SStefano Zampini 
226674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
227674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
228674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
229b8ffe317SStefano Zampini     if (basis_dofs) {
230674ae819SStefano Zampini       ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
231674ae819SStefano Zampini       ierr = VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
232674ae819SStefano Zampini       ierr = VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
233674ae819SStefano Zampini       ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
234674ae819SStefano Zampini       ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
235674ae819SStefano Zampini       ierr = VecResetArray(work1);CHKERRQ(ierr);
236674ae819SStefano Zampini       ierr = VecResetArray(work2);CHKERRQ(ierr);
237674ae819SStefano Zampini     }
238b8ffe317SStefano Zampini   }
239b8ffe317SStefano Zampini 
240b8ffe317SStefano Zampini   if (basis_dofs) {
241674ae819SStefano Zampini     if (nnsp_has_cnst) {
242674ae819SStefano Zampini       ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr);
243674ae819SStefano Zampini       ierr = VecSet(work1,one);CHKERRQ(ierr);
244674ae819SStefano Zampini       ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr);
245674ae819SStefano Zampini       ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
246674ae819SStefano Zampini       ierr = VecResetArray(work1);CHKERRQ(ierr);
247674ae819SStefano Zampini       ierr = VecResetArray(work2);CHKERRQ(ierr);
248674ae819SStefano Zampini     }
249674ae819SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
250674ae819SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
251674ae819SStefano Zampini     ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
252674ae819SStefano Zampini     ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr);
253674ae819SStefano Zampini     ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr);
254674ae819SStefano Zampini 
255674ae819SStefano Zampini     /* Assemble another Mat object in shell context */
256674ae819SStefano Zampini     ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr);
257674ae819SStefano Zampini     ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr);
258674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr);
259674ae819SStefano Zampini     ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr);
260674ae819SStefano Zampini     ierr = ISDestroy(&is_aux);CHKERRQ(ierr);
261785e854fSJed Brown     ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr);
262674ae819SStefano Zampini     for (k=0;k<basis_size;k++) {
263674ae819SStefano Zampini       ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr);
264674ae819SStefano Zampini       ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr);
265674ae819SStefano Zampini       ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr);
266674ae819SStefano Zampini       ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr);
267674ae819SStefano Zampini       ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr);
268674ae819SStefano Zampini       ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
269674ae819SStefano Zampini       for (i=0;i<basis_size;i++) {
270674ae819SStefano Zampini         array_mat[i*basis_size+k]=array[i];
271674ae819SStefano Zampini       }
272674ae819SStefano Zampini       ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr);
273674ae819SStefano Zampini     }
274674ae819SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr);
275674ae819SStefano Zampini     ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr);
276674ae819SStefano Zampini     ierr = PetscFree(array_mat);CHKERRQ(ierr);
277674ae819SStefano Zampini     ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr);
278674ae819SStefano Zampini     ierr = MatDestroy(&small_mat);CHKERRQ(ierr);
279674ae819SStefano Zampini     ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr);
280674ae819SStefano Zampini 
281674ae819SStefano Zampini     /* Rebuild local PC */
282674ae819SStefano Zampini     ierr = KSPGetPC(*local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr);
283674ae819SStefano Zampini     ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr);
284674ae819SStefano Zampini     ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr);
28523ee1639SBarry Smith     ierr = PCSetOperators(newpc,local_mat,local_mat);CHKERRQ(ierr);
286674ae819SStefano Zampini     ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
287674ae819SStefano Zampini     ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr);
288674ae819SStefano Zampini     ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr);
289674ae819SStefano Zampini     ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr);
290674ae819SStefano Zampini     ierr = PCSetUp(newpc);CHKERRQ(ierr);
291674ae819SStefano Zampini     ierr = KSPSetPC(*local_ksp,newpc);CHKERRQ(ierr);
292674ae819SStefano Zampini     ierr = PCDestroy(&newpc);CHKERRQ(ierr);
293674ae819SStefano Zampini     ierr = KSPSetUp(*local_ksp);CHKERRQ(ierr);
294b8ffe317SStefano Zampini   }
295674ae819SStefano Zampini   /* test */
296b8ffe317SStefano Zampini   if (pcbddc->dbg_flag && basis_dofs) {
297674ae819SStefano Zampini     KSP         check_ksp;
298674ae819SStefano Zampini     PC          check_pc;
299674ae819SStefano Zampini     Mat         test_mat;
300674ae819SStefano Zampini     Vec         work3;
301674ae819SStefano Zampini     PetscReal   test_err,lambda_min,lambda_max;
302674ae819SStefano Zampini     PetscBool   setsym,issym=PETSC_FALSE;
303b8ffe317SStefano Zampini     PetscInt    tabs;
304674ae819SStefano Zampini 
305b8ffe317SStefano Zampini     ierr = PetscViewerASCIIGetTab(pcbddc->dbg_viewer,&tabs);CHKERRQ(ierr);
306674ae819SStefano Zampini     ierr = KSPGetPC(*local_ksp,&check_pc);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);
310674ae819SStefano Zampini     ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr);
311674ae819SStefano Zampini     ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr);
312674ae819SStefano Zampini     ierr = VecCopy(work1,work2);CHKERRQ(ierr);
313674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr);
314674ae819SStefano Zampini     ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr);
315674ae819SStefano Zampini     ierr = VecAXPY(work1,m_one,work2);CHKERRQ(ierr);
316674ae819SStefano Zampini     ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr);
317b8ffe317SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank);
318b8ffe317SStefano Zampini     ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_FALSE);CHKERRQ(ierr);
319674ae819SStefano Zampini     if (basis_dofs == n_I) {
320b8ffe317SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Dirichlet ");
321674ae819SStefano Zampini     } else {
322b8ffe317SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Neumann ");
323674ae819SStefano Zampini     }
324b8ffe317SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"solver is :%1.14e\n",test_err);
325b8ffe317SStefano Zampini     ierr = PetscViewerASCIISetTab(pcbddc->dbg_viewer,tabs);CHKERRQ(ierr);
326b8ffe317SStefano Zampini     ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
327674ae819SStefano Zampini 
328674ae819SStefano Zampini     ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr);
329674ae819SStefano Zampini     ierr = MatShift(test_mat,one);CHKERRQ(ierr);
330674ae819SStefano Zampini     ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr);
331674ae819SStefano Zampini     ierr = MatDestroy(&test_mat);CHKERRQ(ierr);
332b8ffe317SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err);
333674ae819SStefano Zampini 
334674ae819SStefano Zampini     /* Create ksp object suitable for extreme eigenvalues' estimation */
335674ae819SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr);
33623ee1639SBarry Smith     ierr = KSPSetOperators(check_ksp,local_mat,local_mat);CHKERRQ(ierr);
337674ae819SStefano Zampini     ierr = KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr);
338674ae819SStefano Zampini     ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr);
339674ae819SStefano Zampini     ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr);
340674ae819SStefano Zampini     if (issym) {
341674ae819SStefano Zampini       ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr);
342674ae819SStefano Zampini     }
343674ae819SStefano Zampini     ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr);
344674ae819SStefano Zampini     ierr = KSPSetUp(check_ksp);CHKERRQ(ierr);
345674ae819SStefano Zampini     ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr);
346674ae819SStefano Zampini     ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr);
347674ae819SStefano Zampini     ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr);
348674ae819SStefano Zampini     ierr = VecAXPY(work2,m_one,work1);CHKERRQ(ierr);
349674ae819SStefano Zampini     ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr);
350674ae819SStefano Zampini     ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr);
351674ae819SStefano Zampini     ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr);
352b8ffe317SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_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);
353674ae819SStefano Zampini     ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr);
354674ae819SStefano Zampini     ierr = VecDestroy(&work1);CHKERRQ(ierr);
355674ae819SStefano Zampini     ierr = VecDestroy(&work2);CHKERRQ(ierr);
356674ae819SStefano Zampini     ierr = VecDestroy(&work3);CHKERRQ(ierr);
357674ae819SStefano Zampini   }
358b8ffe317SStefano Zampini   /* all processes shoud call this, even the void ones */
359b8ffe317SStefano Zampini   if (pcbddc->dbg_flag) {
360b8ffe317SStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
361b8ffe317SStefano Zampini   }
362674ae819SStefano Zampini   PetscFunctionReturn(0);
363674ae819SStefano Zampini }
364674ae819SStefano Zampini 
365674ae819SStefano Zampini #undef __FUNCT__
366674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAdaptGlobal"
367674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc)
368674ae819SStefano Zampini {
369674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
370674ae819SStefano Zampini   KSP            inv_change;
371674ae819SStefano Zampini   const Vec      *nsp_vecs;
372674ae819SStefano Zampini   Vec            *new_nsp_vecs;
373674ae819SStefano Zampini   PetscInt       i,nsp_size,new_nsp_size,start_new;
374674ae819SStefano Zampini   PetscBool      nsp_has_cnst;
375674ae819SStefano Zampini   MatNullSpace   new_nsp;
376674ae819SStefano Zampini   PetscErrorCode ierr;
377674ae819SStefano Zampini 
378674ae819SStefano Zampini   PetscFunctionBegin;
379674ae819SStefano Zampini   /* create KSP for change of basis */
380*906d46d4SStefano Zampini   ierr = MatGetSize(pcbddc->ChangeOfBasisMatrix,&i,NULL);CHKERRQ(ierr);
381*906d46d4SStefano Zampini   ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&inv_change);CHKERRQ(ierr);
38223ee1639SBarry Smith   ierr = KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
383*906d46d4SStefano Zampini   ierr = KSPSetTolerances(inv_change,1.e-8,1.e-8,PETSC_DEFAULT,2*i);CHKERRQ(ierr);
384*906d46d4SStefano Zampini   if (pcbddc->dbg_flag) {
385*906d46d4SStefano Zampini     ierr = KSPMonitorSet(inv_change,KSPMonitorDefault,pcbddc->dbg_viewer,NULL);CHKERRQ(ierr);
386*906d46d4SStefano Zampini   }
387674ae819SStefano Zampini   ierr = KSPSetUp(inv_change);CHKERRQ(ierr);
388*906d46d4SStefano Zampini 
389674ae819SStefano Zampini   /* get nullspace and transform it */
390674ae819SStefano Zampini   ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr);
391674ae819SStefano Zampini   new_nsp_size = nsp_size;
392674ae819SStefano Zampini   if (nsp_has_cnst) {
393674ae819SStefano Zampini     new_nsp_size++;
394674ae819SStefano Zampini   }
395785e854fSJed Brown   ierr = PetscMalloc1(new_nsp_size,&new_nsp_vecs);CHKERRQ(ierr);
396674ae819SStefano Zampini   for (i=0;i<new_nsp_size;i++) {
397*906d46d4SStefano Zampini     ierr = MatGetVecs(pcbddc->ChangeOfBasisMatrix,&new_nsp_vecs[i],NULL);CHKERRQ(ierr);
398674ae819SStefano Zampini   }
399*906d46d4SStefano Zampini 
400674ae819SStefano Zampini   start_new = 0;
401674ae819SStefano Zampini   if (nsp_has_cnst) {
402674ae819SStefano Zampini     start_new = 1;
403674ae819SStefano Zampini     ierr = VecSet(new_nsp_vecs[0],1.0);CHKERRQ(ierr);
404*906d46d4SStefano Zampini     if (pcbddc->dbg_flag) {
405*906d46d4SStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
406*906d46d4SStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Mapping constant in nullspace\n");CHKERRQ(ierr);
407*906d46d4SStefano Zampini     }
408*906d46d4SStefano Zampini     ierr = KSPSolve(inv_change,new_nsp_vecs[0],new_nsp_vecs[0]);CHKERRQ(ierr);
409674ae819SStefano Zampini   }
410674ae819SStefano Zampini   for (i=0;i<nsp_size;i++) {
411*906d46d4SStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
412*906d46d4SStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Mapping %dth vector in nullspace\n",i);CHKERRQ(ierr);
413*906d46d4SStefano Zampini     ierr = KSPSolve(inv_change,nsp_vecs[i],new_nsp_vecs[i+start_new]);CHKERRQ(ierr);
414674ae819SStefano Zampini   }
4159a7d3425SStefano Zampini   ierr = PCBDDCOrthonormalizeVecs(new_nsp_size,new_nsp_vecs);CHKERRQ(ierr);
41698a51de6SStefano Zampini   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)pc),PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);CHKERRQ(ierr);
417674ae819SStefano Zampini   ierr = PCBDDCSetNullSpace(pc,new_nsp);CHKERRQ(ierr);
418674ae819SStefano Zampini 
419*906d46d4SStefano Zampini   /* free */
420*906d46d4SStefano Zampini   ierr = KSPDestroy(&inv_change);CHKERRQ(ierr);
421*906d46d4SStefano Zampini   ierr = MatNullSpaceDestroy(&new_nsp);CHKERRQ(ierr);
422674ae819SStefano Zampini   for (i=0;i<new_nsp_size;i++) {
423674ae819SStefano Zampini     ierr = VecDestroy(&new_nsp_vecs[i]);CHKERRQ(ierr);
424674ae819SStefano Zampini   }
425674ae819SStefano Zampini   ierr = PetscFree(new_nsp_vecs);CHKERRQ(ierr);
426*906d46d4SStefano Zampini 
427*906d46d4SStefano Zampini   /* check */
428*906d46d4SStefano Zampini   if (pcbddc->dbg_flag) {
429*906d46d4SStefano Zampini     PetscBool nsp_t=PETSC_FALSE;
430*906d46d4SStefano Zampini     Mat       temp_mat;
431*906d46d4SStefano Zampini     Mat_IS*   matis = (Mat_IS*)pc->pmat->data;
432*906d46d4SStefano Zampini 
433*906d46d4SStefano Zampini     temp_mat = matis->A;
434*906d46d4SStefano Zampini     matis->A = pcbddc->local_mat;
435*906d46d4SStefano Zampini     pcbddc->local_mat = temp_mat;
436*906d46d4SStefano Zampini     ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr);
437*906d46d4SStefano Zampini     ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"Check nullspace with change of basis: %d\n",nsp_t);CHKERRQ(ierr);
438*906d46d4SStefano Zampini     temp_mat = matis->A;
439*906d46d4SStefano Zampini     matis->A = pcbddc->local_mat;
440*906d46d4SStefano Zampini     pcbddc->local_mat = temp_mat;
441*906d46d4SStefano Zampini   }
442674ae819SStefano Zampini   PetscFunctionReturn(0);
443674ae819SStefano Zampini }
444