xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 7c334f0243eb3f00fc22f2f6efb0b24b1906cf86)
123ce1328SBarry Smith 
2b4319ba4SBarry Smith #include "src/ksp/pc/impls/is/pcis.h"
3b4319ba4SBarry Smith 
4b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
5b4319ba4SBarry Smith /*
6b4319ba4SBarry Smith    PCISSetUp -
7b4319ba4SBarry Smith */
8b4319ba4SBarry Smith #undef __FUNCT__
9b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp"
10dfbe8321SBarry Smith PetscErrorCode PCISSetUp(PC pc)
11b4319ba4SBarry Smith {
12b4319ba4SBarry Smith   PC_IS      *pcis = (PC_IS*)(pc->data);
13b4319ba4SBarry Smith   Mat_IS     *matis = (Mat_IS*)pc->mat->data;
146849ba73SBarry Smith   int        i;
156849ba73SBarry Smith   PetscErrorCode  ierr;
16b4319ba4SBarry Smith   PetscTruth flg;
17b4319ba4SBarry Smith 
18b4319ba4SBarry Smith   PetscFunctionBegin;
19b4319ba4SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr);
20b4319ba4SBarry Smith   if (!flg){
211302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
22b4319ba4SBarry Smith   }
23b4319ba4SBarry Smith 
24b4319ba4SBarry Smith   pcis->pure_neumann = matis->pure_neumann;
25b4319ba4SBarry Smith 
26b4319ba4SBarry Smith   /*
27b4319ba4SBarry Smith     Creating the local vector vec1_N, containing the inverse of the number
28b4319ba4SBarry Smith     of subdomains to which each local node (either owned or ghost)
29b4319ba4SBarry Smith     pertains. To accomplish that, we scatter local vectors of 1's to
30b4319ba4SBarry Smith     a global vector (adding the values); scatter the result back to
31b4319ba4SBarry Smith     local vectors and finally invert the result.
32b4319ba4SBarry Smith   */
33b4319ba4SBarry Smith   {
34b4319ba4SBarry Smith     Vec    counter;
35b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
36b4319ba4SBarry Smith     ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
3723ce1328SBarry Smith     ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
38b4319ba4SBarry Smith     ierr = VecSet(&zero,counter);CHKERRQ(ierr);
39b4319ba4SBarry Smith     ierr = VecSet(&one,pcis->vec1_N);CHKERRQ(ierr);
40b4319ba4SBarry Smith     ierr = VecScatterBegin(pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);CHKERRQ(ierr);
41b4319ba4SBarry Smith     ierr = VecScatterEnd  (pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);CHKERRQ(ierr);
42b4319ba4SBarry Smith     ierr = VecScatterBegin(counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);CHKERRQ(ierr);
43b4319ba4SBarry Smith     ierr = VecScatterEnd  (counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);CHKERRQ(ierr);
44b4319ba4SBarry Smith     ierr = VecDestroy(counter);CHKERRQ(ierr);
45b4319ba4SBarry Smith   }
46b4319ba4SBarry Smith   /*
47b4319ba4SBarry Smith     Creating local and global index sets for interior and
48b4319ba4SBarry Smith     inteface nodes. Notice that interior nodes have D[i]==1.0.
49b4319ba4SBarry Smith   */
50b4319ba4SBarry Smith   {
51b4319ba4SBarry Smith     int     n_I;
52b4319ba4SBarry Smith     int    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
53b4319ba4SBarry Smith     PetscScalar *array;
54b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
55b4319ba4SBarry Smith     ierr = VecGetSize(pcis->vec1_N,&pcis->n);CHKERRQ(ierr);
56b4319ba4SBarry Smith     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
57b4319ba4SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(int),&idx_I_local);CHKERRQ(ierr);
58b4319ba4SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(int),&idx_B_local);CHKERRQ(ierr);
59b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
60b4319ba4SBarry Smith       if (array[i] == 1.0) { idx_I_local[n_I]       = i; n_I++;       }
61b4319ba4SBarry Smith       else                 { idx_B_local[pcis->n_B] = i; pcis->n_B++; }
62b4319ba4SBarry Smith     }
63b4319ba4SBarry Smith     /* Getting the global numbering */
64b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
65b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
66b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
67b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I,      idx_I_local,idx_I_global);CHKERRQ(ierr);
68b4319ba4SBarry Smith     /* Creating the index sets. */
69b4319ba4SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local, &pcis->is_B_local);CHKERRQ(ierr);
70b4319ba4SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,&pcis->is_B_global);CHKERRQ(ierr);
71b4319ba4SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,n_I      ,idx_I_local, &pcis->is_I_local);CHKERRQ(ierr);
72b4319ba4SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,n_I      ,idx_I_global,&pcis->is_I_global);CHKERRQ(ierr);
73b4319ba4SBarry Smith     /* Freeing memory and restoring arrays */
74b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
75b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
76b4319ba4SBarry Smith     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
77b4319ba4SBarry Smith   }
78b4319ba4SBarry Smith 
79b4319ba4SBarry Smith   /*
80b4319ba4SBarry Smith     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
81b4319ba4SBarry Smith     is such that interior nodes come first than the interface ones, we have
82b4319ba4SBarry Smith 
83b4319ba4SBarry Smith     [           |      ]
84b4319ba4SBarry Smith     [    A_II   | A_IB ]
85b4319ba4SBarry Smith     A = [           |      ]
86b4319ba4SBarry Smith     [-----------+------]
87b4319ba4SBarry Smith     [    A_BI   | A_BB ]
88b4319ba4SBarry Smith   */
89b4319ba4SBarry Smith 
90b4319ba4SBarry Smith   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr);
91b4319ba4SBarry Smith   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
92b4319ba4SBarry Smith   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
93b4319ba4SBarry Smith   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
94b4319ba4SBarry Smith 
95b4319ba4SBarry Smith   /*
96b4319ba4SBarry Smith     Creating work vectors and arrays
97b4319ba4SBarry Smith   */
98b4319ba4SBarry Smith   /* pcis->vec1_N has already been created */
99b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
100b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
101b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
102b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
103b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
104b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
105b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
10623ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
107b4319ba4SBarry Smith   ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr);
108b4319ba4SBarry Smith 
109b4319ba4SBarry Smith   /* Creating the scatter contexts */
11023ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
111b4319ba4SBarry Smith   ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
11223ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
113b4319ba4SBarry Smith 
114b4319ba4SBarry Smith   /* Creating scaling "matrix" D, from information in vec1_N */
115b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
116b4319ba4SBarry Smith   ierr = VecScatterBegin(pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr);
117b4319ba4SBarry Smith   ierr = VecScatterEnd  (pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr);
118b4319ba4SBarry Smith   ierr = VecReciprocal(pcis->D);CHKERRQ(ierr);
119b4319ba4SBarry Smith 
120b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
121b4319ba4SBarry Smith 
122b4319ba4SBarry Smith   /*
123b4319ba4SBarry Smith     Creating the KSP contexts for the local Dirichlet and Neumann problems.
124b4319ba4SBarry Smith   */
125b4319ba4SBarry Smith   {
126b4319ba4SBarry Smith     PC  pc_ctx;
127b4319ba4SBarry Smith     /* Dirichlet */
128b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
129b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr);
130b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
131b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
132b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
133b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
134b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
135b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
136b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
137b4319ba4SBarry Smith     /* Neumann */
138b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
139b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr);
140b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
141b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
142b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
143b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
144b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
145b4319ba4SBarry Smith     {
146b4319ba4SBarry Smith       PetscTruth damp_fixed,
147b4319ba4SBarry Smith                  remove_nullspace_fixed,
148b4319ba4SBarry Smith                  set_damping_factor_floating,
149b4319ba4SBarry Smith                  not_damp_floating,
150b4319ba4SBarry Smith                  not_remove_nullspace_floating;
151b4319ba4SBarry Smith       PetscReal  fixed_factor,
152b4319ba4SBarry Smith                  floating_factor;
153b4319ba4SBarry Smith 
154b4319ba4SBarry Smith       ierr = PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
155b4319ba4SBarry Smith       if (!damp_fixed) { fixed_factor = 0.0; }
156b4319ba4SBarry Smith       ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_damp_fixed",&damp_fixed);CHKERRQ(ierr);
157b4319ba4SBarry Smith 
158b4319ba4SBarry Smith       ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed);CHKERRQ(ierr);
159b4319ba4SBarry Smith 
160b4319ba4SBarry Smith       ierr = PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_set_damping_factor_floating",
161b4319ba4SBarry Smith 			      &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
162b4319ba4SBarry Smith       if (!set_damping_factor_floating) { floating_factor = 0.0; }
163b4319ba4SBarry Smith       ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating);CHKERRQ(ierr);
164b4319ba4SBarry Smith       if (!set_damping_factor_floating) { floating_factor = 1.e-12; }
165b4319ba4SBarry Smith 
166b4319ba4SBarry Smith       ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_damp_floating",&not_damp_floating);CHKERRQ(ierr);
167b4319ba4SBarry Smith 
168b4319ba4SBarry Smith       ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating);CHKERRQ(ierr);
169b4319ba4SBarry Smith 
170b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
171b4319ba4SBarry Smith 	if (!(not_damp_floating)) {
172b4319ba4SBarry Smith 	  ierr = PCLUSetDamping (pc_ctx,floating_factor);CHKERRQ(ierr);
173b4319ba4SBarry Smith 	  ierr = PCILUSetDamping(pc_ctx,floating_factor);CHKERRQ(ierr);
174b4319ba4SBarry Smith 	}
175b4319ba4SBarry Smith 	if (!(not_remove_nullspace_floating)){
176b4319ba4SBarry Smith 	  MatNullSpace nullsp;
177ee512c37SSatish Balay 	  ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr);
178d8fd42c4SBarry Smith 	  ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
179b4319ba4SBarry Smith 	  ierr = MatNullSpaceDestroy(nullsp);CHKERRQ(ierr);
180b4319ba4SBarry Smith 	}
181b4319ba4SBarry Smith       } else {  /* fixed subdomain */
182b4319ba4SBarry Smith 	if (damp_fixed) {
183b4319ba4SBarry Smith 	  ierr = PCLUSetDamping (pc_ctx,fixed_factor);CHKERRQ(ierr);
184b4319ba4SBarry Smith 	  ierr = PCILUSetDamping(pc_ctx,fixed_factor);CHKERRQ(ierr);
185b4319ba4SBarry Smith 	}
186b4319ba4SBarry Smith 	if (remove_nullspace_fixed) {
187b4319ba4SBarry Smith 	  MatNullSpace nullsp;
188ee512c37SSatish Balay 	  ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr);
189d8fd42c4SBarry Smith 	  ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
190b4319ba4SBarry Smith 	  ierr = MatNullSpaceDestroy(nullsp);CHKERRQ(ierr);
191b4319ba4SBarry Smith 	}
192b4319ba4SBarry Smith       }
193b4319ba4SBarry Smith     }
194b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
195b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
196b4319ba4SBarry Smith   }
197b4319ba4SBarry Smith 
198*7c334f02SBarry Smith   ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
199b4319ba4SBarry Smith   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE;
200b4319ba4SBarry Smith 
201b4319ba4SBarry Smith   PetscFunctionReturn(0);
202b4319ba4SBarry Smith }
203b4319ba4SBarry Smith 
204b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
205b4319ba4SBarry Smith /*
206b4319ba4SBarry Smith    PCISDestroy -
207b4319ba4SBarry Smith */
208b4319ba4SBarry Smith #undef __FUNCT__
209b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy"
210dfbe8321SBarry Smith PetscErrorCode PCISDestroy(PC pc)
211b4319ba4SBarry Smith {
212b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS*)(pc->data);
213dfbe8321SBarry Smith   PetscErrorCode ierr;
214b4319ba4SBarry Smith 
215b4319ba4SBarry Smith   PetscFunctionBegin;
216b4319ba4SBarry Smith 
217b4319ba4SBarry Smith   if (pcis->is_B_local)  {ierr = ISDestroy(pcis->is_B_local);CHKERRQ(ierr);}
218b4319ba4SBarry Smith   if (pcis->is_I_local)  {ierr = ISDestroy(pcis->is_I_local);CHKERRQ(ierr);}
219b4319ba4SBarry Smith   if (pcis->is_B_global) {ierr = ISDestroy(pcis->is_B_global);CHKERRQ(ierr);}
220b4319ba4SBarry Smith   if (pcis->is_I_global) {ierr = ISDestroy(pcis->is_I_global);CHKERRQ(ierr);}
221b4319ba4SBarry Smith   if (pcis->A_II)        {ierr = MatDestroy(pcis->A_II);CHKERRQ(ierr);}
222b4319ba4SBarry Smith   if (pcis->A_IB)        {ierr = MatDestroy(pcis->A_IB);CHKERRQ(ierr);}
223b4319ba4SBarry Smith   if (pcis->A_BI)        {ierr = MatDestroy(pcis->A_BI);CHKERRQ(ierr);}
224b4319ba4SBarry Smith   if (pcis->A_BB)        {ierr = MatDestroy(pcis->A_BB);CHKERRQ(ierr);}
225b4319ba4SBarry Smith   if (pcis->D)           {ierr = VecDestroy(pcis->D);CHKERRQ(ierr);}
226b4319ba4SBarry Smith   if (pcis->ksp_N)      {ierr = KSPDestroy(pcis->ksp_N);CHKERRQ(ierr);}
227b4319ba4SBarry Smith   if (pcis->ksp_D)      {ierr = KSPDestroy(pcis->ksp_D);CHKERRQ(ierr);}
228b4319ba4SBarry Smith   if (pcis->vec1_N)      {ierr = VecDestroy(pcis->vec1_N);CHKERRQ(ierr);}
229b4319ba4SBarry Smith   if (pcis->vec2_N)      {ierr = VecDestroy(pcis->vec2_N);CHKERRQ(ierr);}
230b4319ba4SBarry Smith   if (pcis->vec1_D)      {ierr = VecDestroy(pcis->vec1_D);CHKERRQ(ierr);}
231b4319ba4SBarry Smith   if (pcis->vec2_D)      {ierr = VecDestroy(pcis->vec2_D);CHKERRQ(ierr);}
232b4319ba4SBarry Smith   if (pcis->vec3_D)      {ierr = VecDestroy(pcis->vec3_D);CHKERRQ(ierr);}
233b4319ba4SBarry Smith   if (pcis->vec1_B)      {ierr = VecDestroy(pcis->vec1_B);CHKERRQ(ierr);}
234b4319ba4SBarry Smith   if (pcis->vec2_B)      {ierr = VecDestroy(pcis->vec2_B);CHKERRQ(ierr);}
235b4319ba4SBarry Smith   if (pcis->vec3_B)      {ierr = VecDestroy(pcis->vec3_B);CHKERRQ(ierr);}
236b4319ba4SBarry Smith   if (pcis->vec1_global) {ierr = VecDestroy(pcis->vec1_global);CHKERRQ(ierr);}
237b4319ba4SBarry Smith   if (pcis->work_N)      {ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);}
238b4319ba4SBarry Smith   if (pcis->global_to_D) {ierr = VecScatterDestroy(pcis->global_to_D);CHKERRQ(ierr);}
239b4319ba4SBarry Smith   if (pcis->N_to_B)      {ierr = VecScatterDestroy(pcis->N_to_B);CHKERRQ(ierr);}
240b4319ba4SBarry Smith   if (pcis->global_to_B) {ierr = VecScatterDestroy(pcis->global_to_B);CHKERRQ(ierr);}
241b4319ba4SBarry Smith   if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) {
242b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
243b4319ba4SBarry Smith   }
244b4319ba4SBarry Smith 
245b4319ba4SBarry Smith   PetscFunctionReturn(0);
246b4319ba4SBarry Smith }
247b4319ba4SBarry Smith 
248b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
249b4319ba4SBarry Smith /*
250b4319ba4SBarry Smith    PCISCreate -
251b4319ba4SBarry Smith */
252b4319ba4SBarry Smith #undef __FUNCT__
253b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate"
2546849ba73SBarry Smith PetscErrorCode PCISCreate(PC pc)
255b4319ba4SBarry Smith {
256b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS*)(pc->data);
257b4319ba4SBarry Smith 
258b4319ba4SBarry Smith   PetscFunctionBegin;
259b4319ba4SBarry Smith 
260b4319ba4SBarry Smith   pcis->is_B_local  = 0;
261b4319ba4SBarry Smith   pcis->is_I_local  = 0;
262b4319ba4SBarry Smith   pcis->is_B_global = 0;
263b4319ba4SBarry Smith   pcis->is_I_global = 0;
264b4319ba4SBarry Smith   pcis->A_II        = 0;
265b4319ba4SBarry Smith   pcis->A_IB        = 0;
266b4319ba4SBarry Smith   pcis->A_BI        = 0;
267b4319ba4SBarry Smith   pcis->A_BB        = 0;
268b4319ba4SBarry Smith   pcis->D           = 0;
269b4319ba4SBarry Smith   pcis->ksp_N      = 0;
270b4319ba4SBarry Smith   pcis->ksp_D      = 0;
271b4319ba4SBarry Smith   pcis->vec1_N      = 0;
272b4319ba4SBarry Smith   pcis->vec2_N      = 0;
273b4319ba4SBarry Smith   pcis->vec1_D      = 0;
274b4319ba4SBarry Smith   pcis->vec2_D      = 0;
275b4319ba4SBarry Smith   pcis->vec3_D      = 0;
276b4319ba4SBarry Smith   pcis->vec1_B      = 0;
277b4319ba4SBarry Smith   pcis->vec2_B      = 0;
278b4319ba4SBarry Smith   pcis->vec3_B      = 0;
279b4319ba4SBarry Smith   pcis->vec1_global = 0;
280b4319ba4SBarry Smith   pcis->work_N      = 0;
281b4319ba4SBarry Smith   pcis->global_to_D = 0;
282b4319ba4SBarry Smith   pcis->N_to_B      = 0;
283b4319ba4SBarry Smith   pcis->global_to_B = 0;
284b4319ba4SBarry Smith   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE;
285b4319ba4SBarry Smith 
286b4319ba4SBarry Smith   PetscFunctionReturn(0);
287b4319ba4SBarry Smith }
288b4319ba4SBarry Smith 
289b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
290b4319ba4SBarry Smith /*
291b4319ba4SBarry Smith    PCISApplySchur -
292b4319ba4SBarry Smith 
293b4319ba4SBarry Smith    Input parameters:
294b4319ba4SBarry Smith .  pc - preconditioner context
295b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
296b4319ba4SBarry Smith 
297b4319ba4SBarry Smith    Output parameters:
298b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
299b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
300b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
301b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
302b4319ba4SBarry Smith 
303b4319ba4SBarry Smith */
304b4319ba4SBarry Smith #undef __FUNCT__
305b4319ba4SBarry Smith #define __FUNCT__ "PCIterSuApplySchur"
3066849ba73SBarry Smith PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
307b4319ba4SBarry Smith {
308dfbe8321SBarry Smith   PetscErrorCode ierr;
309b4319ba4SBarry Smith   PetscScalar m_one = -1.0;
310b4319ba4SBarry Smith   PC_IS       *pcis = (PC_IS*)(pc->data);
311b4319ba4SBarry Smith 
312b4319ba4SBarry Smith   PetscFunctionBegin;
313b4319ba4SBarry Smith 
314b4319ba4SBarry Smith   if (vec2_B == (Vec)0) { vec2_B = v; }
315b4319ba4SBarry Smith 
316b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
317b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
31823ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
319b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
320b4319ba4SBarry Smith   ierr = VecAXPY(&m_one,vec2_B,vec1_B);CHKERRQ(ierr);
321b4319ba4SBarry Smith 
322b4319ba4SBarry Smith   PetscFunctionReturn(0);
323b4319ba4SBarry Smith }
324b4319ba4SBarry Smith 
325b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
326b4319ba4SBarry Smith /*
327b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
328b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
329b4319ba4SBarry Smith    mode.
330b4319ba4SBarry Smith 
331b4319ba4SBarry Smith    Input parameters:
332b4319ba4SBarry Smith .  pc - preconditioner context
333b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
334b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
335b4319ba4SBarry Smith 
336b4319ba4SBarry Smith    Output parameter:
337b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
338b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
339b4319ba4SBarry Smith 
340b4319ba4SBarry Smith    Notes:
341b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
342b4319ba4SBarry Smith */
343b4319ba4SBarry Smith #undef __FUNCT__
344b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB"
3456849ba73SBarry Smith PetscErrorCode PCISScatterArrayNToVecB (PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
346b4319ba4SBarry Smith {
3476849ba73SBarry Smith   int         i, *idex;
3486849ba73SBarry Smith   PetscErrorCode ierr;
349b4319ba4SBarry Smith   PetscScalar *array_B;
350b4319ba4SBarry Smith   PC_IS       *pcis = (PC_IS*)(pc->data);
351b4319ba4SBarry Smith 
352b4319ba4SBarry Smith   PetscFunctionBegin;
353b4319ba4SBarry Smith 
354b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
355b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
356b4319ba4SBarry Smith 
357b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
358b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
359b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_B[i]  = array_N[idex[i]]; }
360b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
361b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_B[i] += array_N[idex[i]]; }
362b4319ba4SBarry Smith     }
363b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
364b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
365b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_N[idex[i]]  = array_B[i]; }
366b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
367b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] += array_B[i]; }
368b4319ba4SBarry Smith     }
369b4319ba4SBarry Smith   }
370b4319ba4SBarry Smith 
371b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
372b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
373b4319ba4SBarry Smith 
374b4319ba4SBarry Smith   PetscFunctionReturn(0);
375b4319ba4SBarry Smith }
376b4319ba4SBarry Smith 
377b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
378b4319ba4SBarry Smith /*
379b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
380b4319ba4SBarry Smith    More precisely, solves the problem:
381b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
382b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
383b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
384b4319ba4SBarry Smith 
385b4319ba4SBarry Smith    Input parameters:
386b4319ba4SBarry Smith .  pc - preconditioner context
387b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
388b4319ba4SBarry Smith 
389b4319ba4SBarry Smith    Output parameters:
390b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
391b4319ba4SBarry Smith        complement to b
392b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
393b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
394b4319ba4SBarry Smith 
395b4319ba4SBarry Smith */
396b4319ba4SBarry Smith #undef __FUNCT__
397b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur"
3986849ba73SBarry Smith PetscErrorCode PCISApplyInvSchur (PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
399b4319ba4SBarry Smith {
400dfbe8321SBarry Smith   PetscErrorCode ierr;
401b4319ba4SBarry Smith   PC_IS       *pcis = (PC_IS*)(pc->data);
402b4319ba4SBarry Smith   PetscScalar zero  = 0.0;
403b4319ba4SBarry Smith 
404b4319ba4SBarry Smith   PetscFunctionBegin;
405b4319ba4SBarry Smith 
406b4319ba4SBarry Smith   /*
407b4319ba4SBarry Smith     Neumann solvers.
408b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
409b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
410b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
411b4319ba4SBarry Smith     is stored in x.
412b4319ba4SBarry Smith   */
413b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
414b4319ba4SBarry Smith   ierr = VecSet(&zero,vec1_N);CHKERRQ(ierr);
415b4319ba4SBarry Smith   ierr = VecScatterBegin(b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);CHKERRQ(ierr);
416b4319ba4SBarry Smith   ierr = VecScatterEnd  (b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);CHKERRQ(ierr);
417b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
418b4319ba4SBarry Smith   {
419b4319ba4SBarry Smith     PetscTruth flg;
420b4319ba4SBarry Smith     ierr = PetscOptionsHasName(PETSC_NULL,"-pc_is_check_consistency",&flg);CHKERRQ(ierr);
421b4319ba4SBarry Smith     if (flg) {
422b4319ba4SBarry Smith       PetscScalar average;
423b4319ba4SBarry Smith       ierr = VecSum(vec1_N,&average);CHKERRQ(ierr);
424b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
425b4319ba4SBarry Smith       if (pcis->pure_neumann) {
426b4319ba4SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is floating. Average = % 1.14e\n",
427b4319ba4SBarry Smith                                              PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
428b4319ba4SBarry Smith       } else {
429b4319ba4SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is fixed.    Average = % 1.14e\n",
430b4319ba4SBarry Smith                                              PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
431b4319ba4SBarry Smith       }
432b4319ba4SBarry Smith       PetscViewerFlush(PETSC_VIEWER_STDOUT_(pc->comm));
433b4319ba4SBarry Smith     }
434b4319ba4SBarry Smith   }
435b4319ba4SBarry Smith   /* Solving the system for vec2_N */
43623ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
437b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
438b4319ba4SBarry Smith   ierr = VecScatterBegin(vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr);
439b4319ba4SBarry Smith   ierr = VecScatterEnd  (vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr);
440b4319ba4SBarry Smith 
441b4319ba4SBarry Smith   PetscFunctionReturn(0);
442b4319ba4SBarry Smith }
443b4319ba4SBarry Smith 
444b4319ba4SBarry Smith 
445b4319ba4SBarry Smith 
446b4319ba4SBarry Smith 
447b4319ba4SBarry Smith 
448b4319ba4SBarry Smith 
449b4319ba4SBarry Smith 
450b4319ba4SBarry Smith 
451b4319ba4SBarry Smith 
452