xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision 5e5bbd0a99ea6ef6b187abb9b4a6405ccdf6db1e)
123ce1328SBarry Smith 
2*5e5bbd0aSStefano Zampini #include <petsc/private/pcisimpl.h> /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
4b965d45eSStefano Zampini static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
5b965d45eSStefano Zampini {
6b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
7b965d45eSStefano Zampini 
8b965d45eSStefano Zampini   PetscFunctionBegin;
9b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
10b965d45eSStefano Zampini   PetscFunctionReturn(0);
11b965d45eSStefano Zampini }
12b965d45eSStefano Zampini 
13b965d45eSStefano Zampini /*@
14b965d45eSStefano Zampini  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
15b965d45eSStefano Zampini                               local matrices' diagonal.
16b965d45eSStefano Zampini 
17b965d45eSStefano Zampini    Not collective
18b965d45eSStefano Zampini 
19b965d45eSStefano Zampini    Input Parameters:
20b965d45eSStefano Zampini +  pc - the preconditioning context
21b965d45eSStefano Zampini -  use - whether or not pcis use matrix diagonal to build partition of unity.
22b965d45eSStefano Zampini 
23b965d45eSStefano Zampini    Level: intermediate
24b965d45eSStefano Zampini 
25b965d45eSStefano Zampini    Notes:
26b965d45eSStefano Zampini 
27db781477SPatrick Sanan .seealso: `PCBDDC`
28b965d45eSStefano Zampini @*/
29b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
30b965d45eSStefano Zampini {
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));
35b965d45eSStefano Zampini   PetscFunctionReturn(0);
36b965d45eSStefano Zampini }
37b965d45eSStefano Zampini 
38ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
39ba1573a8SStefano Zampini {
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   }
58ba1573a8SStefano Zampini   PetscFunctionReturn(0);
59ba1573a8SStefano Zampini }
60ba1573a8SStefano Zampini 
61ba1573a8SStefano Zampini /*@
62ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
63ba1573a8SStefano Zampini 
64ba1573a8SStefano Zampini    Not 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 
72ba1573a8SStefano Zampini    Notes:
73ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
74ba1573a8SStefano Zampini 
75db781477SPatrick Sanan .seealso: `PCBDDC`
76ba1573a8SStefano Zampini @*/
77ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
78ba1573a8SStefano Zampini {
79ba1573a8SStefano Zampini   PetscFunctionBegin;
80ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
81a007db60SStefano Zampini   PetscValidHeaderSpecific(scaling_factors,VEC_CLASSID,2);
82cac4c232SBarry Smith   PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));
83ba1573a8SStefano Zampini   PetscFunctionReturn(0);
84ba1573a8SStefano Zampini }
85ba1573a8SStefano Zampini 
86831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
87831a100dSStefano Zampini {
88831a100dSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
89831a100dSStefano Zampini 
90831a100dSStefano Zampini   PetscFunctionBegin;
91831a100dSStefano Zampini   pcis->scaling_factor = scal;
92a007db60SStefano Zampini   if (pcis->D) {
93a007db60SStefano Zampini 
949566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->D,pcis->scaling_factor));
95a007db60SStefano Zampini   }
96831a100dSStefano Zampini   PetscFunctionReturn(0);
97831a100dSStefano Zampini }
98831a100dSStefano Zampini 
99831a100dSStefano Zampini /*@
100831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
101831a100dSStefano Zampini 
102831a100dSStefano Zampini    Not collective
103831a100dSStefano Zampini 
104831a100dSStefano Zampini    Input Parameters:
105831a100dSStefano Zampini +  pc - the preconditioning context
106831a100dSStefano Zampini -  scal - scaling factor for the subdomain
107831a100dSStefano Zampini 
108831a100dSStefano Zampini    Level: intermediate
109831a100dSStefano Zampini 
110831a100dSStefano Zampini    Notes:
111831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
112831a100dSStefano Zampini 
113db781477SPatrick Sanan .seealso: `PCBDDC`
114831a100dSStefano Zampini @*/
115831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
116831a100dSStefano Zampini {
117831a100dSStefano Zampini   PetscFunctionBegin;
118831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119cac4c232SBarry Smith   PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));
120831a100dSStefano Zampini   PetscFunctionReturn(0);
121831a100dSStefano Zampini }
122831a100dSStefano Zampini 
123b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
124b4319ba4SBarry Smith /*
125b4319ba4SBarry Smith    PCISSetUp -
126b4319ba4SBarry Smith */
127d9869140SStefano Zampini PetscErrorCode  PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers)
128b4319ba4SBarry Smith {
129b4319ba4SBarry Smith   PC_IS          *pcis  = (PC_IS*)(pc->data);
130bf327c11SStefano Zampini   Mat_IS         *matis;
1315e8657edSStefano Zampini   MatReuse       reuse;
13285c21eb1SStefano Zampini   PetscBool      flg,issbaij;
133b4319ba4SBarry Smith 
134b4319ba4SBarry Smith   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg));
13628b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires preconditioning matrix of type MATIS");
13785c21eb1SStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1382f37b69bSStefano Zampini   if (pc->useAmat) {
1399566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg));
14028b400f6SJacob Faibussowitsch     PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Requires linear system matrix of type MATIS");
1412f37b69bSStefano Zampini   }
142b4319ba4SBarry Smith 
1435e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1445e8657edSStefano Zampini   if (!pc->setupcalled) {
1455e8657edSStefano Zampini     PetscInt    n_I;
1465e8657edSStefano Zampini     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
1476f516dd7SStefano Zampini     PetscBT     bt;
1485e8657edSStefano Zampini     PetscInt    i,j;
149b4319ba4SBarry Smith 
150892d8026SStefano Zampini     /* get info on mapping */
1519566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)matis->rmapping));
1529566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
153e432b41dSStefano Zampini     pcis->mapping = matis->rmapping;
1549566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n));
1559566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared)));
156892d8026SStefano Zampini 
157b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
1589566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(pcis->n,&bt));
1594843afd6SStefano Zampini     for (i=0;i<pcis->n_neigh;i++)
1606f516dd7SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
1619566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(bt,pcis->shared[i][j]));
1626f516dd7SStefano Zampini       }
163892d8026SStefano Zampini 
1645e8657edSStefano Zampini     /* Creating local and global index sets for interior and inteface nodes. */
1659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n,&idx_I_local));
1669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n,&idx_B_local));
167b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
1686f516dd7SStefano Zampini       if (!PetscBTLookup(bt,i)) {
1692fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1702fa5cd67SKarl Rupp         n_I++;
1712fa5cd67SKarl Rupp       } else {
1722fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1732fa5cd67SKarl Rupp         pcis->n_B++;
1742fa5cd67SKarl Rupp       }
175b4319ba4SBarry Smith     }
1766f516dd7SStefano Zampini 
177b4319ba4SBarry Smith     /* Getting the global numbering */
178b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
179b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1809566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global));
1819566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global));
1822fa5cd67SKarl Rupp 
1835e8657edSStefano Zampini     /* Creating the index sets */
1849566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local));
1859566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global));
1869566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local));
1879566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global));
1882fa5cd67SKarl Rupp 
1895e8657edSStefano Zampini     /* Freeing memory */
1909566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_B_local));
1919566063dSJacob Faibussowitsch     PetscCall(PetscFree(idx_I_local));
1929566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&bt));
193b4319ba4SBarry Smith 
1945e8657edSStefano Zampini     /* Creating work vectors and arrays */
1959566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(matis->x,&pcis->vec1_N));
1969566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_N,&pcis->vec2_N));
1979566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF,&pcis->vec1_D));
1989566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_D,pcis->n-pcis->n_B,PETSC_DECIDE));
1999566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_D,((PetscObject)pcis->vec1_N)->type_name));
2009566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec2_D));
2019566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec3_D));
2029566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_D,&pcis->vec4_D));
2039566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF,&pcis->vec1_B));
2049566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pcis->vec1_B,pcis->n_B,PETSC_DECIDE));
2059566063dSJacob Faibussowitsch     PetscCall(VecSetType(pcis->vec1_B,((PetscObject)pcis->vec1_N)->type_name));
2069566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B,&pcis->vec2_B));
2079566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_B,&pcis->vec3_B));
2089566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat,&pcis->vec1_global,NULL));
2099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n,&pcis->work_N));
2105e8657edSStefano Zampini     /* scaling vector */
211bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2129566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B,&pcis->D));
2139566063dSJacob Faibussowitsch       PetscCall(VecSet(pcis->D,pcis->scaling_factor));
214bf83c2afSStefano Zampini     }
215b4319ba4SBarry Smith 
216b4319ba4SBarry Smith     /* Creating the scatter contexts */
2179566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D));
2189566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D));
2199566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B));
2209566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B));
221b4319ba4SBarry Smith 
2225e8657edSStefano Zampini     /* map from boundary to local */
2239566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap));
2245e8657edSStefano Zampini   }
2255e8657edSStefano Zampini 
226a007db60SStefano Zampini   {
227a007db60SStefano Zampini     PetscInt sn;
228a007db60SStefano Zampini 
2299566063dSJacob Faibussowitsch     PetscCall(VecGetSize(pcis->D,&sn));
230a007db60SStefano Zampini     if (sn == pcis->n) {
2319566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
2329566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
2339566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pcis->D));
2349566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(pcis->vec1_B,&pcis->D));
2359566063dSJacob Faibussowitsch       PetscCall(VecCopy(pcis->vec1_B,pcis->D));
23663a3b9bcSJacob 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);
237a007db60SStefano Zampini   }
238a007db60SStefano Zampini 
2395e8657edSStefano Zampini   /*
2405e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2415e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2425e8657edSStefano Zampini 
2435e8657edSStefano Zampini         [ A_II | A_IB ]
2445e8657edSStefano Zampini     A = [------+------]
2455e8657edSStefano Zampini         [ A_BI | A_BB ]
2465e8657edSStefano Zampini   */
247d9869140SStefano Zampini   if (computematrices) {
2482f37b69bSStefano Zampini     PetscBool amat = (PetscBool)(pc->mat != pc->pmat && pc->useAmat);
2492f37b69bSStefano Zampini     PetscInt  bs,ibs;
2502f37b69bSStefano Zampini 
2515e8657edSStefano Zampini     reuse = MAT_INITIAL_MATRIX;
2523975b054SStefano Zampini     if (pcis->reusesubmatrices && pc->setupcalled) {
2535e8657edSStefano Zampini       if (pc->flag == SAME_NONZERO_PATTERN) {
2545e8657edSStefano Zampini         reuse = MAT_REUSE_MATRIX;
2555e8657edSStefano Zampini       } else {
2565e8657edSStefano Zampini         reuse = MAT_INITIAL_MATRIX;
2575e8657edSStefano Zampini       }
2585e8657edSStefano Zampini     }
2595e8657edSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2609566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2619566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->pA_II));
2629566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_IB));
2639566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BI));
2649566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_BB));
2655e8657edSStefano Zampini     }
2665e8657edSStefano Zampini 
2679566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(pcis->mapping,&ibs));
2689566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(matis->A,&bs));
2699566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->pA_II));
2702f37b69bSStefano Zampini     if (amat) {
2712f37b69bSStefano Zampini       Mat_IS *amatis = (Mat_IS*)pc->mat->data;
2729566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(amatis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II));
2732f37b69bSStefano Zampini     } else {
2749566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)pcis->pA_II));
2759566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pcis->A_II));
2762f37b69bSStefano Zampini       pcis->A_II = pcis->pA_II;
2772f37b69bSStefano Zampini     }
2789566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_II,bs == ibs ? bs : 1));
2799566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->pA_II,bs == ibs ? bs : 1));
2809566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB));
2819566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij));
2825e8657edSStefano Zampini     if (!issbaij) {
2839566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB));
2849566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI));
2855e8657edSStefano Zampini     } else {
2865e8657edSStefano Zampini       Mat newmat;
287d9869140SStefano Zampini 
2889566063dSJacob Faibussowitsch       PetscCall(MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat));
2899566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB));
2909566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI));
2919566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&newmat));
2925e8657edSStefano Zampini     }
2939566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(pcis->A_BB,bs == ibs ? bs : 1));
294d9869140SStefano Zampini   }
2955e8657edSStefano Zampini 
296bf83c2afSStefano Zampini   /* Creating scaling vector D */
2979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL));
298bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
29915f235b8SStefano Zampini     PetscScalar *a;
30015f235b8SStefano Zampini     PetscInt    i,n;
30115f235b8SStefano Zampini 
302d9869140SStefano Zampini     if (pcis->A_BB) {
3039566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pcis->A_BB,pcis->D));
304d9869140SStefano Zampini     } else {
3059566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matis->A,pcis->vec1_N));
3069566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
3079566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
308d9869140SStefano Zampini     }
3099566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
3109566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D,&n));
3119566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D,&a));
31215f235b8SStefano Zampini     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
3139566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D,&a));
314ba1573a8SStefano Zampini   }
3159566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global,0.0));
3169566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3179566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3189566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
3199566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
3209566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B));
321b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
322b4319ba4SBarry Smith 
3235e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
3245e8657edSStefano Zampini   if (computesolvers) {
325b4319ba4SBarry Smith     PC pc_ctx;
3265e8657edSStefano Zampini 
3275e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
328b4319ba4SBarry Smith     /* Dirichlet */
3299566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D));
3309566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure));
3319566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1));
3329566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II));
3339566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_"));
3349566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_D,&pc_ctx));
3359566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx,PCLU));
3369566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_D,KSPPREONLY));
3379566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_D));
338b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3399566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_D));
340b4319ba4SBarry Smith     /* Neumann */
3419566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N));
3429566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure));
3439566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1));
3449566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(pcis->ksp_N,matis->A,matis->A));
3459566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_"));
3469566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcis->ksp_N,&pc_ctx));
3479566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc_ctx,PCLU));
3489566063dSJacob Faibussowitsch     PetscCall(KSPSetType(pcis->ksp_N,KSPPREONLY));
3499566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(pcis->ksp_N));
350b4319ba4SBarry Smith     {
351ace3abfcSBarry Smith       PetscBool damp_fixed                    = PETSC_FALSE,
35290d69ab7SBarry Smith                 remove_nullspace_fixed        = PETSC_FALSE,
35390d69ab7SBarry Smith                 set_damping_factor_floating   = PETSC_FALSE,
35490d69ab7SBarry Smith                 not_damp_floating             = PETSC_FALSE,
35590d69ab7SBarry Smith                 not_remove_nullspace_floating = PETSC_FALSE;
356b4319ba4SBarry Smith       PetscReal fixed_factor,
357b4319ba4SBarry Smith                 floating_factor;
358b4319ba4SBarry Smith 
3599566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed));
3602fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
3619566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL));
362b4319ba4SBarry Smith 
3639566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL));
364b4319ba4SBarry Smith 
365d0609cedSBarry Smith       PetscCall(PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&floating_factor,&set_damping_factor_floating));
3662fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
3679566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL));
3682fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
369b4319ba4SBarry Smith 
3709566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL));
371b4319ba4SBarry Smith 
3729566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL));
373b4319ba4SBarry Smith 
374b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
375b4319ba4SBarry Smith         if (!(not_damp_floating)) {
3769566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO));
3779566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx,floating_factor));
378b4319ba4SBarry Smith         }
379b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
380b4319ba4SBarry Smith           MatNullSpace nullsp;
3819566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp));
3829566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A,nullsp));
3839566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
384b4319ba4SBarry Smith         }
385b4319ba4SBarry Smith       } else {  /* fixed subdomain */
386b4319ba4SBarry Smith         if (damp_fixed) {
3879566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO));
3889566063dSJacob Faibussowitsch           PetscCall(PCFactorSetShiftAmount(pc_ctx,floating_factor));
389b4319ba4SBarry Smith         }
390b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
391b4319ba4SBarry Smith           MatNullSpace nullsp;
3929566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp));
3939566063dSJacob Faibussowitsch           PetscCall(MatSetNullSpace(matis->A,nullsp));
3949566063dSJacob Faibussowitsch           PetscCall(MatNullSpaceDestroy(&nullsp));
395b4319ba4SBarry Smith         }
396b4319ba4SBarry Smith       }
397b4319ba4SBarry Smith     }
398b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
3999566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(pcis->ksp_N));
400b4319ba4SBarry Smith   }
401b4319ba4SBarry Smith   PetscFunctionReturn(0);
402b4319ba4SBarry Smith }
403b4319ba4SBarry Smith 
404b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
405b4319ba4SBarry Smith /*
406b4319ba4SBarry Smith    PCISDestroy -
407b4319ba4SBarry Smith */
4087087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
409b4319ba4SBarry Smith {
410747c8291SStefano Zampini   PC_IS *pcis;
411b4319ba4SBarry Smith 
412b4319ba4SBarry Smith   PetscFunctionBegin;
413747c8291SStefano Zampini   if (!pc) PetscFunctionReturn(0);
414747c8291SStefano Zampini   pcis = (PC_IS*)(pc->data);
4159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_local));
4169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_local));
4179566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_B_global));
4189566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pcis->is_I_global));
4199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_II));
4209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->pA_II));
4219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_IB));
4229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BI));
4239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pcis->A_BB));
4249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->D));
4259566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_N));
4269566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&pcis->ksp_D));
4279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_N));
4289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_N));
4299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_D));
4309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_D));
4319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_D));
4329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec4_D));
4339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_B));
4349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec2_B));
4359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec3_B));
4369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcis->vec1_global));
4379566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_D));
4389566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_B));
4399566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->N_to_D));
4409566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pcis->global_to_B));
4419566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcis->work_N));
4427dbfca69SStefano Zampini   if (pcis->n_neigh > -1) {
4439566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared)));
4447dbfca69SStefano Zampini   }
4459566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->mapping));
4469566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&pcis->BtoNmap));
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL));
4489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL));
4499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL));
450b4319ba4SBarry Smith   PetscFunctionReturn(0);
451b4319ba4SBarry Smith }
452b4319ba4SBarry Smith 
453b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
454b4319ba4SBarry Smith /*
455b4319ba4SBarry Smith    PCISCreate -
456b4319ba4SBarry Smith */
4577087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
458b4319ba4SBarry Smith {
459b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
460b4319ba4SBarry Smith 
461b4319ba4SBarry Smith   PetscFunctionBegin;
462747c8291SStefano Zampini   if (!pcis) {
463747c8291SStefano Zampini     PetscCall(PetscNewLog(pc,&pcis));
464747c8291SStefano Zampini     pc->data = pcis;
465747c8291SStefano Zampini   }
4667dbfca69SStefano Zampini   pcis->n_neigh          = -1;
467831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4683975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
469831a100dSStefano Zampini   /* composing functions */
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS));
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS));
4729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS));
473b4319ba4SBarry Smith   PetscFunctionReturn(0);
474b4319ba4SBarry Smith }
475b4319ba4SBarry Smith 
476b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
477b4319ba4SBarry Smith /*
478b4319ba4SBarry Smith    PCISApplySchur -
479b4319ba4SBarry Smith 
480b4319ba4SBarry Smith    Input parameters:
481b4319ba4SBarry Smith .  pc - preconditioner context
482b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
483b4319ba4SBarry Smith 
484b4319ba4SBarry Smith    Output parameters:
485b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
486b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
487b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
488b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith */
4917087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
492b4319ba4SBarry Smith {
493b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
494b4319ba4SBarry Smith 
495b4319ba4SBarry Smith   PetscFunctionBegin;
4962fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
497b4319ba4SBarry Smith 
4989566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BB,v,vec1_B));
4999566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_IB,v,vec1_D));
5009566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_D,vec1_D,vec2_D));
5019566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_D,pc,vec2_D));
5029566063dSJacob Faibussowitsch   PetscCall(MatMult(pcis->A_BI,vec2_D,vec2_B));
5039566063dSJacob Faibussowitsch   PetscCall(VecAXPY(vec1_B,-1.0,vec2_B));
504b4319ba4SBarry Smith   PetscFunctionReturn(0);
505b4319ba4SBarry Smith }
506b4319ba4SBarry Smith 
507b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
508b4319ba4SBarry Smith /*
509b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
510b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
511b4319ba4SBarry Smith    mode.
512b4319ba4SBarry Smith 
513b4319ba4SBarry Smith    Input parameters:
514b4319ba4SBarry Smith .  pc - preconditioner context
515b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
516b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
517b4319ba4SBarry Smith 
518b4319ba4SBarry Smith    Output parameter:
519b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
520b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
521b4319ba4SBarry Smith 
522b4319ba4SBarry Smith    Notes:
523b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
524b4319ba4SBarry Smith */
5257087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
526b4319ba4SBarry Smith {
5275d0c19d7SBarry Smith   PetscInt       i;
5285d0c19d7SBarry Smith   const PetscInt *idex;
529b4319ba4SBarry Smith   PetscScalar    *array_B;
530b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
531b4319ba4SBarry Smith 
532b4319ba4SBarry Smith   PetscFunctionBegin;
5339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_B,&array_B));
5349566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pcis->is_B_local,&idex));
535b4319ba4SBarry Smith 
536b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
537b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5382fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
539b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5402fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
541b4319ba4SBarry Smith     }
542b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
543b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5442fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
545b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5462fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
547b4319ba4SBarry Smith     }
548b4319ba4SBarry Smith   }
5499566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pcis->is_B_local,&idex));
5509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_B,&array_B));
551b4319ba4SBarry Smith   PetscFunctionReturn(0);
552b4319ba4SBarry Smith }
553b4319ba4SBarry Smith 
554b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
555b4319ba4SBarry Smith /*
556b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
557b4319ba4SBarry Smith    More precisely, solves the problem:
558b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
559b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
560b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
561b4319ba4SBarry Smith 
562b4319ba4SBarry Smith    Input parameters:
563b4319ba4SBarry Smith .  pc - preconditioner context
564b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
565b4319ba4SBarry Smith 
566b4319ba4SBarry Smith    Output parameters:
567b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
568b4319ba4SBarry Smith        complement to b
569b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
570b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
571b4319ba4SBarry Smith 
572b4319ba4SBarry Smith */
5737087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
574b4319ba4SBarry Smith {
575b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
576b4319ba4SBarry Smith 
577b4319ba4SBarry Smith   PetscFunctionBegin;
578b4319ba4SBarry Smith   /*
579b4319ba4SBarry Smith     Neumann solvers.
580b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
581b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
582b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
583b4319ba4SBarry Smith     is stored in x.
584b4319ba4SBarry Smith   */
585b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
5869566063dSJacob Faibussowitsch   PetscCall(VecSet(vec1_N,0.0));
5879566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE));
5889566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE));
589b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
590b4319ba4SBarry Smith   {
591ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
5929566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL));
593b4319ba4SBarry Smith     if (flg) {
594b4319ba4SBarry Smith       PetscScalar average;
5953050cee2SBarry Smith       PetscViewer viewer;
5969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer));
5973050cee2SBarry Smith 
5989566063dSJacob Faibussowitsch       PetscCall(VecSum(vec1_N,&average));
599b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
6009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
601b4319ba4SBarry Smith       if (pcis->pure_neumann) {
60263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,(double)PetscAbsScalar(average)));
603b4319ba4SBarry Smith       } else {
60463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,(double)PetscAbsScalar(average)));
605b4319ba4SBarry Smith       }
6069566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
6079566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
608b4319ba4SBarry Smith     }
609b4319ba4SBarry Smith   }
610b4319ba4SBarry Smith   /* Solving the system for vec2_N */
6119566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcis->ksp_N,vec1_N,vec2_N));
6129566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(pcis->ksp_N,pc,vec2_N));
613b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
6149566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD));
6159566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD));
616b4319ba4SBarry Smith   PetscFunctionReturn(0);
617b4319ba4SBarry Smith }
618