xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision f0d2bb25642116b571c72fa9797d32d784b1f106)
15e5bbd0aSStefano Zampini #include <petsc/private/pcisimpl.h> /*I "petscpc.h" I*/
2831a100dSStefano Zampini 
PCISSetUseStiffnessScaling_IS(PC pc,PetscBool use)3d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
4d71ae5a4SJacob Faibussowitsch {
5b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
6b965d45eSStefano Zampini 
7b965d45eSStefano Zampini   PetscFunctionBegin;
8b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
93ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10b965d45eSStefano Zampini }
11b965d45eSStefano Zampini 
12b965d45eSStefano Zampini /*@
13f1580f4eSBarry Smith   PCISSetUseStiffnessScaling - Tells `PCIS` to construct partition of unity using
14f1580f4eSBarry Smith   the local matrices' diagonal entries
15b965d45eSStefano Zampini 
16f1580f4eSBarry Smith   Logically Collective
17b965d45eSStefano Zampini 
18b965d45eSStefano Zampini   Input Parameters:
19b965d45eSStefano Zampini + pc  - the preconditioning context
2020f4b53cSBarry Smith - use - whether or not it should use matrix diagonal to build partition of unity.
21b965d45eSStefano Zampini 
22b965d45eSStefano Zampini   Level: intermediate
23b965d45eSStefano Zampini 
24562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
2504c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
2604c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
27b965d45eSStefano Zampini @*/
PCISSetUseStiffnessScaling(PC pc,PetscBool use)28d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
29d71ae5a4SJacob Faibussowitsch {
30b965d45eSStefano Zampini   PetscFunctionBegin;
31b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
32064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(pc, use, 2);
33cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetUseStiffnessScaling_C", (PC, PetscBool), (pc, use));
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35b965d45eSStefano Zampini }
36b965d45eSStefano Zampini 
PCISSetSubdomainDiagonalScaling_IS(PC pc,Vec scaling_factors)37d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
38d71ae5a4SJacob Faibussowitsch {
39ba1573a8SStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
40ba1573a8SStefano Zampini 
41ba1573a8SStefano Zampini   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)scaling_factors));
439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
44ba1573a8SStefano Zampini   pcis->D = scaling_factors;
45a007db60SStefano Zampini   if (pc->setupcalled) {
46a007db60SStefano Zampini     PetscInt sn;
47a007db60SStefano Zampini 
489566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
49a007db60SStefano Zampini     if (sn == pcis->n) {
509566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
519566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
529566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
539566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
549566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
5563a3b9bcSJacob 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);
56a007db60SStefano Zampini   }
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58ba1573a8SStefano Zampini }
59ba1573a8SStefano Zampini 
60ba1573a8SStefano Zampini /*@
61f1580f4eSBarry Smith   PCISSetSubdomainDiagonalScaling - Set diagonal scaling for `PCIS`.
62ba1573a8SStefano Zampini 
63f1580f4eSBarry Smith   Logically Collective
64ba1573a8SStefano Zampini 
65ba1573a8SStefano Zampini   Input Parameters:
66ba1573a8SStefano Zampini + pc              - the preconditioning context
67ba1573a8SStefano Zampini - scaling_factors - scaling factors for the subdomain
68ba1573a8SStefano Zampini 
69ba1573a8SStefano Zampini   Level: intermediate
70ba1573a8SStefano Zampini 
71f1580f4eSBarry Smith   Note:
72f1580f4eSBarry Smith   Intended for use with jumping coefficients cases.
73ba1573a8SStefano Zampini 
74562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
7504c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISSetUseStiffnessScaling()`,
7604c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
77ba1573a8SStefano Zampini @*/
PCISSetSubdomainDiagonalScaling(PC pc,Vec scaling_factors)78d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
79d71ae5a4SJacob Faibussowitsch {
80ba1573a8SStefano Zampini   PetscFunctionBegin;
81ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
82a007db60SStefano Zampini   PetscValidHeaderSpecific(scaling_factors, VEC_CLASSID, 2);
83cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainDiagonalScaling_C", (PC, Vec), (pc, scaling_factors));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85ba1573a8SStefano Zampini }
86ba1573a8SStefano Zampini 
PCISSetSubdomainScalingFactor_IS(PC pc,PetscScalar scal)87d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
88d71ae5a4SJacob Faibussowitsch {
89831a100dSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
90831a100dSStefano Zampini 
91831a100dSStefano Zampini   PetscFunctionBegin;
92831a100dSStefano Zampini   pcis->scaling_factor = scal;
9348a46eb9SPierre Jolivet   if (pcis->D) PetscCall(VecSet(pcis->D, pcis->scaling_factor));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95831a100dSStefano Zampini }
96831a100dSStefano Zampini 
97831a100dSStefano Zampini /*@
98f1580f4eSBarry Smith   PCISSetSubdomainScalingFactor - Set scaling factor for `PCIS`.
99831a100dSStefano Zampini 
10020f4b53cSBarry Smith   Not Collective
101831a100dSStefano Zampini 
102831a100dSStefano Zampini   Input Parameters:
103831a100dSStefano Zampini + pc   - the preconditioning context
104831a100dSStefano Zampini - scal - scaling factor for the subdomain
105831a100dSStefano Zampini 
106831a100dSStefano Zampini   Level: intermediate
107831a100dSStefano Zampini 
108f1580f4eSBarry Smith   Note:
109f1580f4eSBarry Smith   Intended for use with the jumping coefficients cases.
110831a100dSStefano Zampini 
111562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
11204c3f3b8SBarry Smith           `PCISSetSubdomainDiagonalScaling()`, `PCISSetUseStiffnessScaling()`,
11304c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
114831a100dSStefano Zampini @*/
PCISSetSubdomainScalingFactor(PC pc,PetscScalar scal)115d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
116d71ae5a4SJacob Faibussowitsch {
117831a100dSStefano Zampini   PetscFunctionBegin;
118831a100dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
119cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainScalingFactor_C", (PC, PetscScalar), (pc, scal));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121831a100dSStefano Zampini }
122831a100dSStefano Zampini 
12304c3f3b8SBarry Smith /*@
12404c3f3b8SBarry Smith   PCISSetUp - sets up the `PC_IS` portion of `PCNN` and `PCBDDC` preconditioner context as part of their setup process
12504c3f3b8SBarry Smith 
12604c3f3b8SBarry Smith   Input Parameters:
12704c3f3b8SBarry Smith + pc              - the `PC` object, must be of type `PCNN` or `PCBDDC`
12804c3f3b8SBarry Smith . computematrices - Extract the blocks `A_II`, `A_BI`, `A_IB` and `A_BB` from the matrix
12904c3f3b8SBarry Smith - computesolvers  - Create the `KSP` for the local Dirichlet and Neumann problems
13004c3f3b8SBarry Smith 
13104c3f3b8SBarry Smith   Level: advanced
13204c3f3b8SBarry Smith 
133562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
13404c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
13504c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
13604c3f3b8SBarry Smith @*/
PCISSetUp(PC pc,PetscBool computematrices,PetscBool computesolvers)137d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
138d71ae5a4SJacob Faibussowitsch {
139f4f49eeaSPierre Jolivet   PC_IS    *pcis = (PC_IS *)pc->data;
140bf327c11SStefano Zampini   Mat_IS   *matis;
1415e8657edSStefano Zampini   MatReuse  reuse;
14285c21eb1SStefano Zampini   PetscBool flg, issbaij;
143b4319ba4SBarry Smith 
144b4319ba4SBarry Smith   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATIS, &flg));
146*7addb90fSBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires matrix for constructing the preconditioner of type MATIS");
14785c21eb1SStefano Zampini   matis = (Mat_IS *)pc->pmat->data;
1482f37b69bSStefano Zampini   if (pc->useAmat) {
1499566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATIS, &flg));
15028b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires linear system matrix of type MATIS");
1512f37b69bSStefano Zampini   }
152b4319ba4SBarry Smith 
1535e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1545e8657edSStefano Zampini   if (!pc->setupcalled) {
1555e8657edSStefano Zampini     PetscInt  n_I;
1561ca0b81aSStefano Zampini     PetscInt *idx_I_local, *idx_B_local, *idx_I_global, *idx_B_global, *count;
1571ca0b81aSStefano Zampini     PetscInt  i;
158b4319ba4SBarry Smith 
159892d8026SStefano Zampini     /* get info on mapping */
1609566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
1619566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
162e432b41dSStefano Zampini     pcis->mapping = matis->rmapping;
1639566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping, &pcis->n));
164f4f49eeaSPierre Jolivet     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping, &pcis->n_neigh, &pcis->neigh, &pcis->n_shared, &pcis->shared));
165892d8026SStefano Zampini 
166b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1671ca0b81aSStefano Zampini     PetscCall(ISLocalToGlobalMappingGetNodeInfo(pcis->mapping, NULL, &count, NULL));
1689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_I_local));
1699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_B_local));
170b4319ba4SBarry Smith     for (i = 0, pcis->n_B = 0, n_I = 0; i < pcis->n; i++) {
1711ca0b81aSStefano Zampini       if (count[i] < 2) {
1722fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1732fa5cd67SKarl Rupp         n_I++;
1742fa5cd67SKarl Rupp       } else {
1752fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1762fa5cd67SKarl Rupp         pcis->n_B++;
1772fa5cd67SKarl Rupp       }
178b4319ba4SBarry Smith     }
1791ca0b81aSStefano Zampini     PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(pcis->mapping, NULL, &count, NULL));
1806f516dd7SStefano Zampini 
181b4319ba4SBarry Smith     /* Getting the global numbering */
1828e3a54c0SPierre Jolivet     idx_B_global = PetscSafePointerPlusOffset(idx_I_local, n_I); /* Just avoiding allocating extra memory, since we have vacant space */
1838e3a54c0SPierre Jolivet     idx_I_global = PetscSafePointerPlusOffset(idx_B_local, pcis->n_B);
1849566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, pcis->n_B, idx_B_local, idx_B_global));
1859566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, n_I, idx_I_local, idx_I_global));
1862fa5cd67SKarl Rupp 
1875e8657edSStefano Zampini     /* Creating the index sets */
1889566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pcis->n_B, idx_B_local, PETSC_COPY_VALUES, &pcis->is_B_local));
1899566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), pcis->n_B, idx_B_global, PETSC_COPY_VALUES, &pcis->is_B_global));
1909566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n_I, idx_I_local, PETSC_COPY_VALUES, &pcis->is_I_local));
1919566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), n_I, idx_I_global, PETSC_COPY_VALUES, &pcis->is_I_global));
1922fa5cd67SKarl Rupp 
1935e8657edSStefano Zampini     /* Freeing memory */
1949566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_B_local));
1959566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_I_local));
196b4319ba4SBarry Smith 
1975e8657edSStefano Zampini     /* Creating work vectors and arrays */
1989566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(matis->x, &pcis->vec1_N));
1999566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_N, &pcis->vec2_N));
2009566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_D));
2019566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_D, pcis->n - pcis->n_B, PETSC_DECIDE));
2029566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_D, ((PetscObject)pcis->vec1_N)->type_name));
2039566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec2_D));
2049566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec3_D));
2059566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec4_D));
2069566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_B));
2079566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_B, pcis->n_B, PETSC_DECIDE));
2089566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_B, ((PetscObject)pcis->vec1_N)->type_name));
2099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec2_B));
2109566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec3_B));
2119566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &pcis->vec1_global, NULL));
2129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &pcis->work_N));
2135e8657edSStefano Zampini     /* scaling vector */
214bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2159566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2169566063dSJacob Faibussowitsch       PetscCall(VecSet(pcis->D, pcis->scaling_factor));
217bf83c2afSStefano Zampini     }
218b4319ba4SBarry Smith 
219b4319ba4SBarry Smith     /* Creating the scatter contexts */
2209566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_I_local, pcis->vec1_D, (IS)0, &pcis->N_to_D));
2219566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_I_global, pcis->vec1_D, (IS)0, &pcis->global_to_D));
2229566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_B_local, pcis->vec1_B, (IS)0, &pcis->N_to_B));
2239566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_B_global, pcis->vec1_B, (IS)0, &pcis->global_to_B));
224b4319ba4SBarry Smith 
2255e8657edSStefano Zampini     /* map from boundary to local */
2269566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local, &pcis->BtoNmap));
2275e8657edSStefano Zampini   }
2285e8657edSStefano Zampini 
229a007db60SStefano Zampini   {
230a007db60SStefano Zampini     PetscInt sn;
231a007db60SStefano Zampini 
2329566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
233a007db60SStefano Zampini     if (sn == pcis->n) {
2349566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2359566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2369566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
2379566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2389566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
23963a3b9bcSJacob 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);
240a007db60SStefano Zampini   }
241a007db60SStefano Zampini 
2425e8657edSStefano Zampini   /*
2435e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2445e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2455e8657edSStefano Zampini 
2465e8657edSStefano Zampini         [ A_II | A_IB ]
2475e8657edSStefano Zampini     A = [------+------]
2485e8657edSStefano Zampini         [ A_BI | A_BB ]
2495e8657edSStefano Zampini   */
250d9869140SStefano Zampini   if (computematrices) {
2512f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2522f37b69bSStefano Zampini     PetscInt  bs, ibs;
2532f37b69bSStefano Zampini 
2545e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2553975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2565e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2575e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2585e8657edSStefano Zampini       } else {
2595e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2605e8657edSStefano Zampini       }
2615e8657edSStefano Zampini     }
2625e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2639566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2649566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->pA_II));
2659566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_IB));
2669566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BI));
2679566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BB));
2685e8657edSStefano Zampini     }
2695e8657edSStefano Zampini 
2709566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping, &ibs));
2719566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(matis->A, &bs));
2729566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->pA_II));
2732f37b69bSStefano Zampini     if (amat) {
2742f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS *)pc->mat->data;
2759566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(amatis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->A_II));
2762f37b69bSStefano Zampini     } else {
2779566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
2789566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2792f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2802f37b69bSStefano Zampini     }
2819566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_II, bs == ibs ? bs : 1));
2829566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->pA_II, bs == ibs ? bs : 1));
2839566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_B_local, reuse, &pcis->A_BB));
2849566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &issbaij));
2855e8657edSStefano Zampini     if (!issbaij) {
2869566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2879566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2885e8657edSStefano Zampini     } else {
2895e8657edSStefano Zampini       Mat newmat;
290d9869140SStefano Zampini 
2919566063dSJacob Faibussowitsch       PetscCall(MatConvert(matis->A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &newmat));
2929566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2939566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2949566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&newmat));
2955e8657edSStefano Zampini     }
2969566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_BB, bs == ibs ? bs : 1));
297d9869140SStefano Zampini   }
2985e8657edSStefano Zampini 
299bf83c2afSStefano Zampini   /* Creating scaling vector D */
3009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_is_use_stiffness_scaling", &pcis->use_stiffness_scaling, NULL));
301bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
30215f235b8SStefano Zampini     PetscScalar *a;
30315f235b8SStefano Zampini     PetscInt     i, n;
30415f235b8SStefano Zampini 
305d9869140SStefano Zampini     if (pcis->A_BB) {
3069566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pcis->A_BB, pcis->D));
307d9869140SStefano Zampini     } else {
3089566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matis->A, pcis->vec1_N));
3099566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
3109566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
311d9869140SStefano Zampini     }
3129566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
3139566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D, &n));
3149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D, &a));
3159371c9d4SSatish Balay     for (i = 0; i < n; i++)
3169371c9d4SSatish Balay       if (PetscAbsScalar(a[i]) < PETSC_SMALL) a[i] = 1.0;
3179566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D, &a));
318ba1573a8SStefano Zampini   }
3199566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global, 0.0));
3209566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3219566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3229566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3239566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3249566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D, pcis->D, pcis->vec1_B));
325b4319ba4SBarry Smith 
3265e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
3275e8657edSStefano Zampini   if (computesolvers) {
328b4319ba4SBarry Smith     PC pc_ctx;
3295e8657edSStefano Zampini 
3305e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
331b4319ba4SBarry Smith     /* Dirichlet */
3329566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_D));
3333821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_D, pc->kspnestlevel));
3349566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D, pc->erroriffailure));
3359566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D, (PetscObject)pc, 1));
3369566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_D, pcis->A_II, pcis->A_II));
3379566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D, "is_localD_"));
3389566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_D, &pc_ctx));
3399566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3409566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_D, KSPPREONLY));
3419566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_D));
342b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3439566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_D));
344b4319ba4SBarry Smith     /* Neumann */
3459566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_N));
3463821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_N, pc->kspnestlevel));
3479566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N, pc->erroriffailure));
3489566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N, (PetscObject)pc, 1));
3499566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_N, matis->A, matis->A));
3509566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N, "is_localN_"));
3519566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_N, &pc_ctx));
3529566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3539566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_N, KSPPREONLY));
3549566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_N));
355b4319ba4SBarry Smith     {
3569371c9d4SSatish 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;
3579371c9d4SSatish Balay       PetscReal fixed_factor, floating_factor;
358b4319ba4SBarry Smith 
3599566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &fixed_factor, &damp_fixed));
3602fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
3619566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &damp_fixed, NULL));
362b4319ba4SBarry Smith 
3639566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_remove_nullspace_fixed", &remove_nullspace_fixed, NULL));
364b4319ba4SBarry Smith 
365d0609cedSBarry Smith       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &floating_factor, &set_damping_factor_floating));
3662fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
3679566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &set_damping_factor_floating, NULL));
3682fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
369b4319ba4SBarry Smith 
3709566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_damp_floating", &not_damp_floating, NULL));
371b4319ba4SBarry Smith 
3729566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_remove_nullspace_floating", &not_remove_nullspace_floating, NULL));
373b4319ba4SBarry Smith 
374b4319ba4SBarry Smith       if (pcis->pure_neumann) { /* floating subdomain */
375b4319ba4SBarry Smith         if (!(not_damp_floating)) {
3769566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3779566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
378b4319ba4SBarry Smith         }
379b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
380b4319ba4SBarry Smith           MatNullSpace nullsp;
3819566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3829566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3839566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
384b4319ba4SBarry Smith         }
385b4319ba4SBarry Smith       } else { /* fixed subdomain */
386b4319ba4SBarry Smith         if (damp_fixed) {
3879566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3889566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
389b4319ba4SBarry Smith         }
390b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
391b4319ba4SBarry Smith           MatNullSpace nullsp;
3929566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3939566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3949566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
395b4319ba4SBarry Smith         }
396b4319ba4SBarry Smith       }
397b4319ba4SBarry Smith     }
398b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3999566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_N));
400b4319ba4SBarry Smith   }
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
402b4319ba4SBarry Smith }
403b4319ba4SBarry Smith 
40404c3f3b8SBarry Smith /*@
40504c3f3b8SBarry Smith   PCISReset - Removes all the `PC_IS` parts of the `PC` implementation data structure
40604c3f3b8SBarry Smith 
40704c3f3b8SBarry Smith   Input Parameter:
40804c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
40904c3f3b8SBarry Smith 
41004c3f3b8SBarry Smith   Level: advanced
41104c3f3b8SBarry Smith 
412562efe2eSBarry Smith .seealso: [](ch_ksp), `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`, `PCISSetSubdomainScalingFactor()`,
41304c3f3b8SBarry Smith           `PCISInitialize()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
41404c3f3b8SBarry Smith @*/
PCISReset(PC pc)41504c3f3b8SBarry Smith PetscErrorCode PCISReset(PC pc)
416d71ae5a4SJacob Faibussowitsch {
417f4f49eeaSPierre Jolivet   PC_IS    *pcis = (PC_IS *)pc->data;
41804c3f3b8SBarry Smith   PetscBool correcttype;
419b4319ba4SBarry Smith 
420b4319ba4SBarry Smith   PetscFunctionBegin;
4213ba16761SJacob Faibussowitsch   if (!pc) PetscFunctionReturn(PETSC_SUCCESS);
42204c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
42304c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_local));
4259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_local));
4269566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_global));
4279566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_global));
4289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_II));
4299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->pA_II));
4309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_IB));
4319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BI));
4329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BB));
4339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
4349566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_N));
4359566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_D));
4369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_N));
4379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_N));
4389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_D));
4399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_D));
4409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_D));
4419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec4_D));
4429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_B));
4439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_B));
4449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_B));
4459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_global));
4469566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_D));
4479566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_B));
4489566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_D));
4499566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_B));
4509566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcis->work_N));
451f4f49eeaSPierre Jolivet   if (pcis->n_neigh > -1) PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping, &pcis->n_neigh, &pcis->neigh, &pcis->n_shared, &pcis->shared));
4529566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
4539566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
4549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", NULL));
4559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", NULL));
4569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", NULL));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458b4319ba4SBarry Smith }
459b4319ba4SBarry Smith 
46004c3f3b8SBarry Smith /*@
46104c3f3b8SBarry Smith   PCISInitialize - initializes the `PC_IS` portion of `PCNN` and `PCBDDC` preconditioner context
46204c3f3b8SBarry Smith 
46304c3f3b8SBarry Smith   Input Parameter:
46404c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
46504c3f3b8SBarry Smith 
46604c3f3b8SBarry Smith   Level: advanced
46704c3f3b8SBarry Smith 
46804c3f3b8SBarry Smith   Note:
46904c3f3b8SBarry Smith   There is no preconditioner the `PCIS` prefixed routines provide functionality needed by `PCNN` or `PCBDDC`
47004c3f3b8SBarry Smith 
471562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
47204c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
47304c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
47404c3f3b8SBarry Smith @*/
PCISInitialize(PC pc)47504c3f3b8SBarry Smith PetscErrorCode PCISInitialize(PC pc)
476d71ae5a4SJacob Faibussowitsch {
477f4f49eeaSPierre Jolivet   PC_IS    *pcis = (PC_IS *)pc->data;
47804c3f3b8SBarry Smith   PetscBool correcttype;
479b4319ba4SBarry Smith 
480b4319ba4SBarry Smith   PetscFunctionBegin;
48104c3f3b8SBarry Smith   PetscCheck(pcis, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC_IS context must be created by caller");
48204c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
48304c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4847dbfca69SStefano Zampini   pcis->n_neigh          = -1;
485831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4863975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
4879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", PCISSetUseStiffnessScaling_IS));
4889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", PCISSetSubdomainScalingFactor_IS));
4899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", PCISSetSubdomainDiagonalScaling_IS));
4903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
491b4319ba4SBarry Smith }
492b4319ba4SBarry Smith 
49304c3f3b8SBarry Smith /*@
49404c3f3b8SBarry Smith   PCISApplySchur - applies the Schur complement arising from the `MATIS` inside the `PCNN` preconditioner
495b4319ba4SBarry Smith 
49604c3f3b8SBarry Smith   Input Parameters:
49704c3f3b8SBarry Smith + pc     - preconditioner context
498b4319ba4SBarry Smith . v      - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
49904c3f3b8SBarry Smith . vec1_B - location to store the result of Schur complement applied to chunk
50004c3f3b8SBarry Smith . vec2_B - workspace or `NULL`, `v` is used as workspace in that case
50104c3f3b8SBarry Smith . vec1_D - work space
50204c3f3b8SBarry Smith - vec2_D - work space
503b4319ba4SBarry Smith 
50404c3f3b8SBarry Smith   Level: advanced
505b4319ba4SBarry Smith 
506562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
50704c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplyInvSchur()`,
50804c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
50904c3f3b8SBarry Smith @*/
PCISApplySchur(PC pc,Vec v,Vec vec1_B,Vec vec2_B,Vec vec1_D,Vec vec2_D)510d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
511d71ae5a4SJacob Faibussowitsch {
512f4f49eeaSPierre Jolivet   PC_IS *pcis = (PC_IS *)pc->data;
513b4319ba4SBarry Smith 
514b4319ba4SBarry Smith   PetscFunctionBegin;
5152fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
516b4319ba4SBarry Smith 
5179566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BB, v, vec1_B));
5189566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_IB, v, vec1_D));
5199566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_D, vec1_D, vec2_D));
5209566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_D, pc, vec2_D));
5219566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BI, vec2_D, vec2_B));
5229566063dSJacob Faibussowitsch   PetscCall(VecAXPY(vec1_B, -1.0, vec2_B));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
524b4319ba4SBarry Smith }
525b4319ba4SBarry Smith 
52604c3f3b8SBarry Smith /*@
527b4319ba4SBarry Smith   PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
52804c3f3b8SBarry Smith   including ghosts) into an interface vector, when in `SCATTER_FORWARD` mode, or vice-versa, when in `SCATTER_REVERSE`
529b4319ba4SBarry Smith   mode.
530b4319ba4SBarry Smith 
531feefa0e1SJacob Faibussowitsch   Input Parameters:
532feefa0e1SJacob Faibussowitsch + pc      - preconditioner context
53304c3f3b8SBarry Smith . array_N - [when in `SCATTER_FORWARD` mode] Array to be scattered into the vector otherwise output array
53404c3f3b8SBarry Smith . imode   - insert mode, `ADD_VALUES` or `INSERT_VALUES`
53504c3f3b8SBarry Smith . smode   - scatter mode, `SCATTER_FORWARD` or `SCATTER_REVERSE` mode]
53604c3f3b8SBarry Smith - v_B     - [when in `SCATTER_REVERSE` mode] Vector to be scattered into the array, otherwise output vector
537b4319ba4SBarry Smith 
53804c3f3b8SBarry Smith   Level: advanced
539b4319ba4SBarry Smith 
540f1580f4eSBarry Smith   Note:
541b4319ba4SBarry Smith   The entries in the array that do not correspond to interface nodes remain unaltered.
54204c3f3b8SBarry Smith 
543562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`,
54404c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplySchur()`, `PCISApplyInvSchur()`,
54504c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `InsertMode`
54604c3f3b8SBarry Smith @*/
PCISScatterArrayNToVecB(PC pc,PetscScalar * array_N,Vec v_B,InsertMode imode,ScatterMode smode)54704c3f3b8SBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PC pc, PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode)
548d71ae5a4SJacob Faibussowitsch {
5495d0c19d7SBarry Smith   PetscInt        i;
5505d0c19d7SBarry Smith   const PetscInt *idex;
551b4319ba4SBarry Smith   PetscScalar    *array_B;
552f4f49eeaSPierre Jolivet   PC_IS          *pcis = (PC_IS *)pc->data;
553b4319ba4SBarry Smith 
554b4319ba4SBarry Smith   PetscFunctionBegin;
5559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_B, &array_B));
5569566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pcis->is_B_local, &idex));
557b4319ba4SBarry Smith 
558b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
559b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5602fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] = array_N[idex[i]];
561b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5622fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] += array_N[idex[i]];
563b4319ba4SBarry Smith     }
564b4319ba4SBarry Smith   } else { /* SCATTER_REVERSE */
565b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5662fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] = array_B[i];
567b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5682fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] += array_B[i];
569b4319ba4SBarry Smith     }
570b4319ba4SBarry Smith   }
5719566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pcis->is_B_local, &idex));
5729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_B, &array_B));
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
574b4319ba4SBarry Smith }
575b4319ba4SBarry Smith 
57604c3f3b8SBarry Smith /*@
577b4319ba4SBarry Smith   PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
57804c3f3b8SBarry Smith 
57904c3f3b8SBarry Smith   Input Parameters:
58004c3f3b8SBarry Smith + pc     - preconditioner context
58104c3f3b8SBarry Smith . b      - vector of local interface nodes (including ghosts)
58204c3f3b8SBarry Smith . x      - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur complement to `b`
58304c3f3b8SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); used as work space
58404c3f3b8SBarry Smith - vec2_N - vector of local nodes (interior and interface, including ghosts); used as work space
58504c3f3b8SBarry Smith 
58604c3f3b8SBarry Smith   Level: advanced
58704c3f3b8SBarry Smith 
58804c3f3b8SBarry Smith   Note:
58904c3f3b8SBarry Smith   Solves the problem
59004c3f3b8SBarry Smith .vb
591b4319ba4SBarry Smith   [ A_II  A_IB ] [ . ]   [ 0 ]
592b4319ba4SBarry Smith   [            ] [   ] = [   ]
593b4319ba4SBarry Smith   [ A_BI  A_BB ] [ x ]   [ b ]
59404c3f3b8SBarry Smith .ve
595b4319ba4SBarry Smith 
596562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
59704c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
59804c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
59904c3f3b8SBarry Smith @*/
PCISApplyInvSchur(PC pc,Vec b,Vec x,Vec vec1_N,Vec vec2_N)600d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
601d71ae5a4SJacob Faibussowitsch {
602f4f49eeaSPierre Jolivet   PC_IS *pcis = (PC_IS *)pc->data;
603b4319ba4SBarry Smith 
604b4319ba4SBarry Smith   PetscFunctionBegin;
605b4319ba4SBarry Smith   /*
606b4319ba4SBarry Smith     Neumann solvers.
607b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
608b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
609b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
610b4319ba4SBarry Smith     is stored in x.
611b4319ba4SBarry Smith   */
612b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
6139566063dSJacob Faibussowitsch   PetscCall(VecSet(vec1_N, 0.0));
6149566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
6159566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
616b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
617b4319ba4SBarry Smith   {
618ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
6199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_is_check_consistency", &flg, NULL));
620b4319ba4SBarry Smith     if (flg) {
621b4319ba4SBarry Smith       PetscScalar average;
6223050cee2SBarry Smith       PetscViewer viewer;
6239566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
6243050cee2SBarry Smith 
6259566063dSJacob Faibussowitsch       PetscCall(VecSum(vec1_N, &average));
626b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
6279566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
628b4319ba4SBarry Smith       if (pcis->pure_neumann) {
62963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is floating. Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
630b4319ba4SBarry Smith       } else {
63163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is fixed.    Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
632b4319ba4SBarry Smith       }
6339566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
6349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
635b4319ba4SBarry Smith     }
636b4319ba4SBarry Smith   }
637b4319ba4SBarry Smith   /* Solving the system for vec2_N */
6389566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_N, vec1_N, vec2_N));
6399566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_N, pc, vec2_N));
640b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
6419566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6429566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
644b4319ba4SBarry Smith }
645