xref: /petsc/src/ksp/pc/impls/is/pcis.c (revision bf83c2af4b1da2dd6967cb3385e75e2eedcb93dc)
123ce1328SBarry Smith 
2aaa7dc30SBarry Smith #include <../src/ksp/pc/impls/is/pcis.h> /*I "petscpc.h" I*/
3831a100dSStefano Zampini 
4831a100dSStefano Zampini #undef __FUNCT__
5b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling_IS"
6b965d45eSStefano Zampini static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use)
7b965d45eSStefano Zampini {
8b965d45eSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
9b965d45eSStefano Zampini 
10b965d45eSStefano Zampini   PetscFunctionBegin;
11b965d45eSStefano Zampini   pcis->use_stiffness_scaling = use;
12b965d45eSStefano Zampini   PetscFunctionReturn(0);
13b965d45eSStefano Zampini }
14b965d45eSStefano Zampini 
15b965d45eSStefano Zampini #undef __FUNCT__
16b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling"
17b965d45eSStefano Zampini /*@
18b965d45eSStefano Zampini  PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using
19b965d45eSStefano Zampini                               local matrices' diagonal.
20b965d45eSStefano Zampini 
21b965d45eSStefano Zampini    Not collective
22b965d45eSStefano Zampini 
23b965d45eSStefano Zampini    Input Parameters:
24b965d45eSStefano Zampini +  pc - the preconditioning context
25b965d45eSStefano Zampini -  use - whether or not pcis use matrix diagonal to build partition of unity.
26b965d45eSStefano Zampini 
27b965d45eSStefano Zampini    Level: intermediate
28b965d45eSStefano Zampini 
29b965d45eSStefano Zampini    Notes:
30b965d45eSStefano Zampini 
31b965d45eSStefano Zampini .seealso: PCBDDC
32b965d45eSStefano Zampini @*/
33b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use)
34b965d45eSStefano Zampini {
35b965d45eSStefano Zampini   PetscErrorCode ierr;
36b965d45eSStefano Zampini 
37b965d45eSStefano Zampini   PetscFunctionBegin;
38b965d45eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
39b965d45eSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
40b965d45eSStefano Zampini   PetscFunctionReturn(0);
41b965d45eSStefano Zampini }
42b965d45eSStefano Zampini 
43b965d45eSStefano Zampini #undef __FUNCT__
44ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling_IS"
45ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors)
46ba1573a8SStefano Zampini {
47ba1573a8SStefano Zampini   PetscErrorCode ierr;
48ba1573a8SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
49ba1573a8SStefano Zampini 
50ba1573a8SStefano Zampini   PetscFunctionBegin;
51ba1573a8SStefano Zampini   ierr    = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr);
52*bf83c2afSStefano Zampini   ierr    = VecDestroy(&pcis->D);CHKERRQ(ierr);
53ba1573a8SStefano Zampini   pcis->D = scaling_factors;
54ba1573a8SStefano Zampini   PetscFunctionReturn(0);
55ba1573a8SStefano Zampini }
56ba1573a8SStefano Zampini 
57ba1573a8SStefano Zampini #undef __FUNCT__
58ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling"
59ba1573a8SStefano Zampini /*@
60ba1573a8SStefano Zampini  PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS.
61ba1573a8SStefano Zampini 
62ba1573a8SStefano Zampini    Not collective
63ba1573a8SStefano Zampini 
64ba1573a8SStefano Zampini    Input Parameters:
65ba1573a8SStefano Zampini +  pc - the preconditioning context
66ba1573a8SStefano Zampini -  scaling_factors - scaling factors for the subdomain
67ba1573a8SStefano Zampini 
68ba1573a8SStefano Zampini    Level: intermediate
69ba1573a8SStefano Zampini 
70ba1573a8SStefano Zampini    Notes:
71ba1573a8SStefano Zampini    Intended to use with jumping coefficients cases.
72ba1573a8SStefano Zampini 
73ba1573a8SStefano Zampini .seealso: PCBDDC
74ba1573a8SStefano Zampini @*/
75ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors)
76ba1573a8SStefano Zampini {
77ba1573a8SStefano Zampini   PetscErrorCode ierr;
78ba1573a8SStefano Zampini 
79ba1573a8SStefano Zampini   PetscFunctionBegin;
80ba1573a8SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
81ba1573a8SStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr);
82ba1573a8SStefano Zampini   PetscFunctionReturn(0);
83ba1573a8SStefano Zampini }
84ba1573a8SStefano Zampini 
85ba1573a8SStefano Zampini #undef __FUNCT__
86831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS"
87831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal)
88831a100dSStefano Zampini {
89831a100dSStefano Zampini   PC_IS *pcis = (PC_IS*)pc->data;
90831a100dSStefano Zampini 
91831a100dSStefano Zampini   PetscFunctionBegin;
92831a100dSStefano Zampini   pcis->scaling_factor = scal;
93831a100dSStefano Zampini   PetscFunctionReturn(0);
94831a100dSStefano Zampini }
95831a100dSStefano Zampini 
96831a100dSStefano Zampini #undef __FUNCT__
97831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor"
98831a100dSStefano Zampini /*@
99831a100dSStefano Zampini  PCISSetSubdomainScalingFactor - Set scaling factor for PCIS.
100831a100dSStefano Zampini 
101831a100dSStefano Zampini    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 
109831a100dSStefano Zampini    Notes:
110831a100dSStefano Zampini    Intended to use with jumping coefficients cases.
111831a100dSStefano Zampini 
112831a100dSStefano Zampini .seealso: PCBDDC
113831a100dSStefano Zampini @*/
114831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal)
115831a100dSStefano Zampini {
116831a100dSStefano Zampini   PetscErrorCode ierr;
117831a100dSStefano Zampini 
118831a100dSStefano Zampini   PetscFunctionBegin;
119831a100dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120831a100dSStefano Zampini   ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr);
121831a100dSStefano Zampini   PetscFunctionReturn(0);
122831a100dSStefano Zampini }
123831a100dSStefano Zampini 
124b4319ba4SBarry Smith 
125b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
126b4319ba4SBarry Smith /*
127b4319ba4SBarry Smith    PCISSetUp -
128b4319ba4SBarry Smith */
129b4319ba4SBarry Smith #undef __FUNCT__
130b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp"
1315e8657edSStefano Zampini PetscErrorCode  PCISSetUp(PC pc, PetscBool computesolvers)
132b4319ba4SBarry Smith {
133b4319ba4SBarry Smith   PC_IS          *pcis  = (PC_IS*)(pc->data);
134bf327c11SStefano Zampini   Mat_IS         *matis;
1355e8657edSStefano Zampini   MatReuse       reuse;
1366849ba73SBarry Smith   PetscErrorCode ierr;
13785c21eb1SStefano Zampini   PetscBool      flg,issbaij;
138831a100dSStefano Zampini   Vec            counter;
139b4319ba4SBarry Smith 
140b4319ba4SBarry Smith   PetscFunctionBegin;
14185c21eb1SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr);
142ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS");
14385c21eb1SStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
144b4319ba4SBarry Smith 
1455e8657edSStefano Zampini   /* first time creation, get info on substructuring */
1465e8657edSStefano Zampini   if (!pc->setupcalled) {
1475e8657edSStefano Zampini     PetscInt    n_I;
1485e8657edSStefano Zampini     PetscInt    *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global;
1495e8657edSStefano Zampini     PetscInt    *array;
1505e8657edSStefano Zampini     PetscInt    i,j;
151b4319ba4SBarry Smith 
152892d8026SStefano Zampini     /* get info on mapping */
1533bbff08aSStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr);
1548d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
1553bbff08aSStefano Zampini     pcis->mapping = pc->pmat->rmap->mapping;
1568d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr);
1578d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
158892d8026SStefano Zampini 
159b4319ba4SBarry Smith     /* Identifying interior and interface nodes, in local numbering */
160785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&array);CHKERRQ(ierr);
161892d8026SStefano Zampini     ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr);
1624843afd6SStefano Zampini     for (i=0;i<pcis->n_neigh;i++)
163892d8026SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++)
164892d8026SStefano Zampini           array[pcis->shared[i][j]] += 1;
165892d8026SStefano Zampini 
1665e8657edSStefano Zampini     /* Creating local and global index sets for interior and inteface nodes. */
167785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr);
168785e854fSJed Brown     ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr);
169b4319ba4SBarry Smith     for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) {
170892d8026SStefano Zampini       if (!array[i]) {
1712fa5cd67SKarl Rupp         idx_I_local[n_I] = i;
1722fa5cd67SKarl Rupp         n_I++;
1732fa5cd67SKarl Rupp       } else {
1742fa5cd67SKarl Rupp         idx_B_local[pcis->n_B] = i;
1752fa5cd67SKarl Rupp         pcis->n_B++;
1762fa5cd67SKarl Rupp       }
177b4319ba4SBarry Smith     }
178b4319ba4SBarry Smith     /* Getting the global numbering */
179b4319ba4SBarry Smith     idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */
180b4319ba4SBarry Smith     idx_I_global = idx_B_local + pcis->n_B;
1818d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr);
1828d12c799SStefano Zampini     ierr         = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr);
1832fa5cd67SKarl Rupp 
1845e8657edSStefano Zampini     /* Creating the index sets */
185bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr);
186bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr);
187bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr);
188bb943df9SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr);
1892fa5cd67SKarl Rupp 
1905e8657edSStefano Zampini     /* Freeing memory */
191b4319ba4SBarry Smith     ierr = PetscFree(idx_B_local);CHKERRQ(ierr);
192b4319ba4SBarry Smith     ierr = PetscFree(idx_I_local);CHKERRQ(ierr);
193892d8026SStefano Zampini     ierr = PetscFree(array);CHKERRQ(ierr);
194b4319ba4SBarry Smith 
1955e8657edSStefano Zampini     /* Creating work vectors and arrays */
196892d8026SStefano Zampini     ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr);
197b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr);
198b4319ba4SBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr);
199b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr);
200b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr);
201df187020SStefano Zampini     ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr);
202b4319ba4SBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr);
203b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr);
204b4319ba4SBarry Smith     ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr);
2052a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr);
206854ce69bSBarry Smith     ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr);
2075e8657edSStefano Zampini     /* scaling vector */
208*bf83c2afSStefano Zampini     if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */
2095e8657edSStefano Zampini       ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr);
210*bf83c2afSStefano Zampini       ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr);
211*bf83c2afSStefano Zampini     }
212b4319ba4SBarry Smith 
213b4319ba4SBarry Smith     /* Creating the scatter contexts */
21416909a7fSStefano Zampini     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr);
21523ce1328SBarry Smith     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr);
216b4319ba4SBarry Smith     ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr);
21723ce1328SBarry Smith     ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr);
218b4319ba4SBarry Smith 
2195e8657edSStefano Zampini     /* map from boundary to local */
2205e8657edSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr);
2215e8657edSStefano Zampini   }
2225e8657edSStefano Zampini 
2235e8657edSStefano Zampini   /*
2245e8657edSStefano Zampini     Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering
2255e8657edSStefano Zampini     is such that interior nodes come first than the interface ones, we have
2265e8657edSStefano Zampini 
2275e8657edSStefano Zampini         [ A_II | A_IB ]
2285e8657edSStefano Zampini     A = [------+------]
2295e8657edSStefano Zampini         [ A_BI | A_BB ]
2305e8657edSStefano Zampini   */
2315e8657edSStefano Zampini   reuse = MAT_INITIAL_MATRIX;
2323975b054SStefano Zampini   if (pcis->reusesubmatrices && pc->setupcalled) {
2335e8657edSStefano Zampini     if (pc->flag == SAME_NONZERO_PATTERN) {
2345e8657edSStefano Zampini       reuse = MAT_REUSE_MATRIX;
2355e8657edSStefano Zampini     } else {
2365e8657edSStefano Zampini       reuse = MAT_INITIAL_MATRIX;
2375e8657edSStefano Zampini     }
2385e8657edSStefano Zampini   }
2395e8657edSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
2405e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
2415e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
2425e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
2435e8657edSStefano Zampini     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
2445e8657edSStefano Zampini   }
2455e8657edSStefano Zampini 
2465e8657edSStefano Zampini   ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
2475e8657edSStefano Zampini   ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
2485e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
2495e8657edSStefano Zampini   if (!issbaij) {
2505e8657edSStefano Zampini     ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2515e8657edSStefano Zampini     ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
2525e8657edSStefano Zampini   } else {
2535e8657edSStefano Zampini     Mat newmat;
2545e8657edSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr);
2555e8657edSStefano Zampini     ierr = MatGetSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
2565e8657edSStefano Zampini     ierr = MatGetSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
2575e8657edSStefano Zampini     ierr = MatDestroy(&newmat);CHKERRQ(ierr);
2585e8657edSStefano Zampini   }
2595e8657edSStefano Zampini 
260*bf83c2afSStefano Zampini   /* Creating scaling vector D */
261c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr);
262*bf83c2afSStefano Zampini   if (pcis->use_stiffness_scaling) {
263*bf83c2afSStefano Zampini     ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr);
264ba1573a8SStefano Zampini   }
2652a7a6963SBarry Smith   ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */
266831a100dSStefano Zampini   ierr = VecSet(counter,0.0);CHKERRQ(ierr);
267*bf83c2afSStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
268*bf83c2afSStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
269ba1573a8SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
270ba1573a8SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
271ba1573a8SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
272892d8026SStefano Zampini   ierr = VecDestroy(&counter);CHKERRQ(ierr);
273b4319ba4SBarry Smith 
274b4319ba4SBarry Smith   /* See historical note 01, at the bottom of this file. */
275b4319ba4SBarry Smith 
2765e8657edSStefano Zampini   /* Creating the KSP contexts for the local Dirichlet and Neumann problems */
2775e8657edSStefano Zampini   if (computesolvers) {
278b4319ba4SBarry Smith     PC pc_ctx;
2795e8657edSStefano Zampini 
2805e8657edSStefano Zampini     pcis->pure_neumann = matis->pure_neumann;
281b4319ba4SBarry Smith     /* Dirichlet */
282b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr);
283422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr);
2841cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
28523ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr);
286b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr);
287b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr);
288b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
289b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr);
290b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr);
291b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
292b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr);
293b4319ba4SBarry Smith     /* Neumann */
294b4319ba4SBarry Smith     ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr);
295422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr);
2961cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr);
29723ee1639SBarry Smith     ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr);
298b4319ba4SBarry Smith     ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr);
299b4319ba4SBarry Smith     ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr);
300b4319ba4SBarry Smith     ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr);
301b4319ba4SBarry Smith     ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr);
302b4319ba4SBarry Smith     ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr);
303b4319ba4SBarry Smith     {
304ace3abfcSBarry Smith       PetscBool damp_fixed                    = PETSC_FALSE,
30590d69ab7SBarry Smith                 remove_nullspace_fixed        = PETSC_FALSE,
30690d69ab7SBarry Smith                 set_damping_factor_floating   = PETSC_FALSE,
30790d69ab7SBarry Smith                 not_damp_floating             = PETSC_FALSE,
30890d69ab7SBarry Smith                 not_remove_nullspace_floating = PETSC_FALSE;
309b4319ba4SBarry Smith       PetscReal fixed_factor,
310b4319ba4SBarry Smith                 floating_factor;
311b4319ba4SBarry Smith 
312c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr);
3132fa5cd67SKarl Rupp       if (!damp_fixed) fixed_factor = 0.0;
314c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr);
315b4319ba4SBarry Smith 
316c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr);
317b4319ba4SBarry Smith 
318c5929fdfSBarry Smith       ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",
319b4319ba4SBarry Smith                               &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr);
3202fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 0.0;
321c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr);
3222fa5cd67SKarl Rupp       if (!set_damping_factor_floating) floating_factor = 1.e-12;
323b4319ba4SBarry Smith 
324c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",&not_damp_floating,NULL);CHKERRQ(ierr);
325b4319ba4SBarry Smith 
326c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",&not_remove_nullspace_floating,NULL);CHKERRQ(ierr);
327b4319ba4SBarry Smith 
328b4319ba4SBarry Smith       if (pcis->pure_neumann) {  /* floating subdomain */
329b4319ba4SBarry Smith         if (!(not_damp_floating)) {
330c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
331c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
332b4319ba4SBarry Smith         }
333b4319ba4SBarry Smith         if (!(not_remove_nullspace_floating)) {
334b4319ba4SBarry Smith           MatNullSpace nullsp;
3350298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
3365fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
337d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
338b4319ba4SBarry Smith         }
339b4319ba4SBarry Smith       } else {  /* fixed subdomain */
340b4319ba4SBarry Smith         if (damp_fixed) {
341c88783f4SHong Zhang           ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr);
342c88783f4SHong Zhang           ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr);
343b4319ba4SBarry Smith         }
344b4319ba4SBarry Smith         if (remove_nullspace_fixed) {
345b4319ba4SBarry Smith           MatNullSpace nullsp;
3460298fd71SBarry Smith           ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr);
3475fa7ec2dSBarry Smith           ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr);
348d34fcf5fSBarry Smith           ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
349b4319ba4SBarry Smith         }
350b4319ba4SBarry Smith       }
351b4319ba4SBarry Smith     }
352b4319ba4SBarry Smith     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
353b4319ba4SBarry Smith     ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr);
354b4319ba4SBarry Smith   }
355b4319ba4SBarry Smith 
356b4319ba4SBarry Smith   PetscFunctionReturn(0);
357b4319ba4SBarry Smith }
358b4319ba4SBarry Smith 
359b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
360b4319ba4SBarry Smith /*
361b4319ba4SBarry Smith    PCISDestroy -
362b4319ba4SBarry Smith */
363b4319ba4SBarry Smith #undef __FUNCT__
364b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy"
3657087cfbeSBarry Smith PetscErrorCode  PCISDestroy(PC pc)
366b4319ba4SBarry Smith {
367b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
368dfbe8321SBarry Smith   PetscErrorCode ierr;
369b4319ba4SBarry Smith 
370b4319ba4SBarry Smith   PetscFunctionBegin;
371fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr);
372fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr);
373fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr);
374fcfd50ebSBarry Smith   ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr);
375fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
376fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
377fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
378fcfd50ebSBarry Smith   ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
379fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->D);CHKERRQ(ierr);
380fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr);
381fcfd50ebSBarry Smith   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
382fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr);
383fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr);
384fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr);
385fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr);
386fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr);
387df187020SStefano Zampini   ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr);
388fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr);
389fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr);
390fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr);
391fcfd50ebSBarry Smith   ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr);
392fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr);
393fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr);
39416909a7fSStefano Zampini   ierr = VecScatterDestroy(&pcis->N_to_D);CHKERRQ(ierr);
395fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr);
39605b42c5fSBarry Smith   ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);
3977dbfca69SStefano Zampini   if (pcis->n_neigh > -1) {
3988d12c799SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr);
3997dbfca69SStefano Zampini   }
4008d12c799SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr);
4015e8657edSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr);
402bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr);
403bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr);
404bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr);
405b4319ba4SBarry Smith   PetscFunctionReturn(0);
406b4319ba4SBarry Smith }
407b4319ba4SBarry Smith 
408b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
409b4319ba4SBarry Smith /*
410b4319ba4SBarry Smith    PCISCreate -
411b4319ba4SBarry Smith */
412b4319ba4SBarry Smith #undef __FUNCT__
413b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate"
4147087cfbeSBarry Smith PetscErrorCode  PCISCreate(PC pc)
415b4319ba4SBarry Smith {
416b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
417831a100dSStefano Zampini   PetscErrorCode ierr;
418b4319ba4SBarry Smith 
419b4319ba4SBarry Smith   PetscFunctionBegin;
420b4319ba4SBarry Smith   pcis->is_B_local       = 0;
421b4319ba4SBarry Smith   pcis->is_I_local       = 0;
422b4319ba4SBarry Smith   pcis->is_B_global      = 0;
423b4319ba4SBarry Smith   pcis->is_I_global      = 0;
424b4319ba4SBarry Smith   pcis->A_II             = 0;
425b4319ba4SBarry Smith   pcis->A_IB             = 0;
426b4319ba4SBarry Smith   pcis->A_BI             = 0;
427b4319ba4SBarry Smith   pcis->A_BB             = 0;
428b4319ba4SBarry Smith   pcis->D                = 0;
429b4319ba4SBarry Smith   pcis->ksp_N            = 0;
430b4319ba4SBarry Smith   pcis->ksp_D            = 0;
431b4319ba4SBarry Smith   pcis->vec1_N           = 0;
432b4319ba4SBarry Smith   pcis->vec2_N           = 0;
433b4319ba4SBarry Smith   pcis->vec1_D           = 0;
434b4319ba4SBarry Smith   pcis->vec2_D           = 0;
435b4319ba4SBarry Smith   pcis->vec3_D           = 0;
436b4319ba4SBarry Smith   pcis->vec1_B           = 0;
437b4319ba4SBarry Smith   pcis->vec2_B           = 0;
438b4319ba4SBarry Smith   pcis->vec3_B           = 0;
439b4319ba4SBarry Smith   pcis->vec1_global      = 0;
440b4319ba4SBarry Smith   pcis->work_N           = 0;
441b4319ba4SBarry Smith   pcis->global_to_D      = 0;
442b4319ba4SBarry Smith   pcis->N_to_B           = 0;
44316909a7fSStefano Zampini   pcis->N_to_D           = 0;
444b4319ba4SBarry Smith   pcis->global_to_B      = 0;
4458d12c799SStefano Zampini   pcis->mapping          = 0;
4465e8657edSStefano Zampini   pcis->BtoNmap          = 0;
4477dbfca69SStefano Zampini   pcis->n_neigh          = -1;
448831a100dSStefano Zampini   pcis->scaling_factor   = 1.0;
4493975b054SStefano Zampini   pcis->reusesubmatrices = PETSC_TRUE;
450831a100dSStefano Zampini   /* composing functions */
451bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr);
452bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr);
453bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr);
454b4319ba4SBarry Smith   PetscFunctionReturn(0);
455b4319ba4SBarry Smith }
456b4319ba4SBarry Smith 
457b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
458b4319ba4SBarry Smith /*
459b4319ba4SBarry Smith    PCISApplySchur -
460b4319ba4SBarry Smith 
461b4319ba4SBarry Smith    Input parameters:
462b4319ba4SBarry Smith .  pc - preconditioner context
463b4319ba4SBarry Smith .  v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null)
464b4319ba4SBarry Smith 
465b4319ba4SBarry Smith    Output parameters:
466b4319ba4SBarry Smith .  vec1_B - result of Schur complement applied to chunk
467b4319ba4SBarry Smith .  vec2_B - garbage (used as work space), or null (and v is used as workspace)
468b4319ba4SBarry Smith .  vec1_D - garbage (used as work space)
469b4319ba4SBarry Smith .  vec2_D - garbage (used as work space)
470b4319ba4SBarry Smith 
471b4319ba4SBarry Smith */
472b4319ba4SBarry Smith #undef __FUNCT__
473831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur"
4747087cfbeSBarry Smith PetscErrorCode  PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
475b4319ba4SBarry Smith {
476dfbe8321SBarry Smith   PetscErrorCode ierr;
477b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
478b4319ba4SBarry Smith 
479b4319ba4SBarry Smith   PetscFunctionBegin;
4802fa5cd67SKarl Rupp   if (!vec2_B) vec2_B = v;
481b4319ba4SBarry Smith 
482b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
483b4319ba4SBarry Smith   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
48423ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
485b4319ba4SBarry Smith   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
486efb30889SBarry Smith   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
487b4319ba4SBarry Smith   PetscFunctionReturn(0);
488b4319ba4SBarry Smith }
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
491b4319ba4SBarry Smith /*
492b4319ba4SBarry Smith    PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface,
493b4319ba4SBarry Smith    including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE
494b4319ba4SBarry Smith    mode.
495b4319ba4SBarry Smith 
496b4319ba4SBarry Smith    Input parameters:
497b4319ba4SBarry Smith .  pc - preconditioner context
498b4319ba4SBarry Smith .  array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector
499b4319ba4SBarry Smith .  v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array
500b4319ba4SBarry Smith 
501b4319ba4SBarry Smith    Output parameter:
502b4319ba4SBarry Smith .  array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector
503b4319ba4SBarry Smith .  v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array
504b4319ba4SBarry Smith 
505b4319ba4SBarry Smith    Notes:
506b4319ba4SBarry Smith    The entries in the array that do not correspond to interface nodes remain unaltered.
507b4319ba4SBarry Smith */
508b4319ba4SBarry Smith #undef __FUNCT__
509b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB"
5107087cfbeSBarry Smith PetscErrorCode  PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc)
511b4319ba4SBarry Smith {
5125d0c19d7SBarry Smith   PetscInt       i;
5135d0c19d7SBarry Smith   const PetscInt *idex;
5146849ba73SBarry Smith   PetscErrorCode ierr;
515b4319ba4SBarry Smith   PetscScalar    *array_B;
516b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
517b4319ba4SBarry Smith 
518b4319ba4SBarry Smith   PetscFunctionBegin;
519b4319ba4SBarry Smith   ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr);
520b4319ba4SBarry Smith   ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
521b4319ba4SBarry Smith 
522b4319ba4SBarry Smith   if (smode == SCATTER_FORWARD) {
523b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5242fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]];
525b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5262fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]];
527b4319ba4SBarry Smith     }
528b4319ba4SBarry Smith   } else {  /* SCATTER_REVERSE */
529b4319ba4SBarry Smith     if (imode == INSERT_VALUES) {
5302fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i];
531b4319ba4SBarry Smith     } else {  /* ADD_VALUES */
5322fa5cd67SKarl Rupp       for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i];
533b4319ba4SBarry Smith     }
534b4319ba4SBarry Smith   }
535b4319ba4SBarry Smith   ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr);
536b4319ba4SBarry Smith   ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr);
537b4319ba4SBarry Smith   PetscFunctionReturn(0);
538b4319ba4SBarry Smith }
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */
541b4319ba4SBarry Smith /*
542b4319ba4SBarry Smith    PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement.
543b4319ba4SBarry Smith    More precisely, solves the problem:
544b4319ba4SBarry Smith                                         [ A_II  A_IB ] [ . ]   [ 0 ]
545b4319ba4SBarry Smith                                         [            ] [   ] = [   ]
546b4319ba4SBarry Smith                                         [ A_BI  A_BB ] [ x ]   [ b ]
547b4319ba4SBarry Smith 
548b4319ba4SBarry Smith    Input parameters:
549b4319ba4SBarry Smith .  pc - preconditioner context
550b4319ba4SBarry Smith .  b - vector of local interface nodes (including ghosts)
551b4319ba4SBarry Smith 
552b4319ba4SBarry Smith    Output parameters:
553b4319ba4SBarry Smith .  x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur
554b4319ba4SBarry Smith        complement to b
555b4319ba4SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
556b4319ba4SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
557b4319ba4SBarry Smith 
558b4319ba4SBarry Smith */
559b4319ba4SBarry Smith #undef __FUNCT__
560b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur"
5617087cfbeSBarry Smith PetscErrorCode  PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N)
562b4319ba4SBarry Smith {
563dfbe8321SBarry Smith   PetscErrorCode ierr;
564b4319ba4SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
565b4319ba4SBarry Smith 
566b4319ba4SBarry Smith   PetscFunctionBegin;
567b4319ba4SBarry Smith   /*
568b4319ba4SBarry Smith     Neumann solvers.
569b4319ba4SBarry Smith     Applying the inverse of the local Schur complement, i.e, solving a Neumann
570b4319ba4SBarry Smith     Problem with zero at the interior nodes of the RHS and extracting the interface
571b4319ba4SBarry Smith     part of the solution. inverse Schur complement is applied to b and the result
572b4319ba4SBarry Smith     is stored in x.
573b4319ba4SBarry Smith   */
574b4319ba4SBarry Smith   /* Setting the RHS vec1_N */
575efb30889SBarry Smith   ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr);
576ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
577ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
578b4319ba4SBarry Smith   /* Checking for consistency of the RHS */
579b4319ba4SBarry Smith   {
580ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
581c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr);
582b4319ba4SBarry Smith     if (flg) {
583b4319ba4SBarry Smith       PetscScalar average;
5843050cee2SBarry Smith       PetscViewer viewer;
585ce94432eSBarry Smith       ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr);
5863050cee2SBarry Smith 
587b4319ba4SBarry Smith       ierr    = VecSum(vec1_N,&average);CHKERRQ(ierr);
588b4319ba4SBarry Smith       average = average / ((PetscReal)pcis->n);
5891575c14dSBarry Smith       ierr    = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
590b4319ba4SBarry Smith       if (pcis->pure_neumann) {
5919f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
592b4319ba4SBarry Smith       } else {
5939f73f8ecSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed.    Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr);
594b4319ba4SBarry Smith       }
5953050cee2SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5961575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
597b4319ba4SBarry Smith     }
598b4319ba4SBarry Smith   }
599b4319ba4SBarry Smith   /* Solving the system for vec2_N */
60023ce1328SBarry Smith   ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr);
601b4319ba4SBarry Smith   /* Extracting the local interface vector out of the solution */
602ca9f406cSSatish Balay   ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603ca9f406cSSatish Balay   ierr = VecScatterEnd  (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
604b4319ba4SBarry Smith   PetscFunctionReturn(0);
605b4319ba4SBarry Smith }
606b4319ba4SBarry Smith 
607b4319ba4SBarry Smith 
608b4319ba4SBarry Smith 
609b4319ba4SBarry Smith 
610b4319ba4SBarry Smith 
611b4319ba4SBarry Smith 
612b4319ba4SBarry Smith 
613b4319ba4SBarry Smith 
614b4319ba4SBarry Smith 
615