xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 4843afd65bb254a27fba8bce48dfbc5b606537f1)
123ce1328SBarry Smith 
2ccd284c7SBarry Smith #include "../src/ksp/pc/impls/is/pcis.h" /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
4831a100dSStefano Zampini #undef __FUNCT__
5b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling_IS"
6b965d45eSStefano Zampini static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
7b965d45eSStefano Zampini {
8b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
9b965d45eSStefano Zampini 
10b965d45eSStefano Zampini   PetscFunctionBegin;
11b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
12b965d45eSStefano Zampini   PetscFunctionReturn(0);
13b965d45eSStefano Zampini }
14b965d45eSStefano Zampini 
15b965d45eSStefano Zampini #undef __FUNCT__
16b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling"
17b965d45eSStefano Zampini /*@
18b965d45eSStefano Zampini  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
19b965d45eSStefano Zampini                               local matrices' diagonal.
20b965d45eSStefano Zampini 
21b965d45eSStefano Zampini    Not collective
22b965d45eSStefano Zampini 
23b965d45eSStefano Zampini    Input Parameters:
24b965d45eSStefano Zampini +  pc - the preconditioning context
25b965d45eSStefano Zampini -  use - whether or not pcis use matrix diagonal to build partition of unity.
26b965d45eSStefano Zampini 
27b965d45eSStefano Zampini    Level: intermediate
28b965d45eSStefano Zampini 
29b965d45eSStefano Zampini    Notes:
30b965d45eSStefano Zampini 
31b965d45eSStefano Zampini .seealso: PCBDDC
32b965d45eSStefano Zampini @*/
33b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
34b965d45eSStefano Zampini {
35b965d45eSStefano Zampini   PetscErrorCode ierr;
36b965d45eSStefano Zampini 
37b965d45eSStefano Zampini   PetscFunctionBegin;
38b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
39b965d45eSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
40b965d45eSStefano Zampini   PetscFunctionReturn(0);
41b965d45eSStefano Zampini }
42b965d45eSStefano Zampini 
43b965d45eSStefano Zampini #undef __FUNCT__
44ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling_IS"
45ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
46ba1573a8SStefano Zampini {
47ba1573a8SStefano Zampini   PetscErrorCode ierr;
48ba1573a8SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
49ba1573a8SStefano Zampini 
50ba1573a8SStefano Zampini   PetscFunctionBegin;
51ba1573a8SStefano Zampini   ierr    = VecDestroy(&pcis->D);CHKERRQ(ierr);
52ba1573a8SStefano Zampini   ierr    = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr);
53ba1573a8SStefano Zampini   pcis->D = scaling_factors;
54ba1573a8SStefano Zampini   PetscFunctionReturn(0);
55ba1573a8SStefano Zampini }
56ba1573a8SStefano Zampini 
57ba1573a8SStefano Zampini #undef __FUNCT__
58ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling"
59ba1573a8SStefano Zampini /*@
60ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
61ba1573a8SStefano Zampini 
62ba1573a8SStefano Zampini    Not collective
63ba1573a8SStefano Zampini 
64ba1573a8SStefano Zampini    Input Parameters:
65ba1573a8SStefano Zampini +  pc - the preconditioning context
66ba1573a8SStefano Zampini -  scaling_factors - scaling factors for the subdomain
67ba1573a8SStefano Zampini 
68ba1573a8SStefano Zampini    Level: intermediate
69ba1573a8SStefano Zampini 
70ba1573a8SStefano Zampini    Notes:
71ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
72ba1573a8SStefano Zampini 
73ba1573a8SStefano Zampini .seealso: PCBDDC
74ba1573a8SStefano Zampini @*/
75ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
76ba1573a8SStefano Zampini {
77ba1573a8SStefano Zampini   PetscErrorCode ierr;
78ba1573a8SStefano Zampini 
79ba1573a8SStefano Zampini   PetscFunctionBegin;
80ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
81ba1573a8SStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr);
82ba1573a8SStefano Zampini   PetscFunctionReturn(0);
83ba1573a8SStefano Zampini }
84ba1573a8SStefano Zampini 
85ba1573a8SStefano Zampini #undef __FUNCT__
86831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS"
87831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
88831a100dSStefano Zampini {
89831a100dSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
90831a100dSStefano Zampini 
91831a100dSStefano Zampini   PetscFunctionBegin;
92831a100dSStefano Zampini   pcis->scaling_factor = scal;
93831a100dSStefano Zampini   PetscFunctionReturn(0);
94831a100dSStefano Zampini }
95831a100dSStefano Zampini 
96831a100dSStefano Zampini #undef __FUNCT__
97831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor"
98831a100dSStefano Zampini /*@
99831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
100831a100dSStefano Zampini 
101831a100dSStefano Zampini    Not collective
102831a100dSStefano Zampini 
103831a100dSStefano Zampini    Input Parameters:
104831a100dSStefano Zampini +  pc - the preconditioning context
105831a100dSStefano Zampini -  scal - scaling factor for the subdomain
106831a100dSStefano Zampini 
107831a100dSStefano Zampini    Level: intermediate
108831a100dSStefano Zampini 
109831a100dSStefano Zampini    Notes:
110831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
111831a100dSStefano Zampini 
112831a100dSStefano Zampini .seealso: PCBDDC
113831a100dSStefano Zampini @*/
114831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
115831a100dSStefano Zampini {
116831a100dSStefano Zampini   PetscErrorCode ierr;
117831a100dSStefano Zampini 
118831a100dSStefano Zampini   PetscFunctionBegin;
119831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120831a100dSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr);
121831a100dSStefano Zampini   PetscFunctionReturn(0);
122831a100dSStefano Zampini }
123831a100dSStefano Zampini 
124b4319ba4SBarry Smith 
125b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
126b4319ba4SBarry Smith /*
127b4319ba4SBarry Smith    PCISSetUp -
128b4319ba4SBarry Smith */
129b4319ba4SBarry Smith #undef __FUNCT__
130b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp"
1317087cfbeSBarry Smith PetscErrorCode  PCISSetUp(PC pc)
132b4319ba4SBarry Smith {
133b4319ba4SBarry Smith   PC_IS          *pcis  = (PC_IS*)(pc->data);
134b4319ba4SBarry Smith   Mat_IS         *matis = (Mat_IS*)pc->mat->data;
1356849ba73SBarry Smith   PetscErrorCode ierr;
136ace3abfcSBarry Smith   PetscBool      flg;
137831a100dSStefano Zampini   Vec            counter;
138b4319ba4SBarry Smith 
139b4319ba4SBarry Smith   PetscFunctionBegin;
140251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr);
141ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
142b4319ba4SBarry Smith 
143b4319ba4SBarry Smith   pcis->pure_neumann = matis->pure_neumann;
144b4319ba4SBarry Smith 
145892d8026SStefano Zampini   /* get info on mapping */
146892d8026SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
147892d8026SStefano Zampini   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE;
148892d8026SStefano Zampini 
149892d8026SStefano Zampini   /* Creating local and global index sets for interior and inteface nodes. */
150b4319ba4SBarry Smith   {
15113f74950SBarry Smith     PetscInt    n_I;
15213f74950SBarry Smith     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
153892d8026SStefano Zampini     PetscInt    *array;
154892d8026SStefano Zampini     PetscInt    i,j;
155892d8026SStefano Zampini 
156b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
157892d8026SStefano Zampini     ierr = VecGetSize(matis->x,&pcis->n);CHKERRQ(ierr);
158892d8026SStefano Zampini     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array);CHKERRQ(ierr);
159892d8026SStefano Zampini     ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr);
160*4843afd6SStefano Zampini     for (i=0;i<pcis->n_neigh;i++)
161892d8026SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++)
162892d8026SStefano Zampini           array[pcis->shared[i][j]] += 1;
163892d8026SStefano Zampini 
16413f74950SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);CHKERRQ(ierr);
16513f74950SBarry Smith     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);CHKERRQ(ierr);
166b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
167892d8026SStefano Zampini       if (!array[i]) {
1682fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1692fa5cd67SKarl Rupp         n_I++;
1702fa5cd67SKarl Rupp       } else {
1712fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1722fa5cd67SKarl Rupp         pcis->n_B++;
1732fa5cd67SKarl Rupp       }
174b4319ba4SBarry Smith     }
175b4319ba4SBarry Smith     /* Getting the global numbering */
176b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
177b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
178b4319ba4SBarry Smith     ierr         = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
179b4319ba4SBarry Smith     ierr         = ISLocalToGlobalMappingApply(matis->mapping,n_I,      idx_I_local,idx_I_global);CHKERRQ(ierr);
1802fa5cd67SKarl Rupp 
181b4319ba4SBarry Smith     /* Creating the index sets. */
182bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
183bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
184bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
185bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
1862fa5cd67SKarl Rupp 
187b4319ba4SBarry Smith     /* Freeing memory and restoring arrays */
188b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
189b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
190892d8026SStefano Zampini     ierr = PetscFree(array);CHKERRQ(ierr);
191b4319ba4SBarry Smith   }
192b4319ba4SBarry Smith 
193b4319ba4SBarry Smith   /*
194b4319ba4SBarry Smith     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
195b4319ba4SBarry Smith     is such that interior nodes come first than the interface ones, we have
196b4319ba4SBarry Smith 
197b4319ba4SBarry Smith     [           |      ]
198b4319ba4SBarry Smith     [    A_II   | A_IB ]
199b4319ba4SBarry Smith     A = [           |      ]
200b4319ba4SBarry Smith     [-----------+------]
201b4319ba4SBarry Smith     [    A_BI   | A_BB ]
202b4319ba4SBarry Smith   */
203b4319ba4SBarry Smith 
2044aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr);
2054aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
2064aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
2074aa3045dSJed Brown   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
208b4319ba4SBarry Smith 
209b4319ba4SBarry Smith   /*
210b4319ba4SBarry Smith     Creating work vectors and arrays
211b4319ba4SBarry Smith   */
212892d8026SStefano Zampini   ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
213b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
214b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
215b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
216b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
217b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
218b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
219b4319ba4SBarry Smith   ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
22023ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
221b4319ba4SBarry Smith   ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr);
222b4319ba4SBarry Smith 
223b4319ba4SBarry Smith   /* Creating the scatter contexts */
22423ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
225b4319ba4SBarry Smith   ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
22623ce1328SBarry Smith   ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
227b4319ba4SBarry Smith 
228831a100dSStefano Zampini   /* Creating scaling "matrix" D */
2290298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr);
230831a100dSStefano Zampini   if (!pcis->D) {
231b965d45eSStefano Zampini     ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
232b965d45eSStefano Zampini     if (!pcis->use_stiffness_scaling) {
233b965d45eSStefano Zampini       ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
234ba1573a8SStefano Zampini     } else {
235b965d45eSStefano Zampini       ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr);
236b965d45eSStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
237b965d45eSStefano Zampini       ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
238ba1573a8SStefano Zampini     }
239b965d45eSStefano Zampini   }
240b965d45eSStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
241892d8026SStefano Zampini   ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
242831a100dSStefano Zampini   ierr = VecSet(counter,0.0);CHKERRQ(ierr);
243ba1573a8SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
244ba1573a8SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
245ba1573a8SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
246ba1573a8SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
247ba1573a8SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
248892d8026SStefano Zampini   ierr = VecDestroy(&counter);CHKERRQ(ierr);
249b4319ba4SBarry Smith 
250b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
251b4319ba4SBarry Smith 
252b4319ba4SBarry Smith   /*
253b4319ba4SBarry Smith     Creating the KSP contexts for the local Dirichlet and Neumann problems.
254b4319ba4SBarry Smith   */
255b4319ba4SBarry Smith   {
256b4319ba4SBarry Smith     PC pc_ctx;
257b4319ba4SBarry Smith     /* Dirichlet */
258b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
2591cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
260b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr);
261b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
262b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
263b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
264b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
265b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
266b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
267b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
268b4319ba4SBarry Smith     /* Neumann */
269b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
2701cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
271b4319ba4SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr);
272b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
273b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
274b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
275b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
276b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
277b4319ba4SBarry Smith     {
278ace3abfcSBarry Smith       PetscBool damp_fixed                    = PETSC_FALSE,
27990d69ab7SBarry Smith                 remove_nullspace_fixed        = PETSC_FALSE,
28090d69ab7SBarry Smith                 set_damping_factor_floating   = PETSC_FALSE,
28190d69ab7SBarry Smith                 not_damp_floating             = PETSC_FALSE,
28290d69ab7SBarry Smith                 not_remove_nullspace_floating = PETSC_FALSE;
283b4319ba4SBarry Smith       PetscReal fixed_factor,
284b4319ba4SBarry Smith                 floating_factor;
285b4319ba4SBarry Smith 
2867adad957SLisandro Dalcin       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
2872fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
2880298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr);
289b4319ba4SBarry Smith 
2900298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr);
291b4319ba4SBarry Smith 
2927adad957SLisandro Dalcin       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
293b4319ba4SBarry Smith                               &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
2942fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
2950298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr);
2962fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
297b4319ba4SBarry Smith 
2980298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL);CHKERRQ(ierr);
299b4319ba4SBarry Smith 
3000298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL);CHKERRQ(ierr);
301b4319ba4SBarry Smith 
302b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
303b4319ba4SBarry Smith         if (!(not_damp_floating)) {
304c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
305c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
306b4319ba4SBarry Smith         }
307b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
308b4319ba4SBarry Smith           MatNullSpace nullsp;
3090298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
310d8fd42c4SBarry Smith           ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
311d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
312b4319ba4SBarry Smith         }
313b4319ba4SBarry Smith       } else {  /* fixed subdomain */
314b4319ba4SBarry Smith         if (damp_fixed) {
315c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
316c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
317b4319ba4SBarry Smith         }
318b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
319b4319ba4SBarry Smith           MatNullSpace nullsp;
3200298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
321d8fd42c4SBarry Smith           ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr);
322d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
323b4319ba4SBarry Smith         }
324b4319ba4SBarry Smith       }
325b4319ba4SBarry Smith     }
326b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
327b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
328b4319ba4SBarry Smith   }
329b4319ba4SBarry Smith 
330b4319ba4SBarry Smith   PetscFunctionReturn(0);
331b4319ba4SBarry Smith }
332b4319ba4SBarry Smith 
333b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
334b4319ba4SBarry Smith /*
335b4319ba4SBarry Smith    PCISDestroy -
336b4319ba4SBarry Smith */
337b4319ba4SBarry Smith #undef __FUNCT__
338b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy"
3397087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
340b4319ba4SBarry Smith {
341b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
342dfbe8321SBarry Smith   PetscErrorCode ierr;
343b4319ba4SBarry Smith 
344b4319ba4SBarry Smith   PetscFunctionBegin;
345fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
346fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
347fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
348fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
349fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
350fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
351fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
352fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
353fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
354fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
355fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
356fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
357fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
358fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
359fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
360fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
361fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
362fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
363fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
364fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
365fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
366fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
367fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
36805b42c5fSBarry Smith   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
369b4319ba4SBarry Smith   if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) {
370b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
371b4319ba4SBarry Smith   }
372bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr);
373bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr);
374bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr);
375b4319ba4SBarry Smith   PetscFunctionReturn(0);
376b4319ba4SBarry Smith }
377b4319ba4SBarry Smith 
378b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
379b4319ba4SBarry Smith /*
380b4319ba4SBarry Smith    PCISCreate -
381b4319ba4SBarry Smith */
382b4319ba4SBarry Smith #undef __FUNCT__
383b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate"
3847087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
385b4319ba4SBarry Smith {
386b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
387831a100dSStefano Zampini   PetscErrorCode ierr;
388b4319ba4SBarry Smith 
389b4319ba4SBarry Smith   PetscFunctionBegin;
390b4319ba4SBarry Smith   pcis->is_B_local  = 0;
391b4319ba4SBarry Smith   pcis->is_I_local  = 0;
392b4319ba4SBarry Smith   pcis->is_B_global = 0;
393b4319ba4SBarry Smith   pcis->is_I_global = 0;
394b4319ba4SBarry Smith   pcis->A_II        = 0;
395b4319ba4SBarry Smith   pcis->A_IB        = 0;
396b4319ba4SBarry Smith   pcis->A_BI        = 0;
397b4319ba4SBarry Smith   pcis->A_BB        = 0;
398b4319ba4SBarry Smith   pcis->D           = 0;
399b4319ba4SBarry Smith   pcis->ksp_N       = 0;
400b4319ba4SBarry Smith   pcis->ksp_D      = 0;
401b4319ba4SBarry Smith   pcis->vec1_N      = 0;
402b4319ba4SBarry Smith   pcis->vec2_N      = 0;
403b4319ba4SBarry Smith   pcis->vec1_D      = 0;
404b4319ba4SBarry Smith   pcis->vec2_D      = 0;
405b4319ba4SBarry Smith   pcis->vec3_D      = 0;
406b4319ba4SBarry Smith   pcis->vec1_B      = 0;
407b4319ba4SBarry Smith   pcis->vec2_B      = 0;
408b4319ba4SBarry Smith   pcis->vec3_B      = 0;
409b4319ba4SBarry Smith   pcis->vec1_global = 0;
410b4319ba4SBarry Smith   pcis->work_N      = 0;
411b4319ba4SBarry Smith   pcis->global_to_D = 0;
412b4319ba4SBarry Smith   pcis->N_to_B      = 0;
413b4319ba4SBarry Smith   pcis->global_to_B = 0;
4142fa5cd67SKarl Rupp 
415b4319ba4SBarry Smith   pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE;
4162fa5cd67SKarl Rupp 
417831a100dSStefano Zampini   pcis->scaling_factor = 1.0;
418831a100dSStefano Zampini   /* composing functions */
419bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr);
420bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
421bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr);
422b4319ba4SBarry Smith   PetscFunctionReturn(0);
423b4319ba4SBarry Smith }
424b4319ba4SBarry Smith 
425b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
426b4319ba4SBarry Smith /*
427b4319ba4SBarry Smith    PCISApplySchur -
428b4319ba4SBarry Smith 
429b4319ba4SBarry Smith    Input parameters:
430b4319ba4SBarry Smith .  pc - preconditioner context
431b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
432b4319ba4SBarry Smith 
433b4319ba4SBarry Smith    Output parameters:
434b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
435b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
436b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
437b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
438b4319ba4SBarry Smith 
439b4319ba4SBarry Smith */
440b4319ba4SBarry Smith #undef __FUNCT__
441831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur"
4427087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
443b4319ba4SBarry Smith {
444dfbe8321SBarry Smith   PetscErrorCode ierr;
445b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
446b4319ba4SBarry Smith 
447b4319ba4SBarry Smith   PetscFunctionBegin;
4482fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
449b4319ba4SBarry Smith 
450b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
451b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
45223ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
453b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
454efb30889SBarry Smith   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
455b4319ba4SBarry Smith   PetscFunctionReturn(0);
456b4319ba4SBarry Smith }
457b4319ba4SBarry Smith 
458b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
459b4319ba4SBarry Smith /*
460b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
461b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
462b4319ba4SBarry Smith    mode.
463b4319ba4SBarry Smith 
464b4319ba4SBarry Smith    Input parameters:
465b4319ba4SBarry Smith .  pc - preconditioner context
466b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
467b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
468b4319ba4SBarry Smith 
469b4319ba4SBarry Smith    Output parameter:
470b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
471b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
472b4319ba4SBarry Smith 
473b4319ba4SBarry Smith    Notes:
474b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
475b4319ba4SBarry Smith */
476b4319ba4SBarry Smith #undef __FUNCT__
477b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB"
4787087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
479b4319ba4SBarry Smith {
4805d0c19d7SBarry Smith   PetscInt       i;
4815d0c19d7SBarry Smith   const PetscInt *idex;
4826849ba73SBarry Smith   PetscErrorCode ierr;
483b4319ba4SBarry Smith   PetscScalar    *array_B;
484b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
485b4319ba4SBarry Smith 
486b4319ba4SBarry Smith   PetscFunctionBegin;
487b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
488b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
491b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
4922fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
493b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
4942fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
495b4319ba4SBarry Smith     }
496b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
497b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
4982fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
499b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5002fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
501b4319ba4SBarry Smith     }
502b4319ba4SBarry Smith   }
503b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
504b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
505b4319ba4SBarry Smith   PetscFunctionReturn(0);
506b4319ba4SBarry Smith }
507b4319ba4SBarry Smith 
508b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
509b4319ba4SBarry Smith /*
510b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
511b4319ba4SBarry Smith    More precisely, solves the problem:
512b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
513b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
514b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
515b4319ba4SBarry Smith 
516b4319ba4SBarry Smith    Input parameters:
517b4319ba4SBarry Smith .  pc - preconditioner context
518b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
519b4319ba4SBarry Smith 
520b4319ba4SBarry Smith    Output parameters:
521b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
522b4319ba4SBarry Smith        complement to b
523b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
524b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
525b4319ba4SBarry Smith 
526b4319ba4SBarry Smith */
527b4319ba4SBarry Smith #undef __FUNCT__
528b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur"
5297087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
530b4319ba4SBarry Smith {
531dfbe8321SBarry Smith   PetscErrorCode ierr;
532b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
533b4319ba4SBarry Smith 
534b4319ba4SBarry Smith   PetscFunctionBegin;
535b4319ba4SBarry Smith   /*
536b4319ba4SBarry Smith     Neumann solvers.
537b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
538b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
539b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
540b4319ba4SBarry Smith     is stored in x.
541b4319ba4SBarry Smith   */
542b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
543efb30889SBarry Smith   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
544ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
547b4319ba4SBarry Smith   {
548ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
5490298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr);
550b4319ba4SBarry Smith     if (flg) {
551b4319ba4SBarry Smith       PetscScalar average;
5523050cee2SBarry Smith       PetscViewer viewer;
553ce94432eSBarry Smith       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
5543050cee2SBarry Smith 
555b4319ba4SBarry Smith       ierr    = VecSum(vec1_N,&average);CHKERRQ(ierr);
556b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
5577b23a99aSBarry Smith       ierr    = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
558b4319ba4SBarry Smith       if (pcis->pure_neumann) {
5599f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
560b4319ba4SBarry Smith       } else {
5619f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
562b4319ba4SBarry Smith       }
5633050cee2SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5647b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
565b4319ba4SBarry Smith     }
566b4319ba4SBarry Smith   }
567b4319ba4SBarry Smith   /* Solving the system for vec2_N */
56823ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
569b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
570ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
571ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
572b4319ba4SBarry Smith   PetscFunctionReturn(0);
573b4319ba4SBarry Smith }
574b4319ba4SBarry Smith 
575b4319ba4SBarry Smith 
576b4319ba4SBarry Smith 
577b4319ba4SBarry Smith 
578b4319ba4SBarry Smith 
579b4319ba4SBarry Smith 
580b4319ba4SBarry Smith 
581b4319ba4SBarry Smith 
582b4319ba4SBarry Smith 
583