xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
123ce1328SBarry Smith 
25e5bbd0aSStefano Zampini #include <petsc/private/pcisimpl.h> /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
49371c9d4SSatish Balay static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use) {
5b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
6b965d45eSStefano Zampini 
7b965d45eSStefano Zampini   PetscFunctionBegin;
8b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
9b965d45eSStefano Zampini   PetscFunctionReturn(0);
10b965d45eSStefano Zampini }
11b965d45eSStefano Zampini 
12b965d45eSStefano Zampini /*@
13b965d45eSStefano Zampini  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
14b965d45eSStefano Zampini                               local matrices' diagonal.
15b965d45eSStefano Zampini 
16b965d45eSStefano Zampini    Not collective
17b965d45eSStefano Zampini 
18b965d45eSStefano Zampini    Input Parameters:
19b965d45eSStefano Zampini +  pc - the preconditioning context
20b965d45eSStefano Zampini -  use - whether or not pcis use matrix diagonal to build partition of unity.
21b965d45eSStefano Zampini 
22b965d45eSStefano Zampini    Level: intermediate
23b965d45eSStefano Zampini 
24b965d45eSStefano Zampini    Notes:
25b965d45eSStefano Zampini 
26db781477SPatrick Sanan .seealso: `PCBDDC`
27b965d45eSStefano Zampini @*/
289371c9d4SSatish Balay PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use) {
29b965d45eSStefano Zampini   PetscFunctionBegin;
30b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
31064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(pc, use, 2);
32cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetUseStiffnessScaling_C", (PC, PetscBool), (pc, use));
33b965d45eSStefano Zampini   PetscFunctionReturn(0);
34b965d45eSStefano Zampini }
35b965d45eSStefano Zampini 
369371c9d4SSatish Balay static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors) {
37ba1573a8SStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
38ba1573a8SStefano Zampini 
39ba1573a8SStefano Zampini   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)scaling_factors));
419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
42ba1573a8SStefano Zampini   pcis->D = scaling_factors;
43a007db60SStefano Zampini   if (pc->setupcalled) {
44a007db60SStefano Zampini     PetscInt sn;
45a007db60SStefano Zampini 
469566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
47a007db60SStefano Zampini     if (sn == pcis->n) {
489566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
499566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
509566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
519566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
529566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
5363a3b9bcSJacob Faibussowitsch     } else PetscCheck(sn == pcis->n_B, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Invalid size for scaling vector. Expected %" PetscInt_FMT " (or full %" PetscInt_FMT "), found %" PetscInt_FMT, pcis->n_B, pcis->n, sn);
54a007db60SStefano Zampini   }
55ba1573a8SStefano Zampini   PetscFunctionReturn(0);
56ba1573a8SStefano Zampini }
57ba1573a8SStefano Zampini 
58ba1573a8SStefano Zampini /*@
59ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
60ba1573a8SStefano Zampini 
61ba1573a8SStefano Zampini    Not collective
62ba1573a8SStefano Zampini 
63ba1573a8SStefano Zampini    Input Parameters:
64ba1573a8SStefano Zampini +  pc - the preconditioning context
65ba1573a8SStefano Zampini -  scaling_factors - scaling factors for the subdomain
66ba1573a8SStefano Zampini 
67ba1573a8SStefano Zampini    Level: intermediate
68ba1573a8SStefano Zampini 
69ba1573a8SStefano Zampini    Notes:
70ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
71ba1573a8SStefano Zampini 
72db781477SPatrick Sanan .seealso: `PCBDDC`
73ba1573a8SStefano Zampini @*/
749371c9d4SSatish Balay PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors) {
75ba1573a8SStefano Zampini   PetscFunctionBegin;
76ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
77a007db60SStefano Zampini   PetscValidHeaderSpecific(scaling_factors, VEC_CLASSID, 2);
78cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainDiagonalScaling_C", (PC, Vec), (pc, scaling_factors));
79ba1573a8SStefano Zampini   PetscFunctionReturn(0);
80ba1573a8SStefano Zampini }
81ba1573a8SStefano Zampini 
829371c9d4SSatish Balay static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal) {
83831a100dSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
84831a100dSStefano Zampini 
85831a100dSStefano Zampini   PetscFunctionBegin;
86831a100dSStefano Zampini   pcis->scaling_factor = scal;
87*48a46eb9SPierre Jolivet   if (pcis->D) PetscCall(VecSet(pcis->D, pcis->scaling_factor));
88831a100dSStefano Zampini   PetscFunctionReturn(0);
89831a100dSStefano Zampini }
90831a100dSStefano Zampini 
91831a100dSStefano Zampini /*@
92831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
93831a100dSStefano Zampini 
94831a100dSStefano Zampini    Not collective
95831a100dSStefano Zampini 
96831a100dSStefano Zampini    Input Parameters:
97831a100dSStefano Zampini +  pc - the preconditioning context
98831a100dSStefano Zampini -  scal - scaling factor for the subdomain
99831a100dSStefano Zampini 
100831a100dSStefano Zampini    Level: intermediate
101831a100dSStefano Zampini 
102831a100dSStefano Zampini    Notes:
103831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
104831a100dSStefano Zampini 
105db781477SPatrick Sanan .seealso: `PCBDDC`
106831a100dSStefano Zampini @*/
1079371c9d4SSatish Balay PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal) {
108831a100dSStefano Zampini   PetscFunctionBegin;
109831a100dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
110cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainScalingFactor_C", (PC, PetscScalar), (pc, scal));
111831a100dSStefano Zampini   PetscFunctionReturn(0);
112831a100dSStefano Zampini }
113831a100dSStefano Zampini 
114b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
115b4319ba4SBarry Smith /*
116b4319ba4SBarry Smith    PCISSetUp -
117b4319ba4SBarry Smith */
1189371c9d4SSatish Balay PetscErrorCode PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers) {
119b4319ba4SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
120bf327c11SStefano Zampini   Mat_IS   *matis;
1215e8657edSStefano Zampini   MatReuse  reuse;
12285c21eb1SStefano Zampini   PetscBool flg, issbaij;
123b4319ba4SBarry Smith 
124b4319ba4SBarry Smith   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATIS, &flg));
12628b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires preconditioning matrix of type MATIS");
12785c21eb1SStefano Zampini   matis = (Mat_IS *)pc->pmat->data;
1282f37b69bSStefano Zampini   if (pc->useAmat) {
1299566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATIS, &flg));
13028b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires linear system matrix of type MATIS");
1312f37b69bSStefano Zampini   }
132b4319ba4SBarry Smith 
1335e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1345e8657edSStefano Zampini   if (!pc->setupcalled) {
1355e8657edSStefano Zampini     PetscInt  n_I;
1365e8657edSStefano Zampini     PetscInt *idx_I_local, *idx_B_local, *idx_I_global, *idx_B_global;
1376f516dd7SStefano Zampini     PetscBT   bt;
1385e8657edSStefano Zampini     PetscInt  i, j;
139b4319ba4SBarry Smith 
140892d8026SStefano Zampini     /* get info on mapping */
1419566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
1429566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
143e432b41dSStefano Zampini     pcis->mapping = matis->rmapping;
1449566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping, &pcis->n));
1459566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
146892d8026SStefano Zampini 
147b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1489566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(pcis->n, &bt));
1494843afd6SStefano Zampini     for (i = 0; i < pcis->n_neigh; i++)
150*48a46eb9SPierre Jolivet       for (j = 0; j < pcis->n_shared[i]; j++) PetscCall(PetscBTSet(bt, pcis->shared[i][j]));
151892d8026SStefano Zampini 
1525e8657edSStefano Zampini     /* Creating local and global index sets for interior and inteface nodes. */
1539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_I_local));
1549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_B_local));
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;
1689566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, pcis->n_B, idx_B_local, idx_B_global));
1699566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, n_I, idx_I_local, idx_I_global));
1702fa5cd67SKarl Rupp 
1715e8657edSStefano Zampini     /* Creating the index sets */
1729566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pcis->n_B, idx_B_local, PETSC_COPY_VALUES, &pcis->is_B_local));
1739566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), pcis->n_B, idx_B_global, PETSC_COPY_VALUES, &pcis->is_B_global));
1749566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n_I, idx_I_local, PETSC_COPY_VALUES, &pcis->is_I_local));
1759566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), n_I, idx_I_global, PETSC_COPY_VALUES, &pcis->is_I_global));
1762fa5cd67SKarl Rupp 
1775e8657edSStefano Zampini     /* Freeing memory */
1789566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_B_local));
1799566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_I_local));
1809566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&bt));
181b4319ba4SBarry Smith 
1825e8657edSStefano Zampini     /* Creating work vectors and arrays */
1839566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(matis->x, &pcis->vec1_N));
1849566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_N, &pcis->vec2_N));
1859566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_D));
1869566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_D, pcis->n - pcis->n_B, PETSC_DECIDE));
1879566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_D, ((PetscObject)pcis->vec1_N)->type_name));
1889566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec2_D));
1899566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec3_D));
1909566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec4_D));
1919566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_B));
1929566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_B, pcis->n_B, PETSC_DECIDE));
1939566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_B, ((PetscObject)pcis->vec1_N)->type_name));
1949566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec2_B));
1959566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec3_B));
1969566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &pcis->vec1_global, NULL));
1979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &pcis->work_N));
1985e8657edSStefano Zampini     /* scaling vector */
199bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2009566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2019566063dSJacob Faibussowitsch       PetscCall(VecSet(pcis->D, pcis->scaling_factor));
202bf83c2afSStefano Zampini     }
203b4319ba4SBarry Smith 
204b4319ba4SBarry Smith     /* Creating the scatter contexts */
2059566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_I_local, pcis->vec1_D, (IS)0, &pcis->N_to_D));
2069566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_I_global, pcis->vec1_D, (IS)0, &pcis->global_to_D));
2079566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_B_local, pcis->vec1_B, (IS)0, &pcis->N_to_B));
2089566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_B_global, pcis->vec1_B, (IS)0, &pcis->global_to_B));
209b4319ba4SBarry Smith 
2105e8657edSStefano Zampini     /* map from boundary to local */
2119566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local, &pcis->BtoNmap));
2125e8657edSStefano Zampini   }
2135e8657edSStefano Zampini 
214a007db60SStefano Zampini   {
215a007db60SStefano Zampini     PetscInt sn;
216a007db60SStefano Zampini 
2179566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
218a007db60SStefano Zampini     if (sn == pcis->n) {
2199566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2209566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2219566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
2229566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2239566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
22463a3b9bcSJacob Faibussowitsch     } else PetscCheck(sn == pcis->n_B, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Invalid size for scaling vector. Expected %" PetscInt_FMT " (or full %" PetscInt_FMT "), found %" PetscInt_FMT, pcis->n_B, pcis->n, sn);
225a007db60SStefano Zampini   }
226a007db60SStefano Zampini 
2275e8657edSStefano Zampini   /*
2285e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2295e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2305e8657edSStefano Zampini 
2315e8657edSStefano Zampini         [ A_II | A_IB ]
2325e8657edSStefano Zampini     A = [------+------]
2335e8657edSStefano Zampini         [ A_BI | A_BB ]
2345e8657edSStefano Zampini   */
235d9869140SStefano Zampini   if (computematrices) {
2362f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2372f37b69bSStefano Zampini     PetscInt  bs, ibs;
2382f37b69bSStefano Zampini 
2395e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2403975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2415e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2425e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2435e8657edSStefano Zampini       } else {
2445e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2455e8657edSStefano Zampini       }
2465e8657edSStefano Zampini     }
2475e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2489566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2499566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->pA_II));
2509566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_IB));
2519566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BI));
2529566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BB));
2535e8657edSStefano Zampini     }
2545e8657edSStefano Zampini 
2559566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping, &ibs));
2569566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(matis->A, &bs));
2579566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->pA_II));
2582f37b69bSStefano Zampini     if (amat) {
2592f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS *)pc->mat->data;
2609566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(amatis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->A_II));
2612f37b69bSStefano Zampini     } else {
2629566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
2639566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2642f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2652f37b69bSStefano Zampini     }
2669566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_II, bs == ibs ? bs : 1));
2679566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->pA_II, bs == ibs ? bs : 1));
2689566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_B_local, reuse, &pcis->A_BB));
2699566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &issbaij));
2705e8657edSStefano Zampini     if (!issbaij) {
2719566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2729566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2735e8657edSStefano Zampini     } else {
2745e8657edSStefano Zampini       Mat newmat;
275d9869140SStefano Zampini 
2769566063dSJacob Faibussowitsch       PetscCall(MatConvert(matis->A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &newmat));
2779566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2789566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2799566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&newmat));
2805e8657edSStefano Zampini     }
2819566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_BB, bs == ibs ? bs : 1));
282d9869140SStefano Zampini   }
2835e8657edSStefano Zampini 
284bf83c2afSStefano Zampini   /* Creating scaling vector D */
2859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_is_use_stiffness_scaling", &pcis->use_stiffness_scaling, NULL));
286bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
28715f235b8SStefano Zampini     PetscScalar *a;
28815f235b8SStefano Zampini     PetscInt     i, n;
28915f235b8SStefano Zampini 
290d9869140SStefano Zampini     if (pcis->A_BB) {
2919566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pcis->A_BB, pcis->D));
292d9869140SStefano Zampini     } else {
2939566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matis->A, pcis->vec1_N));
2949566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
2959566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
296d9869140SStefano Zampini     }
2979566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
2989566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D, &n));
2999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D, &a));
3009371c9d4SSatish Balay     for (i = 0; i < n; i++)
3019371c9d4SSatish Balay       if (PetscAbsScalar(a[i]) < PETSC_SMALL) a[i] = 1.0;
3029566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D, &a));
303ba1573a8SStefano Zampini   }
3049566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global, 0.0));
3059566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3069566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3079566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3089566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3099566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D, pcis->D, pcis->vec1_B));
310b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
311b4319ba4SBarry Smith 
3125e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
3135e8657edSStefano Zampini   if (computesolvers) {
314b4319ba4SBarry Smith     PC pc_ctx;
3155e8657edSStefano Zampini 
3165e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
317b4319ba4SBarry Smith     /* Dirichlet */
3189566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_D));
3199566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D, pc->erroriffailure));
3209566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D, (PetscObject)pc, 1));
3219566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_D, pcis->A_II, pcis->A_II));
3229566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D, "is_localD_"));
3239566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_D, &pc_ctx));
3249566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3259566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_D, KSPPREONLY));
3269566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_D));
327b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3289566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_D));
329b4319ba4SBarry Smith     /* Neumann */
3309566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_N));
3319566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N, pc->erroriffailure));
3329566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N, (PetscObject)pc, 1));
3339566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_N, matis->A, matis->A));
3349566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N, "is_localN_"));
3359566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_N, &pc_ctx));
3369566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3379566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_N, KSPPREONLY));
3389566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_N));
339b4319ba4SBarry Smith     {
3409371c9d4SSatish Balay       PetscBool damp_fixed = PETSC_FALSE, remove_nullspace_fixed = PETSC_FALSE, set_damping_factor_floating = PETSC_FALSE, not_damp_floating = PETSC_FALSE, not_remove_nullspace_floating = PETSC_FALSE;
3419371c9d4SSatish Balay       PetscReal fixed_factor, floating_factor;
342b4319ba4SBarry Smith 
3439566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &fixed_factor, &damp_fixed));
3442fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
3459566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &damp_fixed, NULL));
346b4319ba4SBarry Smith 
3479566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_remove_nullspace_fixed", &remove_nullspace_fixed, NULL));
348b4319ba4SBarry Smith 
349d0609cedSBarry Smith       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &floating_factor, &set_damping_factor_floating));
3502fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
3519566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &set_damping_factor_floating, NULL));
3522fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
353b4319ba4SBarry Smith 
3549566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_damp_floating", &not_damp_floating, NULL));
355b4319ba4SBarry Smith 
3569566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_remove_nullspace_floating", &not_remove_nullspace_floating, NULL));
357b4319ba4SBarry Smith 
358b4319ba4SBarry Smith       if (pcis->pure_neumann) { /* floating subdomain */
359b4319ba4SBarry Smith         if (!(not_damp_floating)) {
3609566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3619566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
362b4319ba4SBarry Smith         }
363b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
364b4319ba4SBarry Smith           MatNullSpace nullsp;
3659566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3669566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3679566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
368b4319ba4SBarry Smith         }
369b4319ba4SBarry Smith       } else { /* fixed subdomain */
370b4319ba4SBarry Smith         if (damp_fixed) {
3719566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3729566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
373b4319ba4SBarry Smith         }
374b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
375b4319ba4SBarry Smith           MatNullSpace nullsp;
3769566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3779566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3789566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
379b4319ba4SBarry Smith         }
380b4319ba4SBarry Smith       }
381b4319ba4SBarry Smith     }
382b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3839566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_N));
384b4319ba4SBarry Smith   }
385b4319ba4SBarry Smith   PetscFunctionReturn(0);
386b4319ba4SBarry Smith }
387b4319ba4SBarry Smith 
388b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
389b4319ba4SBarry Smith /*
390b4319ba4SBarry Smith    PCISDestroy -
391b4319ba4SBarry Smith */
3929371c9d4SSatish Balay PetscErrorCode PCISDestroy(PC pc) {
393747c8291SStefano Zampini   PC_IS *pcis;
394b4319ba4SBarry Smith 
395b4319ba4SBarry Smith   PetscFunctionBegin;
396747c8291SStefano Zampini   if (!pc) PetscFunctionReturn(0);
397747c8291SStefano Zampini   pcis = (PC_IS *)(pc->data);
3989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_local));
3999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_local));
4009566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_global));
4019566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_global));
4029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_II));
4039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->pA_II));
4049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_IB));
4059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BI));
4069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BB));
4079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
4089566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_N));
4099566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_D));
4109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_N));
4119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_N));
4129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_D));
4139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_D));
4149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_D));
4159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec4_D));
4169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_B));
4179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_B));
4189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_B));
4199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_global));
4209566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_D));
4219566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_B));
4229566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_D));
4239566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_B));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcis->work_N));
425*48a46eb9SPierre Jolivet   if (pcis->n_neigh > -1) PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
4269566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
4279566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", NULL));
4299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", NULL));
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", NULL));
431b4319ba4SBarry Smith   PetscFunctionReturn(0);
432b4319ba4SBarry Smith }
433b4319ba4SBarry Smith 
434b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
435b4319ba4SBarry Smith /*
436b4319ba4SBarry Smith    PCISCreate -
437b4319ba4SBarry Smith */
4389371c9d4SSatish Balay PetscErrorCode PCISCreate(PC pc) {
439b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
440b4319ba4SBarry Smith 
441b4319ba4SBarry Smith   PetscFunctionBegin;
442747c8291SStefano Zampini   if (!pcis) {
443747c8291SStefano Zampini     PetscCall(PetscNewLog(pc, &pcis));
444747c8291SStefano Zampini     pc->data = pcis;
445747c8291SStefano Zampini   }
4467dbfca69SStefano Zampini   pcis->n_neigh          = -1;
447831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4483975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
449831a100dSStefano Zampini   /* composing functions */
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", PCISSetUseStiffnessScaling_IS));
4519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", PCISSetSubdomainScalingFactor_IS));
4529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", PCISSetSubdomainDiagonalScaling_IS));
453b4319ba4SBarry Smith   PetscFunctionReturn(0);
454b4319ba4SBarry Smith }
455b4319ba4SBarry Smith 
456b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
457b4319ba4SBarry Smith /*
458b4319ba4SBarry Smith    PCISApplySchur -
459b4319ba4SBarry Smith 
460b4319ba4SBarry Smith    Input parameters:
461b4319ba4SBarry Smith .  pc - preconditioner context
462b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
463b4319ba4SBarry Smith 
464b4319ba4SBarry Smith    Output parameters:
465b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
466b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
467b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
468b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
469b4319ba4SBarry Smith 
470b4319ba4SBarry Smith */
4719371c9d4SSatish Balay PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) {
472b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
473b4319ba4SBarry Smith 
474b4319ba4SBarry Smith   PetscFunctionBegin;
4752fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
476b4319ba4SBarry Smith 
4779566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BB, v, vec1_B));
4789566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_IB, v, vec1_D));
4799566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_D, vec1_D, vec2_D));
4809566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_D, pc, vec2_D));
4819566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BI, vec2_D, vec2_B));
4829566063dSJacob Faibussowitsch   PetscCall(VecAXPY(vec1_B, -1.0, vec2_B));
483b4319ba4SBarry Smith   PetscFunctionReturn(0);
484b4319ba4SBarry Smith }
485b4319ba4SBarry Smith 
486b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
487b4319ba4SBarry Smith /*
488b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
489b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
490b4319ba4SBarry Smith    mode.
491b4319ba4SBarry Smith 
492b4319ba4SBarry Smith    Input parameters:
493b4319ba4SBarry Smith .  pc - preconditioner context
494b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
495b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
496b4319ba4SBarry Smith 
497b4319ba4SBarry Smith    Output parameter:
498b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
499b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
500b4319ba4SBarry Smith 
501b4319ba4SBarry Smith    Notes:
502b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
503b4319ba4SBarry Smith */
5049371c9d4SSatish Balay PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) {
5055d0c19d7SBarry Smith   PetscInt        i;
5065d0c19d7SBarry Smith   const PetscInt *idex;
507b4319ba4SBarry Smith   PetscScalar    *array_B;
508b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS *)(pc->data);
509b4319ba4SBarry Smith 
510b4319ba4SBarry Smith   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_B, &array_B));
5129566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pcis->is_B_local, &idex));
513b4319ba4SBarry Smith 
514b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
515b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5162fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] = array_N[idex[i]];
517b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5182fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] += array_N[idex[i]];
519b4319ba4SBarry Smith     }
520b4319ba4SBarry Smith   } else { /* SCATTER_REVERSE */
521b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5222fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] = array_B[i];
523b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5242fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] += array_B[i];
525b4319ba4SBarry Smith     }
526b4319ba4SBarry Smith   }
5279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pcis->is_B_local, &idex));
5289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_B, &array_B));
529b4319ba4SBarry Smith   PetscFunctionReturn(0);
530b4319ba4SBarry Smith }
531b4319ba4SBarry Smith 
532b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
533b4319ba4SBarry Smith /*
534b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
535b4319ba4SBarry Smith    More precisely, solves the problem:
536b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
537b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
538b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith    Input parameters:
541b4319ba4SBarry Smith .  pc - preconditioner context
542b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
543b4319ba4SBarry Smith 
544b4319ba4SBarry Smith    Output parameters:
545b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
546b4319ba4SBarry Smith        complement to b
547b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
548b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
549b4319ba4SBarry Smith 
550b4319ba4SBarry Smith */
5519371c9d4SSatish Balay PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) {
552b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
553b4319ba4SBarry Smith 
554b4319ba4SBarry Smith   PetscFunctionBegin;
555b4319ba4SBarry Smith   /*
556b4319ba4SBarry Smith     Neumann solvers.
557b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
558b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
559b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
560b4319ba4SBarry Smith     is stored in x.
561b4319ba4SBarry Smith   */
562b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
5639566063dSJacob Faibussowitsch   PetscCall(VecSet(vec1_N, 0.0));
5649566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
5659566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
566b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
567b4319ba4SBarry Smith   {
568ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
5699566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_is_check_consistency", &flg, NULL));
570b4319ba4SBarry Smith     if (flg) {
571b4319ba4SBarry Smith       PetscScalar average;
5723050cee2SBarry Smith       PetscViewer viewer;
5739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
5743050cee2SBarry Smith 
5759566063dSJacob Faibussowitsch       PetscCall(VecSum(vec1_N, &average));
576b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
5779566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
578b4319ba4SBarry Smith       if (pcis->pure_neumann) {
57963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is floating. Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
580b4319ba4SBarry Smith       } else {
58163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is fixed.    Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
582b4319ba4SBarry Smith       }
5839566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
5849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
585b4319ba4SBarry Smith     }
586b4319ba4SBarry Smith   }
587b4319ba4SBarry Smith   /* Solving the system for vec2_N */
5889566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_N, vec1_N, vec2_N));
5899566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_N, pc, vec2_N));
590b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
5919566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
5929566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
593b4319ba4SBarry Smith   PetscFunctionReturn(0);
594b4319ba4SBarry Smith }
595