xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
15e5bbd0aSStefano Zampini #include <petsc/private/pcisimpl.h> /*I "petscpc.h" I*/
2831a100dSStefano Zampini 
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 
24*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
2504c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
2604c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
27b965d45eSStefano Zampini @*/
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 
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 
74*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
7504c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISSetUseStiffnessScaling()`,
7604c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
77ba1573a8SStefano Zampini @*/
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 
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 
111*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
11204c3f3b8SBarry Smith           `PCISSetSubdomainDiagonalScaling()`, `PCISSetUseStiffnessScaling()`,
11304c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
114831a100dSStefano Zampini @*/
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 
133*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
13404c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
13504c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
13604c3f3b8SBarry Smith @*/
137d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
138d71ae5a4SJacob Faibussowitsch {
139b4319ba4SBarry Smith   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));
14628b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires preconditioning matrix 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;
1565e8657edSStefano Zampini     PetscInt *idx_I_local, *idx_B_local, *idx_I_global, *idx_B_global;
1576f516dd7SStefano Zampini     PetscBT   bt;
1585e8657edSStefano Zampini     PetscInt  i, j;
159b4319ba4SBarry Smith 
160892d8026SStefano Zampini     /* get info on mapping */
1619566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
1629566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
163e432b41dSStefano Zampini     pcis->mapping = matis->rmapping;
1649566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping, &pcis->n));
1659566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
166892d8026SStefano Zampini 
167b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1689566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(pcis->n, &bt));
1694843afd6SStefano Zampini     for (i = 0; i < pcis->n_neigh; i++)
17048a46eb9SPierre Jolivet       for (j = 0; j < pcis->n_shared[i]; j++) PetscCall(PetscBTSet(bt, pcis->shared[i][j]));
171892d8026SStefano Zampini 
172da81f932SPierre Jolivet     /* Creating local and global index sets for interior and interface nodes. */
1739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_I_local));
1749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_B_local));
175b4319ba4SBarry Smith     for (i = 0, pcis->n_B = 0, n_I = 0; i < pcis->n; i++) {
1766f516dd7SStefano Zampini       if (!PetscBTLookup(bt, i)) {
1772fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1782fa5cd67SKarl Rupp         n_I++;
1792fa5cd67SKarl Rupp       } else {
1802fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1812fa5cd67SKarl Rupp         pcis->n_B++;
1822fa5cd67SKarl Rupp       }
183b4319ba4SBarry Smith     }
1846f516dd7SStefano Zampini 
185b4319ba4SBarry Smith     /* Getting the global numbering */
186b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
187b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1889566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, pcis->n_B, idx_B_local, idx_B_global));
1899566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, n_I, idx_I_local, idx_I_global));
1902fa5cd67SKarl Rupp 
1915e8657edSStefano Zampini     /* Creating the index sets */
1929566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pcis->n_B, idx_B_local, PETSC_COPY_VALUES, &pcis->is_B_local));
1939566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), pcis->n_B, idx_B_global, PETSC_COPY_VALUES, &pcis->is_B_global));
1949566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n_I, idx_I_local, PETSC_COPY_VALUES, &pcis->is_I_local));
1959566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), n_I, idx_I_global, PETSC_COPY_VALUES, &pcis->is_I_global));
1962fa5cd67SKarl Rupp 
1975e8657edSStefano Zampini     /* Freeing memory */
1989566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_B_local));
1999566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_I_local));
2009566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&bt));
201b4319ba4SBarry Smith 
2025e8657edSStefano Zampini     /* Creating work vectors and arrays */
2039566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(matis->x, &pcis->vec1_N));
2049566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_N, &pcis->vec2_N));
2059566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_D));
2069566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_D, pcis->n - pcis->n_B, PETSC_DECIDE));
2079566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_D, ((PetscObject)pcis->vec1_N)->type_name));
2089566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec2_D));
2099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec3_D));
2109566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec4_D));
2119566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_B));
2129566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_B, pcis->n_B, PETSC_DECIDE));
2139566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_B, ((PetscObject)pcis->vec1_N)->type_name));
2149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec2_B));
2159566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec3_B));
2169566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &pcis->vec1_global, NULL));
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &pcis->work_N));
2185e8657edSStefano Zampini     /* scaling vector */
219bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2209566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2219566063dSJacob Faibussowitsch       PetscCall(VecSet(pcis->D, pcis->scaling_factor));
222bf83c2afSStefano Zampini     }
223b4319ba4SBarry Smith 
224b4319ba4SBarry Smith     /* Creating the scatter contexts */
2259566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_I_local, pcis->vec1_D, (IS)0, &pcis->N_to_D));
2269566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_I_global, pcis->vec1_D, (IS)0, &pcis->global_to_D));
2279566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_B_local, pcis->vec1_B, (IS)0, &pcis->N_to_B));
2289566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_B_global, pcis->vec1_B, (IS)0, &pcis->global_to_B));
229b4319ba4SBarry Smith 
2305e8657edSStefano Zampini     /* map from boundary to local */
2319566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local, &pcis->BtoNmap));
2325e8657edSStefano Zampini   }
2335e8657edSStefano Zampini 
234a007db60SStefano Zampini   {
235a007db60SStefano Zampini     PetscInt sn;
236a007db60SStefano Zampini 
2379566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
238a007db60SStefano Zampini     if (sn == pcis->n) {
2399566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2409566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2419566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
2429566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2439566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
24463a3b9bcSJacob 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);
245a007db60SStefano Zampini   }
246a007db60SStefano Zampini 
2475e8657edSStefano Zampini   /*
2485e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2495e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2505e8657edSStefano Zampini 
2515e8657edSStefano Zampini         [ A_II | A_IB ]
2525e8657edSStefano Zampini     A = [------+------]
2535e8657edSStefano Zampini         [ A_BI | A_BB ]
2545e8657edSStefano Zampini   */
255d9869140SStefano Zampini   if (computematrices) {
2562f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2572f37b69bSStefano Zampini     PetscInt  bs, ibs;
2582f37b69bSStefano Zampini 
2595e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2603975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2615e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2625e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2635e8657edSStefano Zampini       } else {
2645e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2655e8657edSStefano Zampini       }
2665e8657edSStefano Zampini     }
2675e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2689566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2699566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->pA_II));
2709566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_IB));
2719566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BI));
2729566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BB));
2735e8657edSStefano Zampini     }
2745e8657edSStefano Zampini 
2759566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping, &ibs));
2769566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(matis->A, &bs));
2779566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->pA_II));
2782f37b69bSStefano Zampini     if (amat) {
2792f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS *)pc->mat->data;
2809566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(amatis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->A_II));
2812f37b69bSStefano Zampini     } else {
2829566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
2839566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2842f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2852f37b69bSStefano Zampini     }
2869566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_II, bs == ibs ? bs : 1));
2879566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->pA_II, bs == ibs ? bs : 1));
2889566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_B_local, reuse, &pcis->A_BB));
2899566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &issbaij));
2905e8657edSStefano Zampini     if (!issbaij) {
2919566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2929566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2935e8657edSStefano Zampini     } else {
2945e8657edSStefano Zampini       Mat newmat;
295d9869140SStefano Zampini 
2969566063dSJacob Faibussowitsch       PetscCall(MatConvert(matis->A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &newmat));
2979566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2989566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2999566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&newmat));
3005e8657edSStefano Zampini     }
3019566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_BB, bs == ibs ? bs : 1));
302d9869140SStefano Zampini   }
3035e8657edSStefano Zampini 
304bf83c2afSStefano Zampini   /* Creating scaling vector D */
3059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_is_use_stiffness_scaling", &pcis->use_stiffness_scaling, NULL));
306bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
30715f235b8SStefano Zampini     PetscScalar *a;
30815f235b8SStefano Zampini     PetscInt     i, n;
30915f235b8SStefano Zampini 
310d9869140SStefano Zampini     if (pcis->A_BB) {
3119566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pcis->A_BB, pcis->D));
312d9869140SStefano Zampini     } else {
3139566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matis->A, pcis->vec1_N));
3149566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
3159566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
316d9869140SStefano Zampini     }
3179566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
3189566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D, &n));
3199566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D, &a));
3209371c9d4SSatish Balay     for (i = 0; i < n; i++)
3219371c9d4SSatish Balay       if (PetscAbsScalar(a[i]) < PETSC_SMALL) a[i] = 1.0;
3229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D, &a));
323ba1573a8SStefano Zampini   }
3249566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global, 0.0));
3259566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3269566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3279566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3289566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3299566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D, pcis->D, pcis->vec1_B));
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 */
3389566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_D));
3393821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_D, pc->kspnestlevel));
3409566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D, pc->erroriffailure));
3419566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D, (PetscObject)pc, 1));
3429566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_D, pcis->A_II, pcis->A_II));
3439566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D, "is_localD_"));
3449566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_D, &pc_ctx));
3459566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3469566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_D, KSPPREONLY));
3479566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_D));
348b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3499566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_D));
350b4319ba4SBarry Smith     /* Neumann */
3519566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_N));
3523821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_N, pc->kspnestlevel));
3539566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N, pc->erroriffailure));
3549566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N, (PetscObject)pc, 1));
3559566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_N, matis->A, matis->A));
3569566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N, "is_localN_"));
3579566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_N, &pc_ctx));
3589566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3599566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_N, KSPPREONLY));
3609566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_N));
361b4319ba4SBarry Smith     {
3629371c9d4SSatish 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;
3639371c9d4SSatish Balay       PetscReal fixed_factor, floating_factor;
364b4319ba4SBarry Smith 
3659566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &fixed_factor, &damp_fixed));
3662fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
3679566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &damp_fixed, NULL));
368b4319ba4SBarry Smith 
3699566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_remove_nullspace_fixed", &remove_nullspace_fixed, NULL));
370b4319ba4SBarry Smith 
371d0609cedSBarry Smith       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &floating_factor, &set_damping_factor_floating));
3722fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
3739566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &set_damping_factor_floating, NULL));
3742fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
375b4319ba4SBarry Smith 
3769566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_damp_floating", &not_damp_floating, NULL));
377b4319ba4SBarry Smith 
3789566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_remove_nullspace_floating", &not_remove_nullspace_floating, NULL));
379b4319ba4SBarry Smith 
380b4319ba4SBarry Smith       if (pcis->pure_neumann) { /* floating subdomain */
381b4319ba4SBarry Smith         if (!(not_damp_floating)) {
3829566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3839566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
384b4319ba4SBarry Smith         }
385b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
386b4319ba4SBarry Smith           MatNullSpace nullsp;
3879566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3889566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3899566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
390b4319ba4SBarry Smith         }
391b4319ba4SBarry Smith       } else { /* fixed subdomain */
392b4319ba4SBarry Smith         if (damp_fixed) {
3939566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3949566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
395b4319ba4SBarry Smith         }
396b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
397b4319ba4SBarry Smith           MatNullSpace nullsp;
3989566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3999566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
4009566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
401b4319ba4SBarry Smith         }
402b4319ba4SBarry Smith       }
403b4319ba4SBarry Smith     }
404b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
4059566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_N));
406b4319ba4SBarry Smith   }
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408b4319ba4SBarry Smith }
409b4319ba4SBarry Smith 
41004c3f3b8SBarry Smith /*@
41104c3f3b8SBarry Smith   PCISReset - Removes all the `PC_IS` parts of the `PC` implementation data structure
41204c3f3b8SBarry Smith 
41304c3f3b8SBarry Smith   Input Parameter:
41404c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
41504c3f3b8SBarry Smith 
41604c3f3b8SBarry Smith   Level: advanced
41704c3f3b8SBarry Smith 
418*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`, `PCISSetSubdomainScalingFactor()`,
41904c3f3b8SBarry Smith           `PCISInitialize()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
42004c3f3b8SBarry Smith @*/
42104c3f3b8SBarry Smith PetscErrorCode PCISReset(PC pc)
422d71ae5a4SJacob Faibussowitsch {
42304c3f3b8SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
42404c3f3b8SBarry Smith   PetscBool correcttype;
425b4319ba4SBarry Smith 
426b4319ba4SBarry Smith   PetscFunctionBegin;
4273ba16761SJacob Faibussowitsch   if (!pc) PetscFunctionReturn(PETSC_SUCCESS);
42804c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
42904c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_local));
4319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_local));
4329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_global));
4339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_global));
4349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_II));
4359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->pA_II));
4369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_IB));
4379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BI));
4389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BB));
4399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
4409566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_N));
4419566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_D));
4429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_N));
4439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_N));
4449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_D));
4459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_D));
4469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_D));
4479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec4_D));
4489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_B));
4499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_B));
4509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_B));
4519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_global));
4529566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_D));
4539566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_B));
4549566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_D));
4559566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_B));
4569566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcis->work_N));
45748a46eb9SPierre Jolivet   if (pcis->n_neigh > -1) PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
4589566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
4599566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", NULL));
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", NULL));
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", NULL));
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
464b4319ba4SBarry Smith }
465b4319ba4SBarry Smith 
46604c3f3b8SBarry Smith /*@
46704c3f3b8SBarry Smith   PCISInitialize - initializes the `PC_IS` portion of `PCNN` and `PCBDDC` preconditioner context
46804c3f3b8SBarry Smith 
46904c3f3b8SBarry Smith   Input Parameter:
47004c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
47104c3f3b8SBarry Smith 
47204c3f3b8SBarry Smith   Level: advanced
47304c3f3b8SBarry Smith 
47404c3f3b8SBarry Smith   Note:
47504c3f3b8SBarry Smith   There is no preconditioner the `PCIS` prefixed routines provide functionality needed by `PCNN` or `PCBDDC`
47604c3f3b8SBarry Smith 
477*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
47804c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
47904c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
48004c3f3b8SBarry Smith @*/
48104c3f3b8SBarry Smith PetscErrorCode PCISInitialize(PC pc)
482d71ae5a4SJacob Faibussowitsch {
483b4319ba4SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
48404c3f3b8SBarry Smith   PetscBool correcttype;
485b4319ba4SBarry Smith 
486b4319ba4SBarry Smith   PetscFunctionBegin;
48704c3f3b8SBarry Smith   PetscCheck(pcis, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC_IS context must be created by caller");
48804c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
48904c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4907dbfca69SStefano Zampini   pcis->n_neigh          = -1;
491831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4923975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", PCISSetUseStiffnessScaling_IS));
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", PCISSetSubdomainScalingFactor_IS));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", PCISSetSubdomainDiagonalScaling_IS));
4963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
497b4319ba4SBarry Smith }
498b4319ba4SBarry Smith 
49904c3f3b8SBarry Smith /*@
50004c3f3b8SBarry Smith   PCISApplySchur - applies the Schur complement arising from the `MATIS` inside the `PCNN` preconditioner
501b4319ba4SBarry Smith 
50204c3f3b8SBarry Smith   Input Parameters:
50304c3f3b8SBarry Smith + pc     - preconditioner context
504b4319ba4SBarry Smith . v      - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
50504c3f3b8SBarry Smith . vec1_B - location to store the result of Schur complement applied to chunk
50604c3f3b8SBarry Smith . vec2_B - workspace or `NULL`, `v` is used as workspace in that case
50704c3f3b8SBarry Smith . vec1_D - work space
50804c3f3b8SBarry Smith - vec2_D - work space
509b4319ba4SBarry Smith 
51004c3f3b8SBarry Smith   Level: advanced
511b4319ba4SBarry Smith 
512*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
51304c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplyInvSchur()`,
51404c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
51504c3f3b8SBarry Smith @*/
516d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
517d71ae5a4SJacob Faibussowitsch {
518b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
519b4319ba4SBarry Smith 
520b4319ba4SBarry Smith   PetscFunctionBegin;
5212fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
522b4319ba4SBarry Smith 
5239566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BB, v, vec1_B));
5249566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_IB, v, vec1_D));
5259566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_D, vec1_D, vec2_D));
5269566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_D, pc, vec2_D));
5279566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BI, vec2_D, vec2_B));
5289566063dSJacob Faibussowitsch   PetscCall(VecAXPY(vec1_B, -1.0, vec2_B));
5293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
530b4319ba4SBarry Smith }
531b4319ba4SBarry Smith 
53204c3f3b8SBarry Smith /*@
533b4319ba4SBarry Smith   PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
53404c3f3b8SBarry Smith   including ghosts) into an interface vector, when in `SCATTER_FORWARD` mode, or vice-versa, when in `SCATTER_REVERSE`
535b4319ba4SBarry Smith   mode.
536b4319ba4SBarry Smith 
537feefa0e1SJacob Faibussowitsch   Input Parameters:
538feefa0e1SJacob Faibussowitsch + pc      - preconditioner context
53904c3f3b8SBarry Smith . array_N - [when in `SCATTER_FORWARD` mode] Array to be scattered into the vector otherwise output array
54004c3f3b8SBarry Smith . imode   - insert mode, `ADD_VALUES` or `INSERT_VALUES`
54104c3f3b8SBarry Smith . smode   - scatter mode, `SCATTER_FORWARD` or `SCATTER_REVERSE` mode]
54204c3f3b8SBarry Smith - v_B     - [when in `SCATTER_REVERSE` mode] Vector to be scattered into the array, otherwise output vector
543b4319ba4SBarry Smith 
54404c3f3b8SBarry Smith   Level: advanced
545b4319ba4SBarry Smith 
546f1580f4eSBarry Smith   Note:
547b4319ba4SBarry Smith   The entries in the array that do not correspond to interface nodes remain unaltered.
54804c3f3b8SBarry Smith 
549*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`,
55004c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplySchur()`, `PCISApplyInvSchur()`,
55104c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `InsertMode`
55204c3f3b8SBarry Smith @*/
55304c3f3b8SBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PC pc, PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode)
554d71ae5a4SJacob Faibussowitsch {
5555d0c19d7SBarry Smith   PetscInt        i;
5565d0c19d7SBarry Smith   const PetscInt *idex;
557b4319ba4SBarry Smith   PetscScalar    *array_B;
558b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS *)(pc->data);
559b4319ba4SBarry Smith 
560b4319ba4SBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_B, &array_B));
5629566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pcis->is_B_local, &idex));
563b4319ba4SBarry Smith 
564b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
565b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5662fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] = array_N[idex[i]];
567b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5682fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] += array_N[idex[i]];
569b4319ba4SBarry Smith     }
570b4319ba4SBarry Smith   } else { /* SCATTER_REVERSE */
571b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5722fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] = array_B[i];
573b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5742fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] += array_B[i];
575b4319ba4SBarry Smith     }
576b4319ba4SBarry Smith   }
5779566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pcis->is_B_local, &idex));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_B, &array_B));
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
580b4319ba4SBarry Smith }
581b4319ba4SBarry Smith 
58204c3f3b8SBarry Smith /*@
583b4319ba4SBarry Smith   PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
58404c3f3b8SBarry Smith 
58504c3f3b8SBarry Smith   Input Parameters:
58604c3f3b8SBarry Smith + pc     - preconditioner context
58704c3f3b8SBarry Smith . b      - vector of local interface nodes (including ghosts)
58804c3f3b8SBarry Smith . x      - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur complement to `b`
58904c3f3b8SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); used as work space
59004c3f3b8SBarry Smith - vec2_N - vector of local nodes (interior and interface, including ghosts); used as work space
59104c3f3b8SBarry Smith 
59204c3f3b8SBarry Smith   Level: advanced
59304c3f3b8SBarry Smith 
59404c3f3b8SBarry Smith   Note:
59504c3f3b8SBarry Smith   Solves the problem
59604c3f3b8SBarry Smith .vb
597b4319ba4SBarry Smith   [ A_II  A_IB ] [ . ]   [ 0 ]
598b4319ba4SBarry Smith   [            ] [   ] = [   ]
599b4319ba4SBarry Smith   [ A_BI  A_BB ] [ x ]   [ b ]
60004c3f3b8SBarry Smith .ve
601b4319ba4SBarry Smith 
602*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
60304c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
60404c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
60504c3f3b8SBarry Smith @*/
606d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
607d71ae5a4SJacob Faibussowitsch {
608b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
609b4319ba4SBarry Smith 
610b4319ba4SBarry Smith   PetscFunctionBegin;
611b4319ba4SBarry Smith   /*
612b4319ba4SBarry Smith     Neumann solvers.
613b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
614b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
615b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
616b4319ba4SBarry Smith     is stored in x.
617b4319ba4SBarry Smith   */
618b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
6199566063dSJacob Faibussowitsch   PetscCall(VecSet(vec1_N, 0.0));
6209566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
6219566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
622b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
623b4319ba4SBarry Smith   {
624ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
6259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_is_check_consistency", &flg, NULL));
626b4319ba4SBarry Smith     if (flg) {
627b4319ba4SBarry Smith       PetscScalar average;
6283050cee2SBarry Smith       PetscViewer viewer;
6299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
6303050cee2SBarry Smith 
6319566063dSJacob Faibussowitsch       PetscCall(VecSum(vec1_N, &average));
632b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
6339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
634b4319ba4SBarry Smith       if (pcis->pure_neumann) {
63563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is floating. Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
636b4319ba4SBarry Smith       } else {
63763a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is fixed.    Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
638b4319ba4SBarry Smith       }
6399566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
6409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
641b4319ba4SBarry Smith     }
642b4319ba4SBarry Smith   }
643b4319ba4SBarry Smith   /* Solving the system for vec2_N */
6449566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_N, vec1_N, vec2_N));
6459566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_N, pc, vec2_N));
646b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
6479566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6489566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650b4319ba4SBarry Smith }
651