xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 831a100d645ee317f94a2a60304e7ccc62bfb9ae)
123ce1328SBarry Smith 
2*831a100dSStefano Zampini #include "pcis.h" /*I "petscpc.h" I*/
3*831a100dSStefano Zampini 
4*831a100dSStefano Zampini EXTERN_C_BEGIN
5*831a100dSStefano Zampini #undef __FUNCT__
6*831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS"
7*831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
8*831a100dSStefano Zampini {
9*831a100dSStefano Zampini   PC_IS  *pcis = (PC_IS*)pc->data;
10*831a100dSStefano Zampini 
11*831a100dSStefano Zampini   PetscFunctionBegin;
12*831a100dSStefano Zampini   pcis->scaling_factor = scal;
13*831a100dSStefano Zampini   PetscFunctionReturn(0);
14*831a100dSStefano Zampini }
15*831a100dSStefano Zampini EXTERN_C_END
16*831a100dSStefano Zampini 
17*831a100dSStefano Zampini #undef __FUNCT__
18*831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor"
19*831a100dSStefano Zampini /*@
20*831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
21*831a100dSStefano Zampini 
22*831a100dSStefano Zampini    Not collective
23*831a100dSStefano Zampini 
24*831a100dSStefano Zampini    Input Parameters:
25*831a100dSStefano Zampini +  pc - the preconditioning context
26*831a100dSStefano Zampini -  scal - scaling factor for the subdomain
27*831a100dSStefano Zampini 
28*831a100dSStefano Zampini    Level: intermediate
29*831a100dSStefano Zampini 
30*831a100dSStefano Zampini    Notes:
31*831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
32*831a100dSStefano Zampini 
33*831a100dSStefano Zampini .seealso: PCBDDC
34*831a100dSStefano Zampini @*/
35*831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
36*831a100dSStefano Zampini {
37*831a100dSStefano Zampini   PetscErrorCode ierr;
38*831a100dSStefano Zampini 
39*831a100dSStefano Zampini   PetscFunctionBegin;
40*831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
41*831a100dSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr);
42*831a100dSStefano Zampini   PetscFunctionReturn(0);
43*831a100dSStefano Zampini }
44*831a100dSStefano Zampini 
45b4319ba4SBarry Smith 
46b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
47b4319ba4SBarry Smith /*
48b4319ba4SBarry Smith    PCISSetUp -
49b4319ba4SBarry Smith */
50b4319ba4SBarry Smith #undef __FUNCT__
51b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp"
527087cfbeSBarry Smith PetscErrorCode  PCISSetUp(PC pc)
53b4319ba4SBarry Smith {
54b4319ba4SBarry Smith   PC_IS           *pcis = (PC_IS*)(pc->data);
55b4319ba4SBarry Smith   Mat_IS          *matis = (Mat_IS*)pc->mat->data;
5613f74950SBarry Smith   PetscInt        i;
576849ba73SBarry Smith   PetscErrorCode  ierr;
58ace3abfcSBarry Smith   PetscBool       flg;
59*831a100dSStefano Zampini   Vec    counter;
60b4319ba4SBarry Smith 
61b4319ba4SBarry Smith   PetscFunctionBegin;
62251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr);
63e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
64b4319ba4SBarry Smith 
65b4319ba4SBarry Smith   pcis->pure_neumann = matis->pure_neumann;
66b4319ba4SBarry Smith 
67b4319ba4SBarry Smith   /*
68b4319ba4SBarry Smith     Creating the local vector vec1_N, containing the inverse of the number
69b4319ba4SBarry Smith     of subdomains to which each local node (either owned or ghost)
70b4319ba4SBarry Smith     pertains. To accomplish that, we scatter local vectors of 1's to
71b4319ba4SBarry Smith     a global vector (adding the values); scatter the result back to
72b4319ba4SBarry Smith     local vectors and finally invert the result.
73b4319ba4SBarry Smith   */
74b4319ba4SBarry Smith   ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
7523ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
76efb30889SBarry Smith   ierr = VecSet(counter,0.0);CHKERRQ(ierr);
77efb30889SBarry Smith   ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr);
78ca9f406cSSatish Balay   ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79ca9f406cSSatish Balay   ierr = VecScatterEnd  (matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
80ca9f406cSSatish Balay   ierr = VecScatterBegin(matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81ca9f406cSSatish Balay   ierr = VecScatterEnd  (matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
82b4319ba4SBarry Smith   /*
83b4319ba4SBarry Smith     Creating local and global index sets for interior and
84b4319ba4SBarry Smith     inteface nodes. Notice that interior nodes have D[i]==1.0.
85b4319ba4SBarry Smith   */
86b4319ba4SBarry Smith   {
8713f74950SBarry Smith     PetscInt     n_I;
8813f74950SBarry Smith     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
89b4319ba4SBarry Smith     PetscScalar *array;
90b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
91b4319ba4SBarry Smith     ierr = VecGetSize(pcis->vec1_N,&pcis->n);CHKERRQ(ierr);
92b4319ba4SBarry Smith     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
9313f74950SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);CHKERRQ(ierr);
9413f74950SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);CHKERRQ(ierr);
95b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
96b4319ba4SBarry Smith       if (array[i] == 1.0) { idx_I_local[n_I]       = i; n_I++;       }
97b4319ba4SBarry Smith       else                 { idx_B_local[pcis->n_B] = i; pcis->n_B++; }
98b4319ba4SBarry Smith     }
99b4319ba4SBarry Smith     /* Getting the global numbering */
100b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
101b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
102b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
103b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I,      idx_I_local,idx_I_global);CHKERRQ(ierr);
104b4319ba4SBarry Smith     /* Creating the index sets. */
10570b3c8c7SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
10670b3c8c7SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
10770b3c8c7SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,n_I      ,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
10870b3c8c7SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,n_I      ,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
109b4319ba4SBarry Smith     /* Freeing memory and restoring arrays */
110b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
111b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
112b4319ba4SBarry Smith     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
113b4319ba4SBarry Smith   }
114b4319ba4SBarry Smith 
115b4319ba4SBarry Smith   /*
116b4319ba4SBarry Smith     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
117b4319ba4SBarry Smith     is such that interior nodes come first than the interface ones, we have
118b4319ba4SBarry Smith 
119b4319ba4SBarry Smith     [           |      ]
120b4319ba4SBarry Smith     [    A_II   | A_IB ]
121b4319ba4SBarry Smith     A = [           |      ]
122b4319ba4SBarry Smith     [-----------+------]
123b4319ba4SBarry Smith     [    A_BI   | A_BB ]
124b4319ba4SBarry Smith   */
125b4319ba4SBarry Smith 
1264aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr);
1274aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
1284aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
1294aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
130b4319ba4SBarry Smith 
131b4319ba4SBarry Smith   /*
132b4319ba4SBarry Smith     Creating work vectors and arrays
133b4319ba4SBarry Smith   */
134b4319ba4SBarry Smith   /* pcis->vec1_N has already been created */
135b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
136b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
137b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
138b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
139b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
140b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
141b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
14223ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
143b4319ba4SBarry Smith   ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr);
144b4319ba4SBarry Smith 
145b4319ba4SBarry Smith   /* Creating the scatter contexts */
14623ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
147b4319ba4SBarry Smith   ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
14823ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
149b4319ba4SBarry Smith 
150*831a100dSStefano Zampini   /* Creating scaling "matrix" D */
151*831a100dSStefano Zampini   if( !pcis->D ) {
152*831a100dSStefano Zampini     ierr = VecSet(counter,0.0);CHKERRQ(ierr);
153*831a100dSStefano Zampini     ierr = VecSet(pcis->vec1_N,pcis->scaling_factor);CHKERRQ(ierr);
154*831a100dSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
155*831a100dSStefano Zampini     ierr = VecScatterEnd  (matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
156*831a100dSStefano Zampini     ierr = VecScatterBegin(matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
157*831a100dSStefano Zampini     ierr = VecScatterEnd  (matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
158b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
159ca9f406cSSatish Balay     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
160ca9f406cSSatish Balay     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
161b4319ba4SBarry Smith     ierr = VecReciprocal(pcis->D);CHKERRQ(ierr);
162*831a100dSStefano Zampini     ierr = VecScale(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
163*831a100dSStefano Zampini   }
164b4319ba4SBarry Smith 
165b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
166b4319ba4SBarry Smith 
167b4319ba4SBarry Smith   /*
168b4319ba4SBarry Smith     Creating the KSP contexts for the local Dirichlet and Neumann problems.
169b4319ba4SBarry Smith   */
170b4319ba4SBarry Smith   {
171b4319ba4SBarry Smith     PC  pc_ctx;
172b4319ba4SBarry Smith     /* Dirichlet */
173b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
1741cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
175b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr);
176b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
177b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
178b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
179b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
180b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
181b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
182b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
183b4319ba4SBarry Smith     /* Neumann */
184b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
1851cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
186b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr);
187b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
188b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
189b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
190b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
191b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
192b4319ba4SBarry Smith     {
193ace3abfcSBarry Smith       PetscBool  damp_fixed = PETSC_FALSE,
19490d69ab7SBarry Smith                  remove_nullspace_fixed = PETSC_FALSE,
19590d69ab7SBarry Smith                  set_damping_factor_floating = PETSC_FALSE,
19690d69ab7SBarry Smith                  not_damp_floating = PETSC_FALSE,
19790d69ab7SBarry Smith                  not_remove_nullspace_floating = PETSC_FALSE;
198b4319ba4SBarry Smith       PetscReal  fixed_factor,
199b4319ba4SBarry Smith                  floating_factor;
200b4319ba4SBarry Smith 
2017adad957SLisandro Dalcin       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
202b4319ba4SBarry Smith       if (!damp_fixed) { fixed_factor = 0.0; }
203acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,PETSC_NULL);CHKERRQ(ierr);
204b4319ba4SBarry Smith 
205acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,PETSC_NULL);CHKERRQ(ierr);
206b4319ba4SBarry Smith 
2077adad957SLisandro Dalcin       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
208b4319ba4SBarry Smith 			      &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
209b4319ba4SBarry Smith       if (!set_damping_factor_floating) { floating_factor = 0.0; }
210acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,PETSC_NULL);CHKERRQ(ierr);
211b4319ba4SBarry Smith       if (!set_damping_factor_floating) { floating_factor = 1.e-12; }
212b4319ba4SBarry Smith 
213acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,PETSC_NULL);CHKERRQ(ierr);
214b4319ba4SBarry Smith 
215acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,PETSC_NULL);CHKERRQ(ierr);
216b4319ba4SBarry Smith 
217b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
218b4319ba4SBarry Smith 	if (!(not_damp_floating)) {
219c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
220c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
221b4319ba4SBarry Smith 	}
222b4319ba4SBarry Smith 	if (!(not_remove_nullspace_floating)){
223b4319ba4SBarry Smith 	  MatNullSpace nullsp;
224ee512c37SSatish Balay 	  ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr);
225d8fd42c4SBarry Smith 	  ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
226d34fcf5fSBarry Smith 	  ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
227b4319ba4SBarry Smith 	}
228b4319ba4SBarry Smith       } else {  /* fixed subdomain */
229b4319ba4SBarry Smith 	if (damp_fixed) {
230c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
231c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
232b4319ba4SBarry Smith 	}
233b4319ba4SBarry Smith 	if (remove_nullspace_fixed) {
234b4319ba4SBarry Smith 	  MatNullSpace nullsp;
235ee512c37SSatish Balay 	  ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr);
236d8fd42c4SBarry Smith 	  ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
237d34fcf5fSBarry Smith 	  ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
238b4319ba4SBarry Smith 	}
239b4319ba4SBarry Smith       }
240b4319ba4SBarry Smith     }
241b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
242b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
243b4319ba4SBarry Smith   }
244b4319ba4SBarry Smith 
2457c334f02SBarry Smith   ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
246b4319ba4SBarry Smith   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE;
247*831a100dSStefano Zampini   ierr = VecDestroy(&counter);CHKERRQ(ierr);
248b4319ba4SBarry Smith   PetscFunctionReturn(0);
249b4319ba4SBarry Smith }
250b4319ba4SBarry Smith 
251b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
252b4319ba4SBarry Smith /*
253b4319ba4SBarry Smith    PCISDestroy -
254b4319ba4SBarry Smith */
255b4319ba4SBarry Smith #undef __FUNCT__
256b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy"
2577087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
258b4319ba4SBarry Smith {
259b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
260dfbe8321SBarry Smith   PetscErrorCode ierr;
261b4319ba4SBarry Smith 
262b4319ba4SBarry Smith   PetscFunctionBegin;
263fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
264fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
265fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
266fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
267fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
268fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
269fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
270fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
271fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
272fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
273fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
274fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
275fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
276fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
277fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
278fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
279fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
280fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
281fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
282fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
283fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
284fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
285fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
28605b42c5fSBarry Smith   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
287b4319ba4SBarry Smith   if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) {
288b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
289b4319ba4SBarry Smith   }
290b4319ba4SBarry Smith   PetscFunctionReturn(0);
291b4319ba4SBarry Smith }
292b4319ba4SBarry Smith 
293b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
294b4319ba4SBarry Smith /*
295b4319ba4SBarry Smith    PCISCreate -
296b4319ba4SBarry Smith */
297b4319ba4SBarry Smith #undef __FUNCT__
298b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate"
2997087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
300b4319ba4SBarry Smith {
301b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS*)(pc->data);
302*831a100dSStefano Zampini   PetscErrorCode ierr;
303b4319ba4SBarry Smith 
304b4319ba4SBarry Smith   PetscFunctionBegin;
305b4319ba4SBarry Smith   pcis->is_B_local  = 0;
306b4319ba4SBarry Smith   pcis->is_I_local  = 0;
307b4319ba4SBarry Smith   pcis->is_B_global = 0;
308b4319ba4SBarry Smith   pcis->is_I_global = 0;
309b4319ba4SBarry Smith   pcis->A_II        = 0;
310b4319ba4SBarry Smith   pcis->A_IB        = 0;
311b4319ba4SBarry Smith   pcis->A_BI        = 0;
312b4319ba4SBarry Smith   pcis->A_BB        = 0;
313b4319ba4SBarry Smith   pcis->D           = 0;
314b4319ba4SBarry Smith   pcis->ksp_N      = 0;
315b4319ba4SBarry Smith   pcis->ksp_D      = 0;
316b4319ba4SBarry Smith   pcis->vec1_N      = 0;
317b4319ba4SBarry Smith   pcis->vec2_N      = 0;
318b4319ba4SBarry Smith   pcis->vec1_D      = 0;
319b4319ba4SBarry Smith   pcis->vec2_D      = 0;
320b4319ba4SBarry Smith   pcis->vec3_D      = 0;
321b4319ba4SBarry Smith   pcis->vec1_B      = 0;
322b4319ba4SBarry Smith   pcis->vec2_B      = 0;
323b4319ba4SBarry Smith   pcis->vec3_B      = 0;
324b4319ba4SBarry Smith   pcis->vec1_global = 0;
325b4319ba4SBarry Smith   pcis->work_N      = 0;
326b4319ba4SBarry Smith   pcis->global_to_D = 0;
327b4319ba4SBarry Smith   pcis->N_to_B      = 0;
328b4319ba4SBarry Smith   pcis->global_to_B = 0;
329b4319ba4SBarry Smith   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE;
330*831a100dSStefano Zampini   pcis->scaling_factor = 1.0;
331*831a100dSStefano Zampini   /* composing functions */
332*831a100dSStefano Zampini   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainScalingFactor_C","PCISSetSubdomainScalingFactor_IS",
333*831a100dSStefano Zampini                     PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
334b4319ba4SBarry Smith   PetscFunctionReturn(0);
335b4319ba4SBarry Smith }
336b4319ba4SBarry Smith 
337b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
338b4319ba4SBarry Smith /*
339b4319ba4SBarry Smith    PCISApplySchur -
340b4319ba4SBarry Smith 
341b4319ba4SBarry Smith    Input parameters:
342b4319ba4SBarry Smith .  pc - preconditioner context
343b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
344b4319ba4SBarry Smith 
345b4319ba4SBarry Smith    Output parameters:
346b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
347b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
348b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
349b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
350b4319ba4SBarry Smith 
351b4319ba4SBarry Smith */
352b4319ba4SBarry Smith #undef __FUNCT__
353*831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur"
3547087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
355b4319ba4SBarry Smith {
356dfbe8321SBarry Smith   PetscErrorCode ierr;
357b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
358b4319ba4SBarry Smith 
359b4319ba4SBarry Smith   PetscFunctionBegin;
36013f74950SBarry Smith   if (!vec2_B) { vec2_B = v; }
361b4319ba4SBarry Smith 
362b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
363b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
36423ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
365b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
366efb30889SBarry Smith   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
367b4319ba4SBarry Smith   PetscFunctionReturn(0);
368b4319ba4SBarry Smith }
369b4319ba4SBarry Smith 
370b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
371b4319ba4SBarry Smith /*
372b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
373b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
374b4319ba4SBarry Smith    mode.
375b4319ba4SBarry Smith 
376b4319ba4SBarry Smith    Input parameters:
377b4319ba4SBarry Smith .  pc - preconditioner context
378b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
379b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
380b4319ba4SBarry Smith 
381b4319ba4SBarry Smith    Output parameter:
382b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
383b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
384b4319ba4SBarry Smith 
385b4319ba4SBarry Smith    Notes:
386b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
387b4319ba4SBarry Smith */
388b4319ba4SBarry Smith #undef __FUNCT__
389b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB"
3907087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB (PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
391b4319ba4SBarry Smith {
3925d0c19d7SBarry Smith   PetscInt       i;
3935d0c19d7SBarry Smith   const PetscInt *idex;
3946849ba73SBarry Smith   PetscErrorCode ierr;
395b4319ba4SBarry Smith   PetscScalar    *array_B;
396b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
397b4319ba4SBarry Smith 
398b4319ba4SBarry Smith   PetscFunctionBegin;
399b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
400b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
401b4319ba4SBarry Smith 
402b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
403b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
404b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_B[i]  = array_N[idex[i]]; }
405b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
406b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_B[i] += array_N[idex[i]]; }
407b4319ba4SBarry Smith     }
408b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
409b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
410b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_N[idex[i]]  = array_B[i]; }
411b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
412b4319ba4SBarry Smith       for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] += array_B[i]; }
413b4319ba4SBarry Smith     }
414b4319ba4SBarry Smith   }
415b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
416b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
417b4319ba4SBarry Smith   PetscFunctionReturn(0);
418b4319ba4SBarry Smith }
419b4319ba4SBarry Smith 
420b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
421b4319ba4SBarry Smith /*
422b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
423b4319ba4SBarry Smith    More precisely, solves the problem:
424b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
425b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
426b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
427b4319ba4SBarry Smith 
428b4319ba4SBarry Smith    Input parameters:
429b4319ba4SBarry Smith .  pc - preconditioner context
430b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
431b4319ba4SBarry Smith 
432b4319ba4SBarry Smith    Output parameters:
433b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
434b4319ba4SBarry Smith        complement to b
435b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
436b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
437b4319ba4SBarry Smith 
438b4319ba4SBarry Smith */
439b4319ba4SBarry Smith #undef __FUNCT__
440b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur"
4417087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur (PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
442b4319ba4SBarry Smith {
443dfbe8321SBarry Smith   PetscErrorCode ierr;
444b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
445b4319ba4SBarry Smith 
446b4319ba4SBarry Smith   PetscFunctionBegin;
447b4319ba4SBarry Smith   /*
448b4319ba4SBarry Smith     Neumann solvers.
449b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
450b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
451b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
452b4319ba4SBarry Smith     is stored in x.
453b4319ba4SBarry Smith   */
454b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
455efb30889SBarry Smith   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
456ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
457ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
458b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
459b4319ba4SBarry Smith   {
460ace3abfcSBarry Smith     PetscBool  flg = PETSC_FALSE;
461acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-pc_is_check_consistency",&flg,PETSC_NULL);CHKERRQ(ierr);
462b4319ba4SBarry Smith     if (flg) {
463b4319ba4SBarry Smith       PetscScalar average;
4643050cee2SBarry Smith       PetscViewer viewer;
4657adad957SLisandro Dalcin       ierr = PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);CHKERRQ(ierr);
4663050cee2SBarry Smith 
467b4319ba4SBarry Smith       ierr = VecSum(vec1_N,&average);CHKERRQ(ierr);
468b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
4697b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
470b4319ba4SBarry Smith       if (pcis->pure_neumann) {
4719f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
472b4319ba4SBarry Smith       } else {
4739f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
474b4319ba4SBarry Smith       }
4753050cee2SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
4767b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
477b4319ba4SBarry Smith     }
478b4319ba4SBarry Smith   }
479b4319ba4SBarry Smith   /* Solving the system for vec2_N */
48023ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
481b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
482ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
483ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
484b4319ba4SBarry Smith   PetscFunctionReturn(0);
485b4319ba4SBarry Smith }
486b4319ba4SBarry Smith 
487b4319ba4SBarry Smith 
488b4319ba4SBarry Smith 
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith 
491b4319ba4SBarry Smith 
492b4319ba4SBarry Smith 
493b4319ba4SBarry Smith 
494b4319ba4SBarry Smith 
495