xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 15f235b85394a90bb1e957abbd361308ab093171)
123ce1328SBarry Smith 
2aaa7dc30SBarry Smith #include <../src/ksp/pc/impls/is/pcis.h> /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
4b965d45eSStefano Zampini static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
5b965d45eSStefano Zampini {
6b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
7b965d45eSStefano Zampini 
8b965d45eSStefano Zampini   PetscFunctionBegin;
9b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
10b965d45eSStefano Zampini   PetscFunctionReturn(0);
11b965d45eSStefano Zampini }
12b965d45eSStefano Zampini 
13b965d45eSStefano Zampini /*@
14b965d45eSStefano Zampini  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
15b965d45eSStefano Zampini                               local matrices' diagonal.
16b965d45eSStefano Zampini 
17b965d45eSStefano Zampini    Not collective
18b965d45eSStefano Zampini 
19b965d45eSStefano Zampini    Input Parameters:
20b965d45eSStefano Zampini +  pc - the preconditioning context
21b965d45eSStefano Zampini -  use - whether or not pcis use matrix diagonal to build partition of unity.
22b965d45eSStefano Zampini 
23b965d45eSStefano Zampini    Level: intermediate
24b965d45eSStefano Zampini 
25b965d45eSStefano Zampini    Notes:
26b965d45eSStefano Zampini 
27b965d45eSStefano Zampini .seealso: PCBDDC
28b965d45eSStefano Zampini @*/
29b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
30b965d45eSStefano Zampini {
31b965d45eSStefano Zampini   PetscErrorCode ierr;
32b965d45eSStefano Zampini 
33b965d45eSStefano Zampini   PetscFunctionBegin;
34b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
35b965d45eSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
36b965d45eSStefano Zampini   PetscFunctionReturn(0);
37b965d45eSStefano Zampini }
38b965d45eSStefano Zampini 
39ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
40ba1573a8SStefano Zampini {
41ba1573a8SStefano Zampini   PetscErrorCode ierr;
42ba1573a8SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
43ba1573a8SStefano Zampini 
44ba1573a8SStefano Zampini   PetscFunctionBegin;
45ba1573a8SStefano Zampini   ierr    = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr);
46bf83c2afSStefano Zampini   ierr    = VecDestroy(&pcis->D);CHKERRQ(ierr);
47ba1573a8SStefano Zampini   pcis->D = scaling_factors;
48ba1573a8SStefano Zampini   PetscFunctionReturn(0);
49ba1573a8SStefano Zampini }
50ba1573a8SStefano Zampini 
51ba1573a8SStefano Zampini /*@
52ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
53ba1573a8SStefano Zampini 
54ba1573a8SStefano Zampini    Not collective
55ba1573a8SStefano Zampini 
56ba1573a8SStefano Zampini    Input Parameters:
57ba1573a8SStefano Zampini +  pc - the preconditioning context
58ba1573a8SStefano Zampini -  scaling_factors - scaling factors for the subdomain
59ba1573a8SStefano Zampini 
60ba1573a8SStefano Zampini    Level: intermediate
61ba1573a8SStefano Zampini 
62ba1573a8SStefano Zampini    Notes:
63ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
64ba1573a8SStefano Zampini 
65ba1573a8SStefano Zampini .seealso: PCBDDC
66ba1573a8SStefano Zampini @*/
67ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
68ba1573a8SStefano Zampini {
69ba1573a8SStefano Zampini   PetscErrorCode ierr;
70ba1573a8SStefano Zampini 
71ba1573a8SStefano Zampini   PetscFunctionBegin;
72ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
73ba1573a8SStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr);
74ba1573a8SStefano Zampini   PetscFunctionReturn(0);
75ba1573a8SStefano Zampini }
76ba1573a8SStefano Zampini 
77831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
78831a100dSStefano Zampini {
79831a100dSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
80831a100dSStefano Zampini 
81831a100dSStefano Zampini   PetscFunctionBegin;
82831a100dSStefano Zampini   pcis->scaling_factor = scal;
83831a100dSStefano Zampini   PetscFunctionReturn(0);
84831a100dSStefano Zampini }
85831a100dSStefano Zampini 
86831a100dSStefano Zampini /*@
87831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
88831a100dSStefano Zampini 
89831a100dSStefano Zampini    Not collective
90831a100dSStefano Zampini 
91831a100dSStefano Zampini    Input Parameters:
92831a100dSStefano Zampini +  pc - the preconditioning context
93831a100dSStefano Zampini -  scal - scaling factor for the subdomain
94831a100dSStefano Zampini 
95831a100dSStefano Zampini    Level: intermediate
96831a100dSStefano Zampini 
97831a100dSStefano Zampini    Notes:
98831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
99831a100dSStefano Zampini 
100831a100dSStefano Zampini .seealso: PCBDDC
101831a100dSStefano Zampini @*/
102831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
103831a100dSStefano Zampini {
104831a100dSStefano Zampini   PetscErrorCode ierr;
105831a100dSStefano Zampini 
106831a100dSStefano Zampini   PetscFunctionBegin;
107831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
108831a100dSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr);
109831a100dSStefano Zampini   PetscFunctionReturn(0);
110831a100dSStefano Zampini }
111831a100dSStefano Zampini 
112b4319ba4SBarry Smith 
113b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
114b4319ba4SBarry Smith /*
115b4319ba4SBarry Smith    PCISSetUp -
116b4319ba4SBarry Smith */
1175e8657edSStefano Zampini PetscErrorCode  PCISSetUp(PC pc, PetscBool computesolvers)
118b4319ba4SBarry Smith {
119b4319ba4SBarry Smith   PC_IS          *pcis  = (PC_IS*)(pc->data);
120bf327c11SStefano Zampini   Mat_IS         *matis;
1215e8657edSStefano Zampini   MatReuse       reuse;
1226849ba73SBarry Smith   PetscErrorCode ierr;
12385c21eb1SStefano Zampini   PetscBool      flg,issbaij;
124831a100dSStefano Zampini   Vec            counter;
125b4319ba4SBarry Smith 
126b4319ba4SBarry Smith   PetscFunctionBegin;
12785c21eb1SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr);
128ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
12985c21eb1SStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
130b4319ba4SBarry Smith 
1315e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1325e8657edSStefano Zampini   if (!pc->setupcalled) {
1335e8657edSStefano Zampini     PetscInt    n_I;
1345e8657edSStefano Zampini     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
1356f516dd7SStefano Zampini     PetscBT     bt;
1365e8657edSStefano Zampini     PetscInt    i,j;
137b4319ba4SBarry Smith 
138892d8026SStefano Zampini     /* get info on mapping */
1393bbff08aSStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr);
1408d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
1413bbff08aSStefano Zampini     pcis->mapping = pc->pmat->rmap->mapping;
1428d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr);
1438d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
144892d8026SStefano Zampini 
145b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1466f516dd7SStefano Zampini     ierr = PetscBTCreate(pcis->n,&bt);CHKERRQ(ierr);
1474843afd6SStefano Zampini     for (i=0;i<pcis->n_neigh;i++)
1486f516dd7SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
1496f516dd7SStefano Zampini           ierr = PetscBTSet(bt,pcis->shared[i][j]);CHKERRQ(ierr);
1506f516dd7SStefano Zampini       }
151892d8026SStefano Zampini 
1525e8657edSStefano Zampini     /* Creating local and global index sets for interior and inteface nodes. */
153785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr);
154785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr);
155b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
1566f516dd7SStefano Zampini       if (!PetscBTLookup(bt,i)) {
1572fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1582fa5cd67SKarl Rupp         n_I++;
1592fa5cd67SKarl Rupp       } else {
1602fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1612fa5cd67SKarl Rupp         pcis->n_B++;
1622fa5cd67SKarl Rupp       }
163b4319ba4SBarry Smith     }
1646f516dd7SStefano Zampini 
165b4319ba4SBarry Smith     /* Getting the global numbering */
166b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
167b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1688d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
1698d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr);
1702fa5cd67SKarl Rupp 
1715e8657edSStefano Zampini     /* Creating the index sets */
172bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
173bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
174bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
175bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
1762fa5cd67SKarl Rupp 
1775e8657edSStefano Zampini     /* Freeing memory */
178b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
179b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
1806f516dd7SStefano Zampini     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
181b4319ba4SBarry Smith 
1825e8657edSStefano Zampini     /* Creating work vectors and arrays */
183892d8026SStefano Zampini     ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
184b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
185b4319ba4SBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
186b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
187b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
188df187020SStefano Zampini     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr);
189b4319ba4SBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
190b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
191b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
1922a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
193854ce69bSBarry Smith     ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr);
1945e8657edSStefano Zampini     /* scaling vector */
195bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
1965e8657edSStefano Zampini       ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
197bf83c2afSStefano Zampini       ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
198bf83c2afSStefano Zampini     }
199b4319ba4SBarry Smith 
200b4319ba4SBarry Smith     /* Creating the scatter contexts */
20116909a7fSStefano Zampini     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr);
20223ce1328SBarry Smith     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
203b4319ba4SBarry Smith     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
20423ce1328SBarry Smith     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
205b4319ba4SBarry Smith 
2065e8657edSStefano Zampini     /* map from boundary to local */
2075e8657edSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr);
2085e8657edSStefano Zampini   }
2095e8657edSStefano Zampini 
2105e8657edSStefano Zampini   /*
2115e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2125e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2135e8657edSStefano Zampini 
2145e8657edSStefano Zampini         [ A_II | A_IB ]
2155e8657edSStefano Zampini     A = [------+------]
2165e8657edSStefano Zampini         [ A_BI | A_BB ]
2175e8657edSStefano Zampini   */
2185e8657edSStefano Zampini   reuse = MAT_INITIAL_MATRIX;
2193975b054SStefano Zampini   if (pcis->reusesubmatrices && pc->setupcalled) {
2205e8657edSStefano Zampini     if (pc->flag == SAME_NONZERO_PATTERN) {
2215e8657edSStefano Zampini       reuse = MAT_REUSE_MATRIX;
2225e8657edSStefano Zampini     } else {
2235e8657edSStefano Zampini       reuse = MAT_INITIAL_MATRIX;
2245e8657edSStefano Zampini     }
2255e8657edSStefano Zampini   }
2265e8657edSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
2275e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
2285e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
2295e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
2305e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
2315e8657edSStefano Zampini   }
2325e8657edSStefano Zampini 
2337dae84e0SHong Zhang   ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
2347dae84e0SHong Zhang   ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
2355e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
2365e8657edSStefano Zampini   if (!issbaij) {
2377dae84e0SHong Zhang     ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2387dae84e0SHong Zhang     ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
2395e8657edSStefano Zampini   } else {
2405e8657edSStefano Zampini     Mat newmat;
2415e8657edSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr);
2427dae84e0SHong Zhang     ierr = MatCreateSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2437dae84e0SHong Zhang     ierr = MatCreateSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
2445e8657edSStefano Zampini     ierr = MatDestroy(&newmat);CHKERRQ(ierr);
2455e8657edSStefano Zampini   }
2465e8657edSStefano Zampini 
247bf83c2afSStefano Zampini   /* Creating scaling vector D */
248c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr);
249bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
250*15f235b8SStefano Zampini     PetscScalar *a;
251*15f235b8SStefano Zampini     PetscInt    i,n;
252*15f235b8SStefano Zampini 
253bf83c2afSStefano Zampini     ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr);
254*15f235b8SStefano Zampini     ierr = VecGetLocalSize(pcis->D,&n);CHKERRQ(ierr);
255*15f235b8SStefano Zampini     ierr = VecGetArray(pcis->D,&a);CHKERRQ(ierr);
256*15f235b8SStefano Zampini     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
257*15f235b8SStefano Zampini     ierr = VecRestoreArray(pcis->D,&a);CHKERRQ(ierr);
258ba1573a8SStefano Zampini   }
2592a7a6963SBarry Smith   ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
260831a100dSStefano Zampini   ierr = VecSet(counter,0.0);CHKERRQ(ierr);
261bf83c2afSStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
262bf83c2afSStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
263ba1573a8SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
264ba1573a8SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
265ba1573a8SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
266892d8026SStefano Zampini   ierr = VecDestroy(&counter);CHKERRQ(ierr);
267b4319ba4SBarry Smith 
268b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
269b4319ba4SBarry Smith 
2705e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
2715e8657edSStefano Zampini   if (computesolvers) {
272b4319ba4SBarry Smith     PC pc_ctx;
2735e8657edSStefano Zampini 
2745e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
275b4319ba4SBarry Smith     /* Dirichlet */
276b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
277422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr);
2781cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
27923ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr);
280b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
281b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
282b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
283b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
284b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
285b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
286b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
287b4319ba4SBarry Smith     /* Neumann */
288b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
289422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr);
2901cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
29123ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr);
292b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
293b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
294b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
295b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
296b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
297b4319ba4SBarry Smith     {
298ace3abfcSBarry Smith       PetscBool damp_fixed                    = PETSC_FALSE,
29990d69ab7SBarry Smith                 remove_nullspace_fixed        = PETSC_FALSE,
30090d69ab7SBarry Smith                 set_damping_factor_floating   = PETSC_FALSE,
30190d69ab7SBarry Smith                 not_damp_floating             = PETSC_FALSE,
30290d69ab7SBarry Smith                 not_remove_nullspace_floating = PETSC_FALSE;
303b4319ba4SBarry Smith       PetscReal fixed_factor,
304b4319ba4SBarry Smith                 floating_factor;
305b4319ba4SBarry Smith 
306c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
3072fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
308c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr);
309b4319ba4SBarry Smith 
310c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr);
311b4319ba4SBarry Smith 
312c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
313b4319ba4SBarry Smith                               &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
3142fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
315c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr);
3162fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
317b4319ba4SBarry Smith 
318c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL);CHKERRQ(ierr);
319b4319ba4SBarry Smith 
320c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL);CHKERRQ(ierr);
321b4319ba4SBarry Smith 
322b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
323b4319ba4SBarry Smith         if (!(not_damp_floating)) {
324c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
325c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
326b4319ba4SBarry Smith         }
327b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
328b4319ba4SBarry Smith           MatNullSpace nullsp;
3290298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
3305fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
331d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
332b4319ba4SBarry Smith         }
333b4319ba4SBarry Smith       } else {  /* fixed subdomain */
334b4319ba4SBarry Smith         if (damp_fixed) {
335c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
336c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
337b4319ba4SBarry Smith         }
338b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
339b4319ba4SBarry Smith           MatNullSpace nullsp;
3400298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
3415fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
342d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
343b4319ba4SBarry Smith         }
344b4319ba4SBarry Smith       }
345b4319ba4SBarry Smith     }
346b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
347b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
348b4319ba4SBarry Smith   }
349b4319ba4SBarry Smith 
350b4319ba4SBarry Smith   PetscFunctionReturn(0);
351b4319ba4SBarry Smith }
352b4319ba4SBarry Smith 
353b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
354b4319ba4SBarry Smith /*
355b4319ba4SBarry Smith    PCISDestroy -
356b4319ba4SBarry Smith */
3577087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
358b4319ba4SBarry Smith {
359b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
360dfbe8321SBarry Smith   PetscErrorCode ierr;
361b4319ba4SBarry Smith 
362b4319ba4SBarry Smith   PetscFunctionBegin;
363fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
364fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
365fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
366fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
367fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
368fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
369fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
370fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
371fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
372fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
373fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
374fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
375fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
376fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
377fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
378fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
379df187020SStefano Zampini   ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr);
380fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
381fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
382fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
383fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
384fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
385fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
38616909a7fSStefano Zampini   ierr = VecScatterDestroy(&pcis->N_to_D);CHKERRQ(ierr);
387fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
38805b42c5fSBarry Smith   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
3897dbfca69SStefano Zampini   if (pcis->n_neigh > -1) {
3908d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
3917dbfca69SStefano Zampini   }
3928d12c799SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
3935e8657edSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr);
394bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr);
395bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr);
396bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr);
397b4319ba4SBarry Smith   PetscFunctionReturn(0);
398b4319ba4SBarry Smith }
399b4319ba4SBarry Smith 
400b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
401b4319ba4SBarry Smith /*
402b4319ba4SBarry Smith    PCISCreate -
403b4319ba4SBarry Smith */
4047087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
405b4319ba4SBarry Smith {
406b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
407831a100dSStefano Zampini   PetscErrorCode ierr;
408b4319ba4SBarry Smith 
409b4319ba4SBarry Smith   PetscFunctionBegin;
410b4319ba4SBarry Smith   pcis->is_B_local       = 0;
411b4319ba4SBarry Smith   pcis->is_I_local       = 0;
412b4319ba4SBarry Smith   pcis->is_B_global      = 0;
413b4319ba4SBarry Smith   pcis->is_I_global      = 0;
414b4319ba4SBarry Smith   pcis->A_II             = 0;
415b4319ba4SBarry Smith   pcis->A_IB             = 0;
416b4319ba4SBarry Smith   pcis->A_BI             = 0;
417b4319ba4SBarry Smith   pcis->A_BB             = 0;
418b4319ba4SBarry Smith   pcis->D                = 0;
419b4319ba4SBarry Smith   pcis->ksp_N            = 0;
420b4319ba4SBarry Smith   pcis->ksp_D            = 0;
421b4319ba4SBarry Smith   pcis->vec1_N           = 0;
422b4319ba4SBarry Smith   pcis->vec2_N           = 0;
423b4319ba4SBarry Smith   pcis->vec1_D           = 0;
424b4319ba4SBarry Smith   pcis->vec2_D           = 0;
425b4319ba4SBarry Smith   pcis->vec3_D           = 0;
426b4319ba4SBarry Smith   pcis->vec1_B           = 0;
427b4319ba4SBarry Smith   pcis->vec2_B           = 0;
428b4319ba4SBarry Smith   pcis->vec3_B           = 0;
429b4319ba4SBarry Smith   pcis->vec1_global      = 0;
430b4319ba4SBarry Smith   pcis->work_N           = 0;
431b4319ba4SBarry Smith   pcis->global_to_D      = 0;
432b4319ba4SBarry Smith   pcis->N_to_B           = 0;
43316909a7fSStefano Zampini   pcis->N_to_D           = 0;
434b4319ba4SBarry Smith   pcis->global_to_B      = 0;
4358d12c799SStefano Zampini   pcis->mapping          = 0;
4365e8657edSStefano Zampini   pcis->BtoNmap          = 0;
4377dbfca69SStefano Zampini   pcis->n_neigh          = -1;
438831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4393975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
440831a100dSStefano Zampini   /* composing functions */
441bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr);
442bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
443bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr);
444b4319ba4SBarry Smith   PetscFunctionReturn(0);
445b4319ba4SBarry Smith }
446b4319ba4SBarry Smith 
447b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
448b4319ba4SBarry Smith /*
449b4319ba4SBarry Smith    PCISApplySchur -
450b4319ba4SBarry Smith 
451b4319ba4SBarry Smith    Input parameters:
452b4319ba4SBarry Smith .  pc - preconditioner context
453b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
454b4319ba4SBarry Smith 
455b4319ba4SBarry Smith    Output parameters:
456b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
457b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
458b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
459b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
460b4319ba4SBarry Smith 
461b4319ba4SBarry Smith */
4627087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
463b4319ba4SBarry Smith {
464dfbe8321SBarry Smith   PetscErrorCode ierr;
465b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
466b4319ba4SBarry Smith 
467b4319ba4SBarry Smith   PetscFunctionBegin;
4682fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
469b4319ba4SBarry Smith 
470b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
471b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
47223ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
473b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
474efb30889SBarry Smith   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
475b4319ba4SBarry Smith   PetscFunctionReturn(0);
476b4319ba4SBarry Smith }
477b4319ba4SBarry Smith 
478b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
479b4319ba4SBarry Smith /*
480b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
481b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
482b4319ba4SBarry Smith    mode.
483b4319ba4SBarry Smith 
484b4319ba4SBarry Smith    Input parameters:
485b4319ba4SBarry Smith .  pc - preconditioner context
486b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
487b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
488b4319ba4SBarry Smith 
489b4319ba4SBarry Smith    Output parameter:
490b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
491b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
492b4319ba4SBarry Smith 
493b4319ba4SBarry Smith    Notes:
494b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
495b4319ba4SBarry Smith */
4967087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
497b4319ba4SBarry Smith {
4985d0c19d7SBarry Smith   PetscInt       i;
4995d0c19d7SBarry Smith   const PetscInt *idex;
5006849ba73SBarry Smith   PetscErrorCode ierr;
501b4319ba4SBarry Smith   PetscScalar    *array_B;
502b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
503b4319ba4SBarry Smith 
504b4319ba4SBarry Smith   PetscFunctionBegin;
505b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
506b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
507b4319ba4SBarry Smith 
508b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
509b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5102fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
511b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5122fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
513b4319ba4SBarry Smith     }
514b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
515b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5162fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
517b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5182fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
519b4319ba4SBarry Smith     }
520b4319ba4SBarry Smith   }
521b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
522b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
523b4319ba4SBarry Smith   PetscFunctionReturn(0);
524b4319ba4SBarry Smith }
525b4319ba4SBarry Smith 
526b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
527b4319ba4SBarry Smith /*
528b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
529b4319ba4SBarry Smith    More precisely, solves the problem:
530b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
531b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
532b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
533b4319ba4SBarry Smith 
534b4319ba4SBarry Smith    Input parameters:
535b4319ba4SBarry Smith .  pc - preconditioner context
536b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
537b4319ba4SBarry Smith 
538b4319ba4SBarry Smith    Output parameters:
539b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
540b4319ba4SBarry Smith        complement to b
541b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
542b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
543b4319ba4SBarry Smith 
544b4319ba4SBarry Smith */
5457087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
546b4319ba4SBarry Smith {
547dfbe8321SBarry Smith   PetscErrorCode ierr;
548b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
549b4319ba4SBarry Smith 
550b4319ba4SBarry Smith   PetscFunctionBegin;
551b4319ba4SBarry Smith   /*
552b4319ba4SBarry Smith     Neumann solvers.
553b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
554b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
555b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
556b4319ba4SBarry Smith     is stored in x.
557b4319ba4SBarry Smith   */
558b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
559efb30889SBarry Smith   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
560ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
561ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
562b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
563b4319ba4SBarry Smith   {
564ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
565c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr);
566b4319ba4SBarry Smith     if (flg) {
567b4319ba4SBarry Smith       PetscScalar average;
5683050cee2SBarry Smith       PetscViewer viewer;
569ce94432eSBarry Smith       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
5703050cee2SBarry Smith 
571b4319ba4SBarry Smith       ierr    = VecSum(vec1_N,&average);CHKERRQ(ierr);
572b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
5731575c14dSBarry Smith       ierr    = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
574b4319ba4SBarry Smith       if (pcis->pure_neumann) {
5759f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
576b4319ba4SBarry Smith       } else {
5779f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
578b4319ba4SBarry Smith       }
5793050cee2SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5801575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
581b4319ba4SBarry Smith     }
582b4319ba4SBarry Smith   }
583b4319ba4SBarry Smith   /* Solving the system for vec2_N */
58423ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
585b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
586ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
587ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
588b4319ba4SBarry Smith   PetscFunctionReturn(0);
589b4319ba4SBarry Smith }
590b4319ba4SBarry Smith 
591b4319ba4SBarry Smith 
592b4319ba4SBarry Smith 
593b4319ba4SBarry Smith 
594b4319ba4SBarry Smith 
595b4319ba4SBarry Smith 
596b4319ba4SBarry Smith 
597b4319ba4SBarry Smith 
598b4319ba4SBarry Smith 
599