123ce1328SBarry Smith 2aaa7dc30SBarry Smith #include <../src/ksp/pc/impls/is/pcis.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 27b965d45eSStefano Zampini .seealso: PCBDDC 28b965d45eSStefano Zampini @*/ 29b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use) 30b965d45eSStefano Zampini { 31b965d45eSStefano Zampini PetscErrorCode ierr; 32b965d45eSStefano Zampini 33b965d45eSStefano Zampini PetscFunctionBegin; 34b965d45eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 35b965d45eSStefano Zampini ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr); 36b965d45eSStefano Zampini PetscFunctionReturn(0); 37b965d45eSStefano Zampini } 38b965d45eSStefano Zampini 39ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors) 40ba1573a8SStefano Zampini { 41ba1573a8SStefano Zampini PetscErrorCode ierr; 42ba1573a8SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 43ba1573a8SStefano Zampini 44ba1573a8SStefano Zampini PetscFunctionBegin; 45ba1573a8SStefano Zampini ierr = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr); 46bf83c2afSStefano Zampini ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 47ba1573a8SStefano Zampini pcis->D = scaling_factors; 48ba1573a8SStefano Zampini PetscFunctionReturn(0); 49ba1573a8SStefano Zampini } 50ba1573a8SStefano Zampini 51ba1573a8SStefano Zampini /*@ 52ba1573a8SStefano Zampini PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS. 53ba1573a8SStefano Zampini 54ba1573a8SStefano Zampini Not collective 55ba1573a8SStefano Zampini 56ba1573a8SStefano Zampini Input Parameters: 57ba1573a8SStefano Zampini + pc - the preconditioning context 58ba1573a8SStefano Zampini - scaling_factors - scaling factors for the subdomain 59ba1573a8SStefano Zampini 60ba1573a8SStefano Zampini Level: intermediate 61ba1573a8SStefano Zampini 62ba1573a8SStefano Zampini Notes: 63ba1573a8SStefano Zampini Intended to use with jumping coefficients cases. 64ba1573a8SStefano Zampini 65ba1573a8SStefano Zampini .seealso: PCBDDC 66ba1573a8SStefano Zampini @*/ 67ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors) 68ba1573a8SStefano Zampini { 69ba1573a8SStefano Zampini PetscErrorCode ierr; 70ba1573a8SStefano Zampini 71ba1573a8SStefano Zampini PetscFunctionBegin; 72ba1573a8SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 73ba1573a8SStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr); 74ba1573a8SStefano Zampini PetscFunctionReturn(0); 75ba1573a8SStefano Zampini } 76ba1573a8SStefano Zampini 77831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal) 78831a100dSStefano Zampini { 79831a100dSStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 80831a100dSStefano Zampini 81831a100dSStefano Zampini PetscFunctionBegin; 82831a100dSStefano Zampini pcis->scaling_factor = scal; 83831a100dSStefano Zampini PetscFunctionReturn(0); 84831a100dSStefano Zampini } 85831a100dSStefano Zampini 86831a100dSStefano Zampini /*@ 87831a100dSStefano Zampini PCISSetSubdomainScalingFactor - Set scaling factor for PCIS. 88831a100dSStefano Zampini 89831a100dSStefano Zampini Not collective 90831a100dSStefano Zampini 91831a100dSStefano Zampini Input Parameters: 92831a100dSStefano Zampini + pc - the preconditioning context 93831a100dSStefano Zampini - scal - scaling factor for the subdomain 94831a100dSStefano Zampini 95831a100dSStefano Zampini Level: intermediate 96831a100dSStefano Zampini 97831a100dSStefano Zampini Notes: 98831a100dSStefano Zampini Intended to use with jumping coefficients cases. 99831a100dSStefano Zampini 100831a100dSStefano Zampini .seealso: PCBDDC 101831a100dSStefano Zampini @*/ 102831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal) 103831a100dSStefano Zampini { 104831a100dSStefano Zampini PetscErrorCode ierr; 105831a100dSStefano Zampini 106831a100dSStefano Zampini PetscFunctionBegin; 107831a100dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 108831a100dSStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr); 109831a100dSStefano Zampini PetscFunctionReturn(0); 110831a100dSStefano Zampini } 111831a100dSStefano Zampini 112b4319ba4SBarry Smith 113b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 114b4319ba4SBarry Smith /* 115b4319ba4SBarry Smith PCISSetUp - 116b4319ba4SBarry Smith */ 117d9869140SStefano Zampini PetscErrorCode PCISSetUp(PC pc, PetscBool computematrices, PetscBool computesolvers) 118b4319ba4SBarry Smith { 119b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 120bf327c11SStefano Zampini Mat_IS *matis; 1215e8657edSStefano Zampini MatReuse reuse; 1226849ba73SBarry Smith PetscErrorCode ierr; 12385c21eb1SStefano Zampini PetscBool flg,issbaij; 124b4319ba4SBarry Smith 125b4319ba4SBarry Smith PetscFunctionBegin; 12685c21eb1SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr); 127ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 12885c21eb1SStefano Zampini matis = (Mat_IS*)pc->pmat->data; 129b4319ba4SBarry Smith 1305e8657edSStefano Zampini /* first time creation, get info on substructuring */ 1315e8657edSStefano Zampini if (!pc->setupcalled) { 1325e8657edSStefano Zampini PetscInt n_I; 1335e8657edSStefano Zampini PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 1346f516dd7SStefano Zampini PetscBT bt; 1355e8657edSStefano Zampini PetscInt i,j; 136b4319ba4SBarry Smith 137892d8026SStefano Zampini /* get info on mapping */ 1383bbff08aSStefano Zampini ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr); 1398d12c799SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 1403bbff08aSStefano Zampini pcis->mapping = pc->pmat->rmap->mapping; 1418d12c799SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr); 1428d12c799SStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 143892d8026SStefano Zampini 144b4319ba4SBarry Smith /* Identifying interior and interface nodes, in local numbering */ 1456f516dd7SStefano Zampini ierr = PetscBTCreate(pcis->n,&bt);CHKERRQ(ierr); 1464843afd6SStefano Zampini for (i=0;i<pcis->n_neigh;i++) 1476f516dd7SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 1486f516dd7SStefano Zampini ierr = PetscBTSet(bt,pcis->shared[i][j]);CHKERRQ(ierr); 1496f516dd7SStefano Zampini } 150892d8026SStefano Zampini 1515e8657edSStefano Zampini /* Creating local and global index sets for interior and inteface nodes. */ 152785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr); 153785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr); 154b4319ba4SBarry Smith for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 1556f516dd7SStefano Zampini if (!PetscBTLookup(bt,i)) { 1562fa5cd67SKarl Rupp idx_I_local[n_I] = i; 1572fa5cd67SKarl Rupp n_I++; 1582fa5cd67SKarl Rupp } else { 1592fa5cd67SKarl Rupp idx_B_local[pcis->n_B] = i; 1602fa5cd67SKarl Rupp pcis->n_B++; 1612fa5cd67SKarl Rupp } 162b4319ba4SBarry Smith } 1636f516dd7SStefano Zampini 164b4319ba4SBarry Smith /* Getting the global numbering */ 165b4319ba4SBarry Smith idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 166b4319ba4SBarry Smith idx_I_global = idx_B_local + pcis->n_B; 1678d12c799SStefano Zampini ierr = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 1688d12c799SStefano Zampini ierr = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr); 1692fa5cd67SKarl Rupp 1705e8657edSStefano Zampini /* Creating the index sets */ 171bb943df9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); 172e850f338SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); 173bb943df9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); 174e850f338SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); 1752fa5cd67SKarl Rupp 1765e8657edSStefano Zampini /* Freeing memory */ 177b4319ba4SBarry Smith ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 178b4319ba4SBarry Smith ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 1796f516dd7SStefano Zampini ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 180b4319ba4SBarry Smith 1815e8657edSStefano Zampini /* Creating work vectors and arrays */ 182892d8026SStefano Zampini ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 183b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 184*f08590b7SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&pcis->vec1_D);CHKERRQ(ierr); 185*f08590b7SStefano Zampini ierr = VecSetSizes(pcis->vec1_D,pcis->n-pcis->n_B,PETSC_DECIDE);CHKERRQ(ierr); 186*f08590b7SStefano Zampini ierr = VecSetType(pcis->vec1_D,((PetscObject)pcis->vec1_N)->type_name);CHKERRQ(ierr); 187b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 188b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 189df187020SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr); 190*f08590b7SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&pcis->vec1_B);CHKERRQ(ierr); 191*f08590b7SStefano Zampini ierr = VecSetSizes(pcis->vec1_B,pcis->n_B,PETSC_DECIDE);CHKERRQ(ierr); 192*f08590b7SStefano Zampini ierr = VecSetType(pcis->vec1_B,((PetscObject)pcis->vec1_N)->type_name);CHKERRQ(ierr); 193b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 194b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 1952a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 196854ce69bSBarry Smith ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr); 1975e8657edSStefano Zampini /* scaling vector */ 198bf83c2afSStefano Zampini if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */ 1995e8657edSStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 200bf83c2afSStefano Zampini ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 201bf83c2afSStefano Zampini } 202b4319ba4SBarry Smith 203b4319ba4SBarry Smith /* Creating the scatter contexts */ 20416909a7fSStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr); 20523ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 206b4319ba4SBarry Smith ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 20723ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 208b4319ba4SBarry Smith 2095e8657edSStefano Zampini /* map from boundary to local */ 2105e8657edSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr); 2115e8657edSStefano Zampini } 2125e8657edSStefano Zampini 2135e8657edSStefano Zampini /* 2145e8657edSStefano Zampini Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 2155e8657edSStefano Zampini is such that interior nodes come first than the interface ones, we have 2165e8657edSStefano Zampini 2175e8657edSStefano Zampini [ A_II | A_IB ] 2185e8657edSStefano Zampini A = [------+------] 2195e8657edSStefano Zampini [ A_BI | A_BB ] 2205e8657edSStefano Zampini */ 221d9869140SStefano Zampini if (computematrices) { 2225e8657edSStefano Zampini reuse = MAT_INITIAL_MATRIX; 2233975b054SStefano Zampini if (pcis->reusesubmatrices && pc->setupcalled) { 2245e8657edSStefano Zampini if (pc->flag == SAME_NONZERO_PATTERN) { 2255e8657edSStefano Zampini reuse = MAT_REUSE_MATRIX; 2265e8657edSStefano Zampini } else { 2275e8657edSStefano Zampini reuse = MAT_INITIAL_MATRIX; 2285e8657edSStefano Zampini } 2295e8657edSStefano Zampini } 2305e8657edSStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 2315e8657edSStefano Zampini ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 2325e8657edSStefano Zampini ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 2335e8657edSStefano Zampini ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 2345e8657edSStefano Zampini ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 2355e8657edSStefano Zampini } 2365e8657edSStefano Zampini 2377dae84e0SHong Zhang ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr); 2387dae84e0SHong Zhang ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr); 2395e8657edSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 2405e8657edSStefano Zampini if (!issbaij) { 2417dae84e0SHong Zhang ierr = MatCreateSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 2427dae84e0SHong Zhang ierr = MatCreateSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 2435e8657edSStefano Zampini } else { 2445e8657edSStefano Zampini Mat newmat; 245d9869140SStefano Zampini 2465e8657edSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr); 2477dae84e0SHong Zhang ierr = MatCreateSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 2487dae84e0SHong Zhang ierr = MatCreateSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 2495e8657edSStefano Zampini ierr = MatDestroy(&newmat);CHKERRQ(ierr); 2505e8657edSStefano Zampini } 251d9869140SStefano Zampini } 2525e8657edSStefano Zampini 253bf83c2afSStefano Zampini /* Creating scaling vector D */ 254c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); 255bf83c2afSStefano Zampini if (pcis->use_stiffness_scaling) { 25615f235b8SStefano Zampini PetscScalar *a; 25715f235b8SStefano Zampini PetscInt i,n; 25815f235b8SStefano Zampini 259d9869140SStefano Zampini if (pcis->A_BB) { 260bf83c2afSStefano Zampini ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr); 261d9869140SStefano Zampini } else { 262d9869140SStefano Zampini ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr); 263d9869140SStefano Zampini ierr = VecScatterBegin(pcis->N_to_D,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 264d9869140SStefano Zampini ierr = VecScatterEnd(pcis->N_to_D,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 265d9869140SStefano Zampini } 26615f235b8SStefano Zampini ierr = VecGetLocalSize(pcis->D,&n);CHKERRQ(ierr); 26715f235b8SStefano Zampini ierr = VecGetArray(pcis->D,&a);CHKERRQ(ierr); 26815f235b8SStefano Zampini for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0; 26915f235b8SStefano Zampini ierr = VecRestoreArray(pcis->D,&a);CHKERRQ(ierr); 270ba1573a8SStefano Zampini } 271d9869140SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,matis->counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 272d9869140SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,matis->counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 273ba1573a8SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 274b4319ba4SBarry Smith 275b4319ba4SBarry Smith /* See historical note 01, at the bottom of this file. */ 276b4319ba4SBarry Smith 2775e8657edSStefano Zampini /* Creating the KSP contexts for the local Dirichlet and Neumann problems */ 2785e8657edSStefano Zampini if (computesolvers) { 279b4319ba4SBarry Smith PC pc_ctx; 2805e8657edSStefano Zampini 2815e8657edSStefano Zampini pcis->pure_neumann = matis->pure_neumann; 282b4319ba4SBarry Smith /* Dirichlet */ 283b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 284422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr); 2851cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 28623ee1639SBarry Smith ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr); 287b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 288b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 289b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 290b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 291b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 292b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 293b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 294b4319ba4SBarry Smith /* Neumann */ 295b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 296422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr); 2971cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); 29823ee1639SBarry Smith ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr); 299b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 300b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 301b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 302b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 303b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 304b4319ba4SBarry Smith { 305ace3abfcSBarry Smith PetscBool damp_fixed = PETSC_FALSE, 30690d69ab7SBarry Smith remove_nullspace_fixed = PETSC_FALSE, 30790d69ab7SBarry Smith set_damping_factor_floating = PETSC_FALSE, 30890d69ab7SBarry Smith not_damp_floating = PETSC_FALSE, 30990d69ab7SBarry Smith not_remove_nullspace_floating = PETSC_FALSE; 310b4319ba4SBarry Smith PetscReal fixed_factor, 311b4319ba4SBarry Smith floating_factor; 312b4319ba4SBarry Smith 313c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 3142fa5cd67SKarl Rupp if (!damp_fixed) fixed_factor = 0.0; 315c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr); 316b4319ba4SBarry Smith 317c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr); 318b4319ba4SBarry Smith 319c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", 320b4319ba4SBarry Smith &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 3212fa5cd67SKarl Rupp if (!set_damping_factor_floating) floating_factor = 0.0; 322c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr); 3232fa5cd67SKarl Rupp if (!set_damping_factor_floating) floating_factor = 1.e-12; 324b4319ba4SBarry Smith 325c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,NULL);CHKERRQ(ierr); 326b4319ba4SBarry Smith 327c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,NULL);CHKERRQ(ierr); 328b4319ba4SBarry Smith 329b4319ba4SBarry Smith if (pcis->pure_neumann) { /* floating subdomain */ 330b4319ba4SBarry Smith if (!(not_damp_floating)) { 331c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 332c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 333b4319ba4SBarry Smith } 334b4319ba4SBarry Smith if (!(not_remove_nullspace_floating)) { 335b4319ba4SBarry Smith MatNullSpace nullsp; 3360298fd71SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 3375fa7ec2dSBarry Smith ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 338d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 339b4319ba4SBarry Smith } 340b4319ba4SBarry Smith } else { /* fixed subdomain */ 341b4319ba4SBarry Smith if (damp_fixed) { 342c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 343c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 344b4319ba4SBarry Smith } 345b4319ba4SBarry Smith if (remove_nullspace_fixed) { 346b4319ba4SBarry Smith MatNullSpace nullsp; 3470298fd71SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 3485fa7ec2dSBarry Smith ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 349d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 350b4319ba4SBarry Smith } 351b4319ba4SBarry Smith } 352b4319ba4SBarry Smith } 353b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 354b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 355b4319ba4SBarry Smith } 356b4319ba4SBarry Smith PetscFunctionReturn(0); 357b4319ba4SBarry Smith } 358b4319ba4SBarry Smith 359b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 360b4319ba4SBarry Smith /* 361b4319ba4SBarry Smith PCISDestroy - 362b4319ba4SBarry Smith */ 3637087cfbeSBarry Smith PetscErrorCode PCISDestroy(PC pc) 364b4319ba4SBarry Smith { 365b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 366dfbe8321SBarry Smith PetscErrorCode ierr; 367b4319ba4SBarry Smith 368b4319ba4SBarry Smith PetscFunctionBegin; 369fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 370fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 371fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 372fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 373fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 374fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 375fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 376fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 377fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 378fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 379fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 380fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 381fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 382fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 383fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 384fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 385df187020SStefano Zampini ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr); 386fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 387fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 388fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 389fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 390fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 391fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 39216909a7fSStefano Zampini ierr = VecScatterDestroy(&pcis->N_to_D);CHKERRQ(ierr); 393fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 39405b42c5fSBarry Smith ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 3957dbfca69SStefano Zampini if (pcis->n_neigh > -1) { 3968d12c799SStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 3977dbfca69SStefano Zampini } 3988d12c799SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 3995e8657edSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr); 400bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr); 401bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr); 402bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr); 403b4319ba4SBarry Smith PetscFunctionReturn(0); 404b4319ba4SBarry Smith } 405b4319ba4SBarry Smith 406b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 407b4319ba4SBarry Smith /* 408b4319ba4SBarry Smith PCISCreate - 409b4319ba4SBarry Smith */ 4107087cfbeSBarry Smith PetscErrorCode PCISCreate(PC pc) 411b4319ba4SBarry Smith { 412b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 413831a100dSStefano Zampini PetscErrorCode ierr; 414b4319ba4SBarry Smith 415b4319ba4SBarry Smith PetscFunctionBegin; 4167dbfca69SStefano Zampini pcis->n_neigh = -1; 417831a100dSStefano Zampini pcis->scaling_factor = 1.0; 4183975b054SStefano Zampini pcis->reusesubmatrices = PETSC_TRUE; 419831a100dSStefano Zampini /* composing functions */ 420bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr); 421bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 422bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 423b4319ba4SBarry Smith PetscFunctionReturn(0); 424b4319ba4SBarry Smith } 425b4319ba4SBarry Smith 426b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 427b4319ba4SBarry Smith /* 428b4319ba4SBarry Smith PCISApplySchur - 429b4319ba4SBarry Smith 430b4319ba4SBarry Smith Input parameters: 431b4319ba4SBarry Smith . pc - preconditioner context 432b4319ba4SBarry Smith . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 433b4319ba4SBarry Smith 434b4319ba4SBarry Smith Output parameters: 435b4319ba4SBarry Smith . vec1_B - result of Schur complement applied to chunk 436b4319ba4SBarry Smith . vec2_B - garbage (used as work space), or null (and v is used as workspace) 437b4319ba4SBarry Smith . vec1_D - garbage (used as work space) 438b4319ba4SBarry Smith . vec2_D - garbage (used as work space) 439b4319ba4SBarry Smith 440b4319ba4SBarry Smith */ 4417087cfbeSBarry Smith PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 442b4319ba4SBarry Smith { 443dfbe8321SBarry Smith PetscErrorCode ierr; 444b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 445b4319ba4SBarry Smith 446b4319ba4SBarry Smith PetscFunctionBegin; 4472fa5cd67SKarl Rupp if (!vec2_B) vec2_B = v; 448b4319ba4SBarry Smith 449b4319ba4SBarry Smith ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 450b4319ba4SBarry Smith ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 45123ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 452b4319ba4SBarry Smith ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 453efb30889SBarry Smith ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 454b4319ba4SBarry Smith PetscFunctionReturn(0); 455b4319ba4SBarry Smith } 456b4319ba4SBarry Smith 457b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 458b4319ba4SBarry Smith /* 459b4319ba4SBarry Smith PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 460b4319ba4SBarry Smith including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 461b4319ba4SBarry Smith mode. 462b4319ba4SBarry Smith 463b4319ba4SBarry Smith Input parameters: 464b4319ba4SBarry Smith . pc - preconditioner context 465b4319ba4SBarry Smith . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 466b4319ba4SBarry Smith . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 467b4319ba4SBarry Smith 468b4319ba4SBarry Smith Output parameter: 469b4319ba4SBarry Smith . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 470b4319ba4SBarry Smith . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 471b4319ba4SBarry Smith 472b4319ba4SBarry Smith Notes: 473b4319ba4SBarry Smith The entries in the array that do not correspond to interface nodes remain unaltered. 474b4319ba4SBarry Smith */ 4757087cfbeSBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 476b4319ba4SBarry Smith { 4775d0c19d7SBarry Smith PetscInt i; 4785d0c19d7SBarry Smith const PetscInt *idex; 4796849ba73SBarry Smith PetscErrorCode ierr; 480b4319ba4SBarry Smith PetscScalar *array_B; 481b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 482b4319ba4SBarry Smith 483b4319ba4SBarry Smith PetscFunctionBegin; 484b4319ba4SBarry Smith ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 485b4319ba4SBarry Smith ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 486b4319ba4SBarry Smith 487b4319ba4SBarry Smith if (smode == SCATTER_FORWARD) { 488b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 4892fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]]; 490b4319ba4SBarry Smith } else { /* ADD_VALUES */ 4912fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]]; 492b4319ba4SBarry Smith } 493b4319ba4SBarry Smith } else { /* SCATTER_REVERSE */ 494b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 4952fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i]; 496b4319ba4SBarry Smith } else { /* ADD_VALUES */ 4972fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i]; 498b4319ba4SBarry Smith } 499b4319ba4SBarry Smith } 500b4319ba4SBarry Smith ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 501b4319ba4SBarry Smith ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 502b4319ba4SBarry Smith PetscFunctionReturn(0); 503b4319ba4SBarry Smith } 504b4319ba4SBarry Smith 505b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 506b4319ba4SBarry Smith /* 507b4319ba4SBarry Smith PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 508b4319ba4SBarry Smith More precisely, solves the problem: 509b4319ba4SBarry Smith [ A_II A_IB ] [ . ] [ 0 ] 510b4319ba4SBarry Smith [ ] [ ] = [ ] 511b4319ba4SBarry Smith [ A_BI A_BB ] [ x ] [ b ] 512b4319ba4SBarry Smith 513b4319ba4SBarry Smith Input parameters: 514b4319ba4SBarry Smith . pc - preconditioner context 515b4319ba4SBarry Smith . b - vector of local interface nodes (including ghosts) 516b4319ba4SBarry Smith 517b4319ba4SBarry Smith Output parameters: 518b4319ba4SBarry Smith . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 519b4319ba4SBarry Smith complement to b 520b4319ba4SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 521b4319ba4SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 522b4319ba4SBarry Smith 523b4319ba4SBarry Smith */ 5247087cfbeSBarry Smith PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 525b4319ba4SBarry Smith { 526dfbe8321SBarry Smith PetscErrorCode ierr; 527b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 528b4319ba4SBarry Smith 529b4319ba4SBarry Smith PetscFunctionBegin; 530b4319ba4SBarry Smith /* 531b4319ba4SBarry Smith Neumann solvers. 532b4319ba4SBarry Smith Applying the inverse of the local Schur complement, i.e, solving a Neumann 533b4319ba4SBarry Smith Problem with zero at the interior nodes of the RHS and extracting the interface 534b4319ba4SBarry Smith part of the solution. inverse Schur complement is applied to b and the result 535b4319ba4SBarry Smith is stored in x. 536b4319ba4SBarry Smith */ 537b4319ba4SBarry Smith /* Setting the RHS vec1_N */ 538efb30889SBarry Smith ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 539ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 540ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 541b4319ba4SBarry Smith /* Checking for consistency of the RHS */ 542b4319ba4SBarry Smith { 543ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 544c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr); 545b4319ba4SBarry Smith if (flg) { 546b4319ba4SBarry Smith PetscScalar average; 5473050cee2SBarry Smith PetscViewer viewer; 548ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 5493050cee2SBarry Smith 550b4319ba4SBarry Smith ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 551b4319ba4SBarry Smith average = average / ((PetscReal)pcis->n); 5521575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 553b4319ba4SBarry Smith if (pcis->pure_neumann) { 5549f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 555b4319ba4SBarry Smith } else { 5569f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 557b4319ba4SBarry Smith } 5583050cee2SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5591575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 560b4319ba4SBarry Smith } 561b4319ba4SBarry Smith } 562b4319ba4SBarry Smith /* Solving the system for vec2_N */ 56323ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 564b4319ba4SBarry Smith /* Extracting the local interface vector out of the solution */ 565ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 567b4319ba4SBarry Smith PetscFunctionReturn(0); 568b4319ba4SBarry Smith } 569