xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 04c3f3b8f563eff258c1de90d4e02d41e6027585)
123ce1328SBarry Smith 
25e5bbd0aSStefano Zampini #include <petsc/private/pcisimpl.h> /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
5d71ae5a4SJacob Faibussowitsch {
6b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
7b965d45eSStefano Zampini 
8b965d45eSStefano Zampini   PetscFunctionBegin;
9b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11b965d45eSStefano Zampini }
12b965d45eSStefano Zampini 
13b965d45eSStefano Zampini /*@
14f1580f4eSBarry Smith   PCISSetUseStiffnessScaling - Tells `PCIS` to construct partition of unity using
15f1580f4eSBarry Smith   the local matrices' diagonal entries
16b965d45eSStefano Zampini 
17f1580f4eSBarry Smith   Logically Collective
18b965d45eSStefano Zampini 
19b965d45eSStefano Zampini   Input Parameters:
20b965d45eSStefano Zampini + pc  - the preconditioning context
2120f4b53cSBarry Smith - use - whether or not it should use matrix diagonal to build partition of unity.
22b965d45eSStefano Zampini 
23b965d45eSStefano Zampini   Level: intermediate
24b965d45eSStefano Zampini 
25*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
26*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
27*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
28b965d45eSStefano Zampini @*/
29d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
30d71ae5a4SJacob Faibussowitsch {
31b965d45eSStefano Zampini   PetscFunctionBegin;
32b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
33064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(pc, use, 2);
34cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetUseStiffnessScaling_C", (PC, PetscBool), (pc, use));
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36b965d45eSStefano Zampini }
37b965d45eSStefano Zampini 
38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
39d71ae5a4SJacob Faibussowitsch {
40ba1573a8SStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
41ba1573a8SStefano Zampini 
42ba1573a8SStefano Zampini   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)scaling_factors));
449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
45ba1573a8SStefano Zampini   pcis->D = scaling_factors;
46a007db60SStefano Zampini   if (pc->setupcalled) {
47a007db60SStefano Zampini     PetscInt sn;
48a007db60SStefano Zampini 
499566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
50a007db60SStefano Zampini     if (sn == pcis->n) {
519566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
529566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
539566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
549566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
559566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
5663a3b9bcSJacob 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);
57a007db60SStefano Zampini   }
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59ba1573a8SStefano Zampini }
60ba1573a8SStefano Zampini 
61ba1573a8SStefano Zampini /*@
62f1580f4eSBarry Smith   PCISSetSubdomainDiagonalScaling - Set diagonal scaling for `PCIS`.
63ba1573a8SStefano Zampini 
64f1580f4eSBarry Smith   Logically Collective
65ba1573a8SStefano Zampini 
66ba1573a8SStefano Zampini   Input Parameters:
67ba1573a8SStefano Zampini + pc              - the preconditioning context
68ba1573a8SStefano Zampini - scaling_factors - scaling factors for the subdomain
69ba1573a8SStefano Zampini 
70ba1573a8SStefano Zampini   Level: intermediate
71ba1573a8SStefano Zampini 
72f1580f4eSBarry Smith   Note:
73f1580f4eSBarry Smith   Intended for use with jumping coefficients cases.
74ba1573a8SStefano Zampini 
75*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
76*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISSetUseStiffnessScaling()`,
77*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
78ba1573a8SStefano Zampini @*/
79d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
80d71ae5a4SJacob Faibussowitsch {
81ba1573a8SStefano Zampini   PetscFunctionBegin;
82ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
83a007db60SStefano Zampini   PetscValidHeaderSpecific(scaling_factors, VEC_CLASSID, 2);
84cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainDiagonalScaling_C", (PC, Vec), (pc, scaling_factors));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86ba1573a8SStefano Zampini }
87ba1573a8SStefano Zampini 
88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
89d71ae5a4SJacob Faibussowitsch {
90831a100dSStefano Zampini   PC_IS *pcis = (PC_IS *)pc->data;
91831a100dSStefano Zampini 
92831a100dSStefano Zampini   PetscFunctionBegin;
93831a100dSStefano Zampini   pcis->scaling_factor = scal;
9448a46eb9SPierre Jolivet   if (pcis->D) PetscCall(VecSet(pcis->D, pcis->scaling_factor));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96831a100dSStefano Zampini }
97831a100dSStefano Zampini 
98831a100dSStefano Zampini /*@
99f1580f4eSBarry Smith   PCISSetSubdomainScalingFactor - Set scaling factor for `PCIS`.
100831a100dSStefano Zampini 
10120f4b53cSBarry Smith   Not Collective
102831a100dSStefano Zampini 
103831a100dSStefano Zampini   Input Parameters:
104831a100dSStefano Zampini + pc   - the preconditioning context
105831a100dSStefano Zampini - scal - scaling factor for the subdomain
106831a100dSStefano Zampini 
107831a100dSStefano Zampini   Level: intermediate
108831a100dSStefano Zampini 
109f1580f4eSBarry Smith   Note:
110f1580f4eSBarry Smith   Intended for use with the jumping coefficients cases.
111831a100dSStefano Zampini 
112*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISScatterArrayNToVecB()`,
113*04c3f3b8SBarry Smith           `PCISSetSubdomainDiagonalScaling()`, `PCISSetUseStiffnessScaling()`,
114*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `PCISApplyInvSchur()`, `PCISApplySchur()`
115831a100dSStefano Zampini @*/
116d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
117d71ae5a4SJacob Faibussowitsch {
118831a100dSStefano Zampini   PetscFunctionBegin;
119831a100dSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
120cac4c232SBarry Smith   PetscTryMethod(pc, "PCISSetSubdomainScalingFactor_C", (PC, PetscScalar), (pc, scal));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122831a100dSStefano Zampini }
123831a100dSStefano Zampini 
124*04c3f3b8SBarry Smith /*@
125*04c3f3b8SBarry Smith   PCISSetUp - sets up the `PC_IS` portion of `PCNN` and `PCBDDC` preconditioner context as part of their setup process
126*04c3f3b8SBarry Smith 
127*04c3f3b8SBarry Smith   Input Parameters:
128*04c3f3b8SBarry Smith + pc              - the `PC` object, must be of type `PCNN` or `PCBDDC`
129*04c3f3b8SBarry Smith . computematrices - Extract the blocks `A_II`, `A_BI`, `A_IB` and `A_BB` from the matrix
130*04c3f3b8SBarry Smith - computesolvers  - Create the `KSP` for the local Dirichlet and Neumann problems
131*04c3f3b8SBarry Smith 
132*04c3f3b8SBarry Smith   Level: advanced
133*04c3f3b8SBarry Smith 
134*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
135*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
136*04c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
137*04c3f3b8SBarry Smith @*/
138d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
139d71ae5a4SJacob Faibussowitsch {
140b4319ba4SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
141bf327c11SStefano Zampini   Mat_IS   *matis;
1425e8657edSStefano Zampini   MatReuse  reuse;
14385c21eb1SStefano Zampini   PetscBool flg, issbaij;
144b4319ba4SBarry Smith 
145b4319ba4SBarry Smith   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATIS, &flg));
14728b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires preconditioning matrix of type MATIS");
14885c21eb1SStefano Zampini   matis = (Mat_IS *)pc->pmat->data;
1492f37b69bSStefano Zampini   if (pc->useAmat) {
1509566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat, MATIS, &flg));
15128b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Requires linear system matrix of type MATIS");
1522f37b69bSStefano Zampini   }
153b4319ba4SBarry Smith 
1545e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1555e8657edSStefano Zampini   if (!pc->setupcalled) {
1565e8657edSStefano Zampini     PetscInt  n_I;
1575e8657edSStefano Zampini     PetscInt *idx_I_local, *idx_B_local, *idx_I_global, *idx_B_global;
1586f516dd7SStefano Zampini     PetscBT   bt;
1595e8657edSStefano Zampini     PetscInt  i, j;
160b4319ba4SBarry Smith 
161892d8026SStefano Zampini     /* get info on mapping */
1629566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
1639566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
164e432b41dSStefano Zampini     pcis->mapping = matis->rmapping;
1659566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping, &pcis->n));
1669566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
167892d8026SStefano Zampini 
168b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1699566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(pcis->n, &bt));
1704843afd6SStefano Zampini     for (i = 0; i < pcis->n_neigh; i++)
17148a46eb9SPierre Jolivet       for (j = 0; j < pcis->n_shared[i]; j++) PetscCall(PetscBTSet(bt, pcis->shared[i][j]));
172892d8026SStefano Zampini 
173da81f932SPierre Jolivet     /* Creating local and global index sets for interior and interface nodes. */
1749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_I_local));
1759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &idx_B_local));
176b4319ba4SBarry Smith     for (i = 0, pcis->n_B = 0, n_I = 0; i < pcis->n; i++) {
1776f516dd7SStefano Zampini       if (!PetscBTLookup(bt, i)) {
1782fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1792fa5cd67SKarl Rupp         n_I++;
1802fa5cd67SKarl Rupp       } else {
1812fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1822fa5cd67SKarl Rupp         pcis->n_B++;
1832fa5cd67SKarl Rupp       }
184b4319ba4SBarry Smith     }
1856f516dd7SStefano Zampini 
186b4319ba4SBarry Smith     /* Getting the global numbering */
187b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
188b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1899566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, pcis->n_B, idx_B_local, idx_B_global));
1909566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping, n_I, idx_I_local, idx_I_global));
1912fa5cd67SKarl Rupp 
1925e8657edSStefano Zampini     /* Creating the index sets */
1939566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pcis->n_B, idx_B_local, PETSC_COPY_VALUES, &pcis->is_B_local));
1949566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), pcis->n_B, idx_B_global, PETSC_COPY_VALUES, &pcis->is_B_global));
1959566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n_I, idx_I_local, PETSC_COPY_VALUES, &pcis->is_I_local));
1969566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), n_I, idx_I_global, PETSC_COPY_VALUES, &pcis->is_I_global));
1972fa5cd67SKarl Rupp 
1985e8657edSStefano Zampini     /* Freeing memory */
1999566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_B_local));
2009566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_I_local));
2019566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&bt));
202b4319ba4SBarry Smith 
2035e8657edSStefano Zampini     /* Creating work vectors and arrays */
2049566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(matis->x, &pcis->vec1_N));
2059566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_N, &pcis->vec2_N));
2069566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_D));
2079566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_D, pcis->n - pcis->n_B, PETSC_DECIDE));
2089566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_D, ((PetscObject)pcis->vec1_N)->type_name));
2099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec2_D));
2109566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec3_D));
2119566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D, &pcis->vec4_D));
2129566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &pcis->vec1_B));
2139566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_B, pcis->n_B, PETSC_DECIDE));
2149566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_B, ((PetscObject)pcis->vec1_N)->type_name));
2159566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec2_B));
2169566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B, &pcis->vec3_B));
2179566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &pcis->vec1_global, NULL));
2189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n, &pcis->work_N));
2195e8657edSStefano Zampini     /* scaling vector */
220bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2219566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2229566063dSJacob Faibussowitsch       PetscCall(VecSet(pcis->D, pcis->scaling_factor));
223bf83c2afSStefano Zampini     }
224b4319ba4SBarry Smith 
225b4319ba4SBarry Smith     /* Creating the scatter contexts */
2269566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_I_local, pcis->vec1_D, (IS)0, &pcis->N_to_D));
2279566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_I_global, pcis->vec1_D, (IS)0, &pcis->global_to_D));
2289566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N, pcis->is_B_local, pcis->vec1_B, (IS)0, &pcis->N_to_B));
2299566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pcis->is_B_global, pcis->vec1_B, (IS)0, &pcis->global_to_B));
230b4319ba4SBarry Smith 
2315e8657edSStefano Zampini     /* map from boundary to local */
2329566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local, &pcis->BtoNmap));
2335e8657edSStefano Zampini   }
2345e8657edSStefano Zampini 
235a007db60SStefano Zampini   {
236a007db60SStefano Zampini     PetscInt sn;
237a007db60SStefano Zampini 
2389566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D, &sn));
239a007db60SStefano Zampini     if (sn == pcis->n) {
2409566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2419566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
2429566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
2439566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B, &pcis->D));
2449566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B, pcis->D));
24563a3b9bcSJacob 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);
246a007db60SStefano Zampini   }
247a007db60SStefano Zampini 
2485e8657edSStefano Zampini   /*
2495e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2505e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2515e8657edSStefano Zampini 
2525e8657edSStefano Zampini         [ A_II | A_IB ]
2535e8657edSStefano Zampini     A = [------+------]
2545e8657edSStefano Zampini         [ A_BI | A_BB ]
2555e8657edSStefano Zampini   */
256d9869140SStefano Zampini   if (computematrices) {
2572f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2582f37b69bSStefano Zampini     PetscInt  bs, ibs;
2592f37b69bSStefano Zampini 
2605e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2613975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2625e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2635e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2645e8657edSStefano Zampini       } else {
2655e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2665e8657edSStefano Zampini       }
2675e8657edSStefano Zampini     }
2685e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2699566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2709566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->pA_II));
2719566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_IB));
2729566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BI));
2739566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BB));
2745e8657edSStefano Zampini     }
2755e8657edSStefano Zampini 
2769566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping, &ibs));
2779566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(matis->A, &bs));
2789566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->pA_II));
2792f37b69bSStefano Zampini     if (amat) {
2802f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS *)pc->mat->data;
2819566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(amatis->A, pcis->is_I_local, pcis->is_I_local, reuse, &pcis->A_II));
2822f37b69bSStefano Zampini     } else {
2839566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
2849566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2852f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2862f37b69bSStefano Zampini     }
2879566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_II, bs == ibs ? bs : 1));
2889566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->pA_II, bs == ibs ? bs : 1));
2899566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_B_local, reuse, &pcis->A_BB));
2909566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A, MATSEQSBAIJ, &issbaij));
2915e8657edSStefano Zampini     if (!issbaij) {
2929566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2939566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
2945e8657edSStefano Zampini     } else {
2955e8657edSStefano Zampini       Mat newmat;
296d9869140SStefano Zampini 
2979566063dSJacob Faibussowitsch       PetscCall(MatConvert(matis->A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &newmat));
2989566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_I_local, pcis->is_B_local, reuse, &pcis->A_IB));
2999566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat, pcis->is_B_local, pcis->is_I_local, reuse, &pcis->A_BI));
3009566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&newmat));
3015e8657edSStefano Zampini     }
3029566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_BB, bs == ibs ? bs : 1));
303d9869140SStefano Zampini   }
3045e8657edSStefano Zampini 
305bf83c2afSStefano Zampini   /* Creating scaling vector D */
3069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_is_use_stiffness_scaling", &pcis->use_stiffness_scaling, NULL));
307bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
30815f235b8SStefano Zampini     PetscScalar *a;
30915f235b8SStefano Zampini     PetscInt     i, n;
31015f235b8SStefano Zampini 
311d9869140SStefano Zampini     if (pcis->A_BB) {
3129566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pcis->A_BB, pcis->D));
313d9869140SStefano Zampini     } else {
3149566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matis->A, pcis->vec1_N));
3159566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
3169566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->D, INSERT_VALUES, SCATTER_FORWARD));
317d9869140SStefano Zampini     }
3189566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
3199566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D, &n));
3209566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D, &a));
3219371c9d4SSatish Balay     for (i = 0; i < n; i++)
3229371c9d4SSatish Balay       if (PetscAbsScalar(a[i]) < PETSC_SMALL) a[i] = 1.0;
3239566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D, &a));
324ba1573a8SStefano Zampini   }
3259566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global, 0.0));
3269566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3279566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->D, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3289566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3299566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
3309566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D, pcis->D, pcis->vec1_B));
331b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
332b4319ba4SBarry Smith 
3335e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
3345e8657edSStefano Zampini   if (computesolvers) {
335b4319ba4SBarry Smith     PC pc_ctx;
3365e8657edSStefano Zampini 
3375e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
338b4319ba4SBarry Smith     /* Dirichlet */
3399566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_D));
3403821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_D, pc->kspnestlevel));
3419566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D, pc->erroriffailure));
3429566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D, (PetscObject)pc, 1));
3439566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_D, pcis->A_II, pcis->A_II));
3449566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D, "is_localD_"));
3459566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_D, &pc_ctx));
3469566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3479566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_D, KSPPREONLY));
3489566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_D));
349b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3509566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_D));
351b4319ba4SBarry Smith     /* Neumann */
3529566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &pcis->ksp_N));
3533821be0aSBarry Smith     PetscCall(KSPSetNestLevel(pcis->ksp_N, pc->kspnestlevel));
3549566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N, pc->erroriffailure));
3559566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N, (PetscObject)pc, 1));
3569566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_N, matis->A, matis->A));
3579566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N, "is_localN_"));
3589566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_N, &pc_ctx));
3599566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx, PCLU));
3609566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_N, KSPPREONLY));
3619566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_N));
362b4319ba4SBarry Smith     {
3639371c9d4SSatish 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;
3649371c9d4SSatish Balay       PetscReal fixed_factor, floating_factor;
365b4319ba4SBarry Smith 
3669566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &fixed_factor, &damp_fixed));
3672fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
3689566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_damp_fixed", &damp_fixed, NULL));
369b4319ba4SBarry Smith 
3709566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_remove_nullspace_fixed", &remove_nullspace_fixed, NULL));
371b4319ba4SBarry Smith 
372d0609cedSBarry Smith       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &floating_factor, &set_damping_factor_floating));
3732fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
3749566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_set_damping_factor_floating", &set_damping_factor_floating, NULL));
3752fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
376b4319ba4SBarry Smith 
3779566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_damp_floating", &not_damp_floating, NULL));
378b4319ba4SBarry Smith 
3799566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options, ((PetscObject)pc_ctx)->prefix, "-pc_is_not_remove_nullspace_floating", &not_remove_nullspace_floating, NULL));
380b4319ba4SBarry Smith 
381b4319ba4SBarry Smith       if (pcis->pure_neumann) { /* floating subdomain */
382b4319ba4SBarry Smith         if (!(not_damp_floating)) {
3839566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3849566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
385b4319ba4SBarry Smith         }
386b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
387b4319ba4SBarry Smith           MatNullSpace nullsp;
3889566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
3899566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
3909566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
391b4319ba4SBarry Smith         }
392b4319ba4SBarry Smith       } else { /* fixed subdomain */
393b4319ba4SBarry Smith         if (damp_fixed) {
3949566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx, MAT_SHIFT_NONZERO));
3959566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx, floating_factor));
396b4319ba4SBarry Smith         }
397b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
398b4319ba4SBarry Smith           MatNullSpace nullsp;
3999566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullsp));
4009566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A, nullsp));
4019566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
402b4319ba4SBarry Smith         }
403b4319ba4SBarry Smith       }
404b4319ba4SBarry Smith     }
405b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
4069566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_N));
407b4319ba4SBarry Smith   }
4083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
409b4319ba4SBarry Smith }
410b4319ba4SBarry Smith 
411*04c3f3b8SBarry Smith /*@
412*04c3f3b8SBarry Smith   PCISReset - Removes all the `PC_IS` parts of the `PC` implementation data structure
413*04c3f3b8SBarry Smith 
414*04c3f3b8SBarry Smith   Input Parameter:
415*04c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
416*04c3f3b8SBarry Smith 
417*04c3f3b8SBarry Smith   Level: advanced
418*04c3f3b8SBarry Smith 
419*04c3f3b8SBarry Smith .seealso: `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`, `PCISSetSubdomainScalingFactor()`,
420*04c3f3b8SBarry Smith           `PCISInitialize()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
421*04c3f3b8SBarry Smith @*/
422*04c3f3b8SBarry Smith PetscErrorCode PCISReset(PC pc)
423d71ae5a4SJacob Faibussowitsch {
424*04c3f3b8SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
425*04c3f3b8SBarry Smith   PetscBool correcttype;
426b4319ba4SBarry Smith 
427b4319ba4SBarry Smith   PetscFunctionBegin;
4283ba16761SJacob Faibussowitsch   if (!pc) PetscFunctionReturn(PETSC_SUCCESS);
429*04c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
430*04c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_local));
4329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_local));
4339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_global));
4349566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_global));
4359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_II));
4369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->pA_II));
4379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_IB));
4389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BI));
4399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BB));
4409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
4419566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_N));
4429566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_D));
4439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_N));
4449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_N));
4459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_D));
4469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_D));
4479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_D));
4489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec4_D));
4499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_B));
4509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_B));
4519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_B));
4529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_global));
4539566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_D));
4549566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_B));
4559566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_D));
4569566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_B));
4579566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcis->work_N));
45848a46eb9SPierre Jolivet   if (pcis->n_neigh > -1) PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping, &(pcis->n_neigh), &(pcis->neigh), &(pcis->n_shared), &(pcis->shared)));
4599566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
4609566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", NULL));
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", NULL));
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", NULL));
4643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
465b4319ba4SBarry Smith }
466b4319ba4SBarry Smith 
467*04c3f3b8SBarry Smith /*@
468*04c3f3b8SBarry Smith   PCISInitialize - initializes the `PC_IS` portion of `PCNN` and `PCBDDC` preconditioner context
469*04c3f3b8SBarry Smith 
470*04c3f3b8SBarry Smith   Input Parameter:
471*04c3f3b8SBarry Smith . pc - the `PC` object, must be of type `PCNN` or `PCBDDC`
472*04c3f3b8SBarry Smith 
473*04c3f3b8SBarry Smith   Level: advanced
474*04c3f3b8SBarry Smith 
475*04c3f3b8SBarry Smith   Note:
476*04c3f3b8SBarry Smith   There is no preconditioner the `PCIS` prefixed routines provide functionality needed by `PCNN` or `PCBDDC`
477*04c3f3b8SBarry Smith 
478*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
479*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
480*04c3f3b8SBarry Smith           `PCISReset()`, `PCISApplySchur()`, `PCISApplyInvSchur()`
481*04c3f3b8SBarry Smith @*/
482*04c3f3b8SBarry Smith PetscErrorCode PCISInitialize(PC pc)
483d71ae5a4SJacob Faibussowitsch {
484b4319ba4SBarry Smith   PC_IS    *pcis = (PC_IS *)(pc->data);
485*04c3f3b8SBarry Smith   PetscBool correcttype;
486b4319ba4SBarry Smith 
487b4319ba4SBarry Smith   PetscFunctionBegin;
488*04c3f3b8SBarry Smith   PetscCheck(pcis, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC_IS context must be created by caller");
489*04c3f3b8SBarry Smith   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc, &correcttype, PCBDDC, PCNN, ""));
490*04c3f3b8SBarry Smith   PetscCheck(correcttype, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PC must be of type PCNN or PCBDDC");
4917dbfca69SStefano Zampini   pcis->n_neigh          = -1;
492831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4933975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetUseStiffnessScaling_C", PCISSetUseStiffnessScaling_IS));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainScalingFactor_C", PCISSetSubdomainScalingFactor_IS));
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCISSetSubdomainDiagonalScaling_C", PCISSetSubdomainDiagonalScaling_IS));
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498b4319ba4SBarry Smith }
499b4319ba4SBarry Smith 
500*04c3f3b8SBarry Smith /*@
501*04c3f3b8SBarry Smith   PCISApplySchur - applies the Schur complement arising from the `MATIS` inside the `PCNN` preconditioner
502b4319ba4SBarry Smith 
503*04c3f3b8SBarry Smith   Input Parameters:
504*04c3f3b8SBarry Smith + pc     - preconditioner context
505b4319ba4SBarry Smith . v      - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
506*04c3f3b8SBarry Smith . vec1_B - location to store the result of Schur complement applied to chunk
507*04c3f3b8SBarry Smith . vec2_B - workspace or `NULL`, `v` is used as workspace in that case
508*04c3f3b8SBarry Smith . vec1_D - work space
509*04c3f3b8SBarry Smith - vec2_D - work space
510b4319ba4SBarry Smith 
511*04c3f3b8SBarry Smith   Level: advanced
512b4319ba4SBarry Smith 
513*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
514*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplyInvSchur()`,
515*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
516*04c3f3b8SBarry Smith @*/
517d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
518d71ae5a4SJacob Faibussowitsch {
519b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
520b4319ba4SBarry Smith 
521b4319ba4SBarry Smith   PetscFunctionBegin;
5222fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
523b4319ba4SBarry Smith 
5249566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BB, v, vec1_B));
5259566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_IB, v, vec1_D));
5269566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_D, vec1_D, vec2_D));
5279566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_D, pc, vec2_D));
5289566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BI, vec2_D, vec2_B));
5299566063dSJacob Faibussowitsch   PetscCall(VecAXPY(vec1_B, -1.0, vec2_B));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531b4319ba4SBarry Smith }
532b4319ba4SBarry Smith 
533*04c3f3b8SBarry Smith /*@
534b4319ba4SBarry Smith   PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
535*04c3f3b8SBarry Smith   including ghosts) into an interface vector, when in `SCATTER_FORWARD` mode, or vice-versa, when in `SCATTER_REVERSE`
536b4319ba4SBarry Smith   mode.
537b4319ba4SBarry Smith 
538feefa0e1SJacob Faibussowitsch   Input Parameters:
539feefa0e1SJacob Faibussowitsch + pc      - preconditioner context
540*04c3f3b8SBarry Smith . array_N - [when in `SCATTER_FORWARD` mode] Array to be scattered into the vector otherwise output array
541*04c3f3b8SBarry Smith . imode   - insert mode, `ADD_VALUES` or `INSERT_VALUES`
542*04c3f3b8SBarry Smith . smode   - scatter mode, `SCATTER_FORWARD` or `SCATTER_REVERSE` mode]
543*04c3f3b8SBarry Smith - v_B     - [when in `SCATTER_REVERSE` mode] Vector to be scattered into the array, otherwise output vector
544b4319ba4SBarry Smith 
545*04c3f3b8SBarry Smith   Level: advanced
546b4319ba4SBarry Smith 
547f1580f4eSBarry Smith   Note:
548b4319ba4SBarry Smith   The entries in the array that do not correspond to interface nodes remain unaltered.
549*04c3f3b8SBarry Smith 
550*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`,
551*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`, `PCISApplySchur()`, `PCISApplyInvSchur()`,
552*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`, `InsertMode`
553*04c3f3b8SBarry Smith @*/
554*04c3f3b8SBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PC pc, PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode)
555d71ae5a4SJacob Faibussowitsch {
5565d0c19d7SBarry Smith   PetscInt        i;
5575d0c19d7SBarry Smith   const PetscInt *idex;
558b4319ba4SBarry Smith   PetscScalar    *array_B;
559b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS *)(pc->data);
560b4319ba4SBarry Smith 
561b4319ba4SBarry Smith   PetscFunctionBegin;
5629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_B, &array_B));
5639566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pcis->is_B_local, &idex));
564b4319ba4SBarry Smith 
565b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
566b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5672fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] = array_N[idex[i]];
568b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5692fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_B[i] += array_N[idex[i]];
570b4319ba4SBarry Smith     }
571b4319ba4SBarry Smith   } else { /* SCATTER_REVERSE */
572b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5732fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] = array_B[i];
574b4319ba4SBarry Smith     } else { /* ADD_VALUES */
5752fa5cd67SKarl Rupp       for (i = 0; i < pcis->n_B; i++) array_N[idex[i]] += array_B[i];
576b4319ba4SBarry Smith     }
577b4319ba4SBarry Smith   }
5789566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pcis->is_B_local, &idex));
5799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_B, &array_B));
5803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
581b4319ba4SBarry Smith }
582b4319ba4SBarry Smith 
583*04c3f3b8SBarry Smith /*@
584b4319ba4SBarry Smith   PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
585*04c3f3b8SBarry Smith 
586*04c3f3b8SBarry Smith   Input Parameters:
587*04c3f3b8SBarry Smith + pc     - preconditioner context
588*04c3f3b8SBarry Smith . b      - vector of local interface nodes (including ghosts)
589*04c3f3b8SBarry Smith . x      - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur complement to `b`
590*04c3f3b8SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); used as work space
591*04c3f3b8SBarry Smith - vec2_N - vector of local nodes (interior and interface, including ghosts); used as work space
592*04c3f3b8SBarry Smith 
593*04c3f3b8SBarry Smith   Level: advanced
594*04c3f3b8SBarry Smith 
595*04c3f3b8SBarry Smith   Note:
596*04c3f3b8SBarry Smith   Solves the problem
597*04c3f3b8SBarry Smith .vb
598b4319ba4SBarry Smith   [ A_II  A_IB ] [ . ]   [ 0 ]
599b4319ba4SBarry Smith   [            ] [   ] = [   ]
600b4319ba4SBarry Smith   [ A_BI  A_BB ] [ x ]   [ b ]
601*04c3f3b8SBarry Smith .ve
602b4319ba4SBarry Smith 
603*04c3f3b8SBarry Smith .seealso: `PCBDDC`, `PCNN`, `PCISSetUseStiffnessScaling()`, `PCISSetSubdomainDiagonalScaling()`, `PCISScatterArrayNToVecB()`,
604*04c3f3b8SBarry Smith           `PCISSetSubdomainScalingFactor()`,
605*04c3f3b8SBarry Smith           `PCISReset()`, `PCISInitialize()`
606*04c3f3b8SBarry Smith @*/
607d71ae5a4SJacob Faibussowitsch PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
608d71ae5a4SJacob Faibussowitsch {
609b4319ba4SBarry Smith   PC_IS *pcis = (PC_IS *)(pc->data);
610b4319ba4SBarry Smith 
611b4319ba4SBarry Smith   PetscFunctionBegin;
612b4319ba4SBarry Smith   /*
613b4319ba4SBarry Smith     Neumann solvers.
614b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
615b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
616b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
617b4319ba4SBarry Smith     is stored in x.
618b4319ba4SBarry Smith   */
619b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
6209566063dSJacob Faibussowitsch   PetscCall(VecSet(vec1_N, 0.0));
6219566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
6229566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, b, vec1_N, INSERT_VALUES, SCATTER_REVERSE));
623b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
624b4319ba4SBarry Smith   {
625ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
6269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_is_check_consistency", &flg, NULL));
627b4319ba4SBarry Smith     if (flg) {
628b4319ba4SBarry Smith       PetscScalar average;
6293050cee2SBarry Smith       PetscViewer viewer;
6309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
6313050cee2SBarry Smith 
6329566063dSJacob Faibussowitsch       PetscCall(VecSum(vec1_N, &average));
633b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
6349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
635b4319ba4SBarry Smith       if (pcis->pure_neumann) {
63663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is floating. Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
637b4319ba4SBarry Smith       } else {
63863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain %04d is fixed.    Average = % 1.14e\n", PetscGlobalRank, (double)PetscAbsScalar(average)));
639b4319ba4SBarry Smith       }
6409566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
6419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
642b4319ba4SBarry Smith     }
643b4319ba4SBarry Smith   }
644b4319ba4SBarry Smith   /* Solving the system for vec2_N */
6459566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_N, vec1_N, vec2_N));
6469566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_N, pc, vec2_N));
647b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
6489566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6499566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->N_to_B, vec2_N, x, INSERT_VALUES, SCATTER_FORWARD));
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
651b4319ba4SBarry Smith }
652