xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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);
35064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(pc,use,2);
36b965d45eSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
37b965d45eSStefano Zampini   PetscFunctionReturn(0);
38b965d45eSStefano Zampini }
39b965d45eSStefano Zampini 
40ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
41ba1573a8SStefano Zampini {
42ba1573a8SStefano Zampini   PetscErrorCode ierr;
43ba1573a8SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
44ba1573a8SStefano Zampini 
45ba1573a8SStefano Zampini   PetscFunctionBegin;
46ba1573a8SStefano Zampini   ierr    = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr);
47bf83c2afSStefano Zampini   ierr    = VecDestroy(&pcis->D);CHKERRQ(ierr);
48ba1573a8SStefano Zampini   pcis->D = scaling_factors;
49a007db60SStefano Zampini   if (pc->setupcalled) {
50a007db60SStefano Zampini     PetscInt sn;
51a007db60SStefano Zampini 
52a007db60SStefano Zampini     ierr = VecGetSize(pcis->D,&sn);CHKERRQ(ierr);
53a007db60SStefano Zampini     if (sn == pcis->n) {
54a007db60SStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
55a007db60SStefano Zampini       ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
56a007db60SStefano Zampini       ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
57a007db60SStefano Zampini       ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
58a007db60SStefano Zampini       ierr = VecCopy(pcis->vec1_B,pcis->D);CHKERRQ(ierr);
59*2c71b3e2SJacob Faibussowitsch     } else PetscCheckFalse(sn != pcis->n_B,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Invalid size for scaling vector. Expected %D (or full %D), found %D",pcis->n_B,pcis->n,sn);
60a007db60SStefano Zampini   }
61ba1573a8SStefano Zampini   PetscFunctionReturn(0);
62ba1573a8SStefano Zampini }
63ba1573a8SStefano Zampini 
64ba1573a8SStefano Zampini /*@
65ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
66ba1573a8SStefano Zampini 
67ba1573a8SStefano Zampini    Not collective
68ba1573a8SStefano Zampini 
69ba1573a8SStefano Zampini    Input Parameters:
70ba1573a8SStefano Zampini +  pc - the preconditioning context
71ba1573a8SStefano Zampini -  scaling_factors - scaling factors for the subdomain
72ba1573a8SStefano Zampini 
73ba1573a8SStefano Zampini    Level: intermediate
74ba1573a8SStefano Zampini 
75ba1573a8SStefano Zampini    Notes:
76ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
77ba1573a8SStefano Zampini 
78ba1573a8SStefano Zampini .seealso: PCBDDC
79ba1573a8SStefano Zampini @*/
80ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
81ba1573a8SStefano Zampini {
82ba1573a8SStefano Zampini   PetscErrorCode ierr;
83ba1573a8SStefano Zampini 
84ba1573a8SStefano Zampini   PetscFunctionBegin;
85ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86a007db60SStefano Zampini   PetscValidHeaderSpecific(scaling_factors,VEC_CLASSID,2);
87ba1573a8SStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr);
88ba1573a8SStefano Zampini   PetscFunctionReturn(0);
89ba1573a8SStefano Zampini }
90ba1573a8SStefano Zampini 
91831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
92831a100dSStefano Zampini {
93831a100dSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
94831a100dSStefano Zampini 
95831a100dSStefano Zampini   PetscFunctionBegin;
96831a100dSStefano Zampini   pcis->scaling_factor = scal;
97a007db60SStefano Zampini   if (pcis->D) {
98a007db60SStefano Zampini     PetscErrorCode ierr;
99a007db60SStefano Zampini 
100a007db60SStefano Zampini     ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
101a007db60SStefano Zampini   }
102831a100dSStefano Zampini   PetscFunctionReturn(0);
103831a100dSStefano Zampini }
104831a100dSStefano Zampini 
105831a100dSStefano Zampini /*@
106831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
107831a100dSStefano Zampini 
108831a100dSStefano Zampini    Not collective
109831a100dSStefano Zampini 
110831a100dSStefano Zampini    Input Parameters:
111831a100dSStefano Zampini +  pc - the preconditioning context
112831a100dSStefano Zampini -  scal - scaling factor for the subdomain
113831a100dSStefano Zampini 
114831a100dSStefano Zampini    Level: intermediate
115831a100dSStefano Zampini 
116831a100dSStefano Zampini    Notes:
117831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
118831a100dSStefano Zampini 
119831a100dSStefano Zampini .seealso: PCBDDC
120831a100dSStefano Zampini @*/
121831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
122831a100dSStefano Zampini {
123831a100dSStefano Zampini   PetscErrorCode ierr;
124831a100dSStefano Zampini 
125831a100dSStefano Zampini   PetscFunctionBegin;
126831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
127831a100dSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr);
128831a100dSStefano Zampini   PetscFunctionReturn(0);
129831a100dSStefano Zampini }
130831a100dSStefano Zampini 
131b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
132b4319ba4SBarry Smith /*
133b4319ba4SBarry Smith    PCISSetUp -
134b4319ba4SBarry Smith */
135d9869140SStefano Zampini PetscErrorCode  PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
136b4319ba4SBarry Smith {
137b4319ba4SBarry Smith   PC_IS          *pcis  = (PC_IS*)(pc->data);
138bf327c11SStefano Zampini   Mat_IS         *matis;
1395e8657edSStefano Zampini   MatReuse       reuse;
1406849ba73SBarry Smith   PetscErrorCode ierr;
14185c21eb1SStefano Zampini   PetscBool      flg,issbaij;
142b4319ba4SBarry Smith 
143b4319ba4SBarry Smith   PetscFunctionBegin;
14485c21eb1SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr);
145*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires preconditioning matrix of type MATIS");
14685c21eb1SStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1472f37b69bSStefano Zampini   if (pc->useAmat) {
1482f37b69bSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr);
149*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires linear system matrix of type MATIS");
1502f37b69bSStefano Zampini   }
151b4319ba4SBarry Smith 
1525e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1535e8657edSStefano Zampini   if (!pc->setupcalled) {
1545e8657edSStefano Zampini     PetscInt    n_I;
1555e8657edSStefano Zampini     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
1566f516dd7SStefano Zampini     PetscBT     bt;
1575e8657edSStefano Zampini     PetscInt    i,j;
158b4319ba4SBarry Smith 
159892d8026SStefano Zampini     /* get info on mapping */
1603bbff08aSStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr);
1618d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
1623bbff08aSStefano Zampini     pcis->mapping = pc->pmat->rmap->mapping;
1638d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr);
1648d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
165892d8026SStefano Zampini 
166b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1676f516dd7SStefano Zampini     ierr = PetscBTCreate(pcis->n,&bt);CHKERRQ(ierr);
1684843afd6SStefano Zampini     for (i=0;i<pcis->n_neigh;i++)
1696f516dd7SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
1706f516dd7SStefano Zampini         ierr = PetscBTSet(bt,pcis->shared[i][j]);CHKERRQ(ierr);
1716f516dd7SStefano Zampini       }
172892d8026SStefano Zampini 
1735e8657edSStefano Zampini     /* Creating local and global index sets for interior and inteface nodes. */
174785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr);
175785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr);
176b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
1776f516dd7SStefano Zampini       if (!PetscBTLookup(bt,i)) {
1782fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1792fa5cd67SKarl Rupp         n_I++;
1802fa5cd67SKarl Rupp       } else {
1812fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1822fa5cd67SKarl Rupp         pcis->n_B++;
1832fa5cd67SKarl Rupp       }
184b4319ba4SBarry Smith     }
1856f516dd7SStefano Zampini 
186b4319ba4SBarry Smith     /* Getting the global numbering */
187b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
188b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1898d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
1908d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr);
1912fa5cd67SKarl Rupp 
1925e8657edSStefano Zampini     /* Creating the index sets */
193bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
194e850f338SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
195bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
196e850f338SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
1972fa5cd67SKarl Rupp 
1985e8657edSStefano Zampini     /* Freeing memory */
199b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
200b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
2016f516dd7SStefano Zampini     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
202b4319ba4SBarry Smith 
2035e8657edSStefano Zampini     /* Creating work vectors and arrays */
204892d8026SStefano Zampini     ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
205b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
206f08590b7SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&pcis->vec1_D);CHKERRQ(ierr);
207f08590b7SStefano Zampini     ierr = VecSetSizes(pcis->vec1_D,pcis->n-pcis->n_B,PETSC_DECIDE);CHKERRQ(ierr);
208f08590b7SStefano Zampini     ierr = VecSetType(pcis->vec1_D,((PetscObject)pcis->vec1_N)->type_name);CHKERRQ(ierr);
209b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
210b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
211df187020SStefano Zampini     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr);
212f08590b7SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&pcis->vec1_B);CHKERRQ(ierr);
213f08590b7SStefano Zampini     ierr = VecSetSizes(pcis->vec1_B,pcis->n_B,PETSC_DECIDE);CHKERRQ(ierr);
214f08590b7SStefano Zampini     ierr = VecSetType(pcis->vec1_B,((PetscObject)pcis->vec1_N)->type_name);CHKERRQ(ierr);
215b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
216b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
2170a545947SLisandro Dalcin     ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,NULL);CHKERRQ(ierr);
218854ce69bSBarry Smith     ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr);
2195e8657edSStefano Zampini     /* scaling vector */
220bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2215e8657edSStefano Zampini       ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
222bf83c2afSStefano Zampini       ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
223bf83c2afSStefano Zampini     }
224b4319ba4SBarry Smith 
225b4319ba4SBarry Smith     /* Creating the scatter contexts */
2269448b7f1SJunchao Zhang     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr);
2279448b7f1SJunchao Zhang     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
2289448b7f1SJunchao Zhang     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
2299448b7f1SJunchao Zhang     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
230b4319ba4SBarry Smith 
2315e8657edSStefano Zampini     /* map from boundary to local */
2325e8657edSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr);
2335e8657edSStefano Zampini   }
2345e8657edSStefano Zampini 
235a007db60SStefano Zampini   {
236a007db60SStefano Zampini     PetscInt sn;
237a007db60SStefano Zampini 
238a007db60SStefano Zampini     ierr = VecGetSize(pcis->D,&sn);CHKERRQ(ierr);
239a007db60SStefano Zampini     if (sn == pcis->n) {
240a007db60SStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
241a007db60SStefano Zampini       ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
242a007db60SStefano Zampini       ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
243a007db60SStefano Zampini       ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
244a007db60SStefano Zampini       ierr = VecCopy(pcis->vec1_B,pcis->D);CHKERRQ(ierr);
245*2c71b3e2SJacob Faibussowitsch     } else PetscCheckFalse(sn != pcis->n_B,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Invalid size for scaling vector. Expected %D (or full %D), found %D",pcis->n_B,pcis->n,sn);
246a007db60SStefano Zampini   }
247a007db60SStefano Zampini 
2485e8657edSStefano Zampini   /*
2495e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2505e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2515e8657edSStefano Zampini 
2525e8657edSStefano Zampini         [ A_II | A_IB ]
2535e8657edSStefano Zampini     A = [------+------]
2545e8657edSStefano Zampini         [ A_BI | A_BB ]
2555e8657edSStefano Zampini   */
256d9869140SStefano Zampini   if (computematrices) {
2572f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2582f37b69bSStefano Zampini     PetscInt  bs,ibs;
2592f37b69bSStefano Zampini 
2605e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2613975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2625e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2635e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2645e8657edSStefano Zampini       } else {
2655e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2665e8657edSStefano Zampini       }
2675e8657edSStefano Zampini     }
2685e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2695e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
2702f37b69bSStefano Zampini       ierr = MatDestroy(&pcis->pA_II);CHKERRQ(ierr);
2715e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
2725e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
2735e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
2745e8657edSStefano Zampini     }
2755e8657edSStefano Zampini 
2762f37b69bSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockSize(pcis->mapping,&ibs);CHKERRQ(ierr);
2772f37b69bSStefano Zampini     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
2782f37b69bSStefano Zampini     ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->pA_II);CHKERRQ(ierr);
2792f37b69bSStefano Zampini     if (amat) {
2802f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS*)pc->mat->data;
2812f37b69bSStefano Zampini       ierr = MatCreateSubMatrix(amatis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
2822f37b69bSStefano Zampini     } else {
2832f37b69bSStefano Zampini       ierr = PetscObjectReference((PetscObject)pcis->pA_II);CHKERRQ(ierr);
2842f37b69bSStefano Zampini       ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
2852f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2862f37b69bSStefano Zampini     }
2872f37b69bSStefano Zampini     ierr = MatSetBlockSize(pcis->A_II,bs == ibs ? bs : 1);CHKERRQ(ierr);
2882f37b69bSStefano Zampini     ierr = MatSetBlockSize(pcis->pA_II,bs == ibs ? bs : 1);CHKERRQ(ierr);
2897dae84e0SHong Zhang     ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
2905e8657edSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
2915e8657edSStefano Zampini     if (!issbaij) {
2927dae84e0SHong Zhang       ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2937dae84e0SHong Zhang       ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
2945e8657edSStefano Zampini     } else {
2955e8657edSStefano Zampini       Mat newmat;
296d9869140SStefano Zampini 
2975e8657edSStefano Zampini       ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr);
2987dae84e0SHong Zhang       ierr = MatCreateSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2997dae84e0SHong Zhang       ierr = MatCreateSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
3005e8657edSStefano Zampini       ierr = MatDestroy(&newmat);CHKERRQ(ierr);
3015e8657edSStefano Zampini     }
3022f37b69bSStefano Zampini     ierr = MatSetBlockSize(pcis->A_BB,bs == ibs ? bs : 1);CHKERRQ(ierr);
303d9869140SStefano Zampini   }
3045e8657edSStefano Zampini 
305bf83c2afSStefano Zampini   /* Creating scaling vector D */
306c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr);
307bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
30815f235b8SStefano Zampini     PetscScalar *a;
30915f235b8SStefano Zampini     PetscInt    i,n;
31015f235b8SStefano Zampini 
311d9869140SStefano Zampini     if (pcis->A_BB) {
312bf83c2afSStefano Zampini       ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr);
313d9869140SStefano Zampini     } else {
314d9869140SStefano Zampini       ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr);
3158d25b668SStefano Zampini       ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3168d25b668SStefano Zampini       ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
317d9869140SStefano Zampini     }
3188d25b668SStefano Zampini     ierr = VecAbs(pcis->D);CHKERRQ(ierr);
31915f235b8SStefano Zampini     ierr = VecGetLocalSize(pcis->D,&n);CHKERRQ(ierr);
32015f235b8SStefano Zampini     ierr = VecGetArray(pcis->D,&a);CHKERRQ(ierr);
32115f235b8SStefano Zampini     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
32215f235b8SStefano Zampini     ierr = VecRestoreArray(pcis->D,&a);CHKERRQ(ierr);
323ba1573a8SStefano Zampini   }
3248d25b668SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
3258d25b668SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3268d25b668SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3278d25b668SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3288d25b668SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
329ba1573a8SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
330b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
331b4319ba4SBarry Smith 
3325e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
3335e8657edSStefano Zampini   if (computesolvers) {
334b4319ba4SBarry Smith     PC pc_ctx;
3355e8657edSStefano Zampini 
3365e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
337b4319ba4SBarry Smith     /* Dirichlet */
338b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
339422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr);
3401cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
34123ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr);
342b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
343b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
344b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
345b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
346b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
347b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
348b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
349b4319ba4SBarry Smith     /* Neumann */
350b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
351422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr);
3521cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
35323ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr);
354b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
355b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
356b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
357b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
358b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
359b4319ba4SBarry Smith     {
360ace3abfcSBarry Smith       PetscBool damp_fixed                    = PETSC_FALSE,
36190d69ab7SBarry Smith                 remove_nullspace_fixed        = PETSC_FALSE,
36290d69ab7SBarry Smith                 set_damping_factor_floating   = PETSC_FALSE,
36390d69ab7SBarry Smith                 not_damp_floating             = PETSC_FALSE,
36490d69ab7SBarry Smith                 not_remove_nullspace_floating = PETSC_FALSE;
365b4319ba4SBarry Smith       PetscReal fixed_factor,
366b4319ba4SBarry Smith                 floating_factor;
367b4319ba4SBarry Smith 
368c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
3692fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
370c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr);
371b4319ba4SBarry Smith 
372c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr);
373b4319ba4SBarry Smith 
374c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
375b4319ba4SBarry Smith                               &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
3762fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
377c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr);
3782fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
379b4319ba4SBarry Smith 
380c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL);CHKERRQ(ierr);
381b4319ba4SBarry Smith 
382c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL);CHKERRQ(ierr);
383b4319ba4SBarry Smith 
384b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
385b4319ba4SBarry Smith         if (!(not_damp_floating)) {
386c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
387c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
388b4319ba4SBarry Smith         }
389b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
390b4319ba4SBarry Smith           MatNullSpace nullsp;
3910298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
3925fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
393d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
394b4319ba4SBarry Smith         }
395b4319ba4SBarry Smith       } else {  /* fixed subdomain */
396b4319ba4SBarry Smith         if (damp_fixed) {
397c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
398c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
399b4319ba4SBarry Smith         }
400b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
401b4319ba4SBarry Smith           MatNullSpace nullsp;
4020298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
4035fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
404d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
405b4319ba4SBarry Smith         }
406b4319ba4SBarry Smith       }
407b4319ba4SBarry Smith     }
408b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
409b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
410b4319ba4SBarry Smith   }
411b4319ba4SBarry Smith   PetscFunctionReturn(0);
412b4319ba4SBarry Smith }
413b4319ba4SBarry Smith 
414b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
415b4319ba4SBarry Smith /*
416b4319ba4SBarry Smith    PCISDestroy -
417b4319ba4SBarry Smith */
4187087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
419b4319ba4SBarry Smith {
420b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
421dfbe8321SBarry Smith   PetscErrorCode ierr;
422b4319ba4SBarry Smith 
423b4319ba4SBarry Smith   PetscFunctionBegin;
424fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
425fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
426fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
427fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
428fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
4292f37b69bSStefano Zampini   ierr = MatDestroy(&pcis->pA_II);CHKERRQ(ierr);
430fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
431fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
432fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
433fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
434fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
435fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
436fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
437fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
438fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
439fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
440fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
441df187020SStefano Zampini   ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr);
442fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
443fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
444fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
445fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
446fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
447fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
44816909a7fSStefano Zampini   ierr = VecScatterDestroy(&pcis->N_to_D);CHKERRQ(ierr);
449fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
45005b42c5fSBarry Smith   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
4517dbfca69SStefano Zampini   if (pcis->n_neigh > -1) {
4528d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
4537dbfca69SStefano Zampini   }
4548d12c799SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
4555e8657edSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr);
456bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr);
457bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr);
458bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr);
459b4319ba4SBarry Smith   PetscFunctionReturn(0);
460b4319ba4SBarry Smith }
461b4319ba4SBarry Smith 
462b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
463b4319ba4SBarry Smith /*
464b4319ba4SBarry Smith    PCISCreate -
465b4319ba4SBarry Smith */
4667087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
467b4319ba4SBarry Smith {
468b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
469831a100dSStefano Zampini   PetscErrorCode ierr;
470b4319ba4SBarry Smith 
471b4319ba4SBarry Smith   PetscFunctionBegin;
4727dbfca69SStefano Zampini   pcis->n_neigh          = -1;
473831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4743975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
475831a100dSStefano Zampini   /* composing functions */
476bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr);
477bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
478bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr);
479b4319ba4SBarry Smith   PetscFunctionReturn(0);
480b4319ba4SBarry Smith }
481b4319ba4SBarry Smith 
482b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
483b4319ba4SBarry Smith /*
484b4319ba4SBarry Smith    PCISApplySchur -
485b4319ba4SBarry Smith 
486b4319ba4SBarry Smith    Input parameters:
487b4319ba4SBarry Smith .  pc - preconditioner context
488b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith    Output parameters:
491b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
492b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
493b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
494b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
495b4319ba4SBarry Smith 
496b4319ba4SBarry Smith */
4977087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
498b4319ba4SBarry Smith {
499dfbe8321SBarry Smith   PetscErrorCode ierr;
500b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
501b4319ba4SBarry Smith 
502b4319ba4SBarry Smith   PetscFunctionBegin;
5032fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
504b4319ba4SBarry Smith 
505b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
506b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
50723ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
508c0decd05SBarry Smith   ierr = KSPCheckSolve(pcis->ksp_D,pc,vec2_D);CHKERRQ(ierr);
509b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
510efb30889SBarry Smith   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
511b4319ba4SBarry Smith   PetscFunctionReturn(0);
512b4319ba4SBarry Smith }
513b4319ba4SBarry Smith 
514b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
515b4319ba4SBarry Smith /*
516b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
517b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
518b4319ba4SBarry Smith    mode.
519b4319ba4SBarry Smith 
520b4319ba4SBarry Smith    Input parameters:
521b4319ba4SBarry Smith .  pc - preconditioner context
522b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
523b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
524b4319ba4SBarry Smith 
525b4319ba4SBarry Smith    Output parameter:
526b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
527b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
528b4319ba4SBarry Smith 
529b4319ba4SBarry Smith    Notes:
530b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
531b4319ba4SBarry Smith */
5327087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
533b4319ba4SBarry Smith {
5345d0c19d7SBarry Smith   PetscInt       i;
5355d0c19d7SBarry Smith   const PetscInt *idex;
5366849ba73SBarry Smith   PetscErrorCode ierr;
537b4319ba4SBarry Smith   PetscScalar    *array_B;
538b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith   PetscFunctionBegin;
541b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
542b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
543b4319ba4SBarry Smith 
544b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
545b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5462fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
547b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5482fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
549b4319ba4SBarry Smith     }
550b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
551b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5522fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
553b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5542fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
555b4319ba4SBarry Smith     }
556b4319ba4SBarry Smith   }
557b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
558b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
559b4319ba4SBarry Smith   PetscFunctionReturn(0);
560b4319ba4SBarry Smith }
561b4319ba4SBarry Smith 
562b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
563b4319ba4SBarry Smith /*
564b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
565b4319ba4SBarry Smith    More precisely, solves the problem:
566b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
567b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
568b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
569b4319ba4SBarry Smith 
570b4319ba4SBarry Smith    Input parameters:
571b4319ba4SBarry Smith .  pc - preconditioner context
572b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
573b4319ba4SBarry Smith 
574b4319ba4SBarry Smith    Output parameters:
575b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
576b4319ba4SBarry Smith        complement to b
577b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
578b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
579b4319ba4SBarry Smith 
580b4319ba4SBarry Smith */
5817087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
582b4319ba4SBarry Smith {
583dfbe8321SBarry Smith   PetscErrorCode ierr;
584b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
585b4319ba4SBarry Smith 
586b4319ba4SBarry Smith   PetscFunctionBegin;
587b4319ba4SBarry Smith   /*
588b4319ba4SBarry Smith     Neumann solvers.
589b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
590b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
591b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
592b4319ba4SBarry Smith     is stored in x.
593b4319ba4SBarry Smith   */
594b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
595efb30889SBarry Smith   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
596ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
597ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
598b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
599b4319ba4SBarry Smith   {
600ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
601c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr);
602b4319ba4SBarry Smith     if (flg) {
603b4319ba4SBarry Smith       PetscScalar average;
6043050cee2SBarry Smith       PetscViewer viewer;
605ce94432eSBarry Smith       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
6063050cee2SBarry Smith 
607b4319ba4SBarry Smith       ierr    = VecSum(vec1_N,&average);CHKERRQ(ierr);
608b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
6091575c14dSBarry Smith       ierr    = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
610b4319ba4SBarry Smith       if (pcis->pure_neumann) {
6119f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
612b4319ba4SBarry Smith       } else {
6139f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
614b4319ba4SBarry Smith       }
6153050cee2SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6161575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
617b4319ba4SBarry Smith     }
618b4319ba4SBarry Smith   }
619b4319ba4SBarry Smith   /* Solving the system for vec2_N */
62023ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
621c0decd05SBarry Smith   ierr = KSPCheckSolve(pcis->ksp_N,pc,vec2_N);CHKERRQ(ierr);
622b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
623ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
624ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
625b4319ba4SBarry Smith   PetscFunctionReturn(0);
626b4319ba4SBarry Smith }
627