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 = VecDestroy(&pcis->D);CHKERRQ(ierr); 52ba1573a8SStefano Zampini ierr = PetscObjectReference((PetscObject)scaling_factors);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" 1317087cfbeSBarry Smith PetscErrorCode PCISSetUp(PC pc) 132b4319ba4SBarry Smith { 133b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 134bf327c11SStefano Zampini Mat_IS *matis; 1356849ba73SBarry Smith PetscErrorCode ierr; 13685c21eb1SStefano Zampini PetscBool flg,issbaij; 137831a100dSStefano Zampini Vec counter; 138b4319ba4SBarry Smith 139b4319ba4SBarry Smith PetscFunctionBegin; 14085c21eb1SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr); 141ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 14285c21eb1SStefano Zampini matis = (Mat_IS*)pc->pmat->data; 143b4319ba4SBarry Smith 144b4319ba4SBarry Smith pcis->pure_neumann = matis->pure_neumann; 145b4319ba4SBarry Smith 146892d8026SStefano Zampini /* get info on mapping */ 1478d12c799SStefano Zampini ierr = PetscObjectReference((PetscObject)matis->mapping);CHKERRQ(ierr); 1488d12c799SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 1498d12c799SStefano Zampini pcis->mapping = matis->mapping; 1508d12c799SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr); 1518d12c799SStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 152892d8026SStefano Zampini 153892d8026SStefano Zampini /* Creating local and global index sets for interior and inteface nodes. */ 154b4319ba4SBarry Smith { 15513f74950SBarry Smith PetscInt n_I; 15613f74950SBarry Smith PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 157892d8026SStefano Zampini PetscInt *array; 158892d8026SStefano Zampini PetscInt i,j; 159892d8026SStefano Zampini 160b4319ba4SBarry Smith /* Identifying interior and interface nodes, in local numbering */ 161785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&array);CHKERRQ(ierr); 162892d8026SStefano Zampini ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr); 1634843afd6SStefano Zampini for (i=0;i<pcis->n_neigh;i++) 164892d8026SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) 165892d8026SStefano Zampini array[pcis->shared[i][j]] += 1; 166892d8026SStefano Zampini 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 184b4319ba4SBarry Smith /* 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 190b4319ba4SBarry Smith /* Freeing memory and restoring arrays */ 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 } 195b4319ba4SBarry Smith 196b4319ba4SBarry Smith /* 197b4319ba4SBarry Smith Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 198b4319ba4SBarry Smith is such that interior nodes come first than the interface ones, we have 199b4319ba4SBarry Smith 200b4319ba4SBarry Smith [ | ] 201b4319ba4SBarry Smith [ A_II | A_IB ] 202b4319ba4SBarry Smith A = [ | ] 203b4319ba4SBarry Smith [-----------+------] 204b4319ba4SBarry Smith [ A_BI | A_BB ] 205b4319ba4SBarry Smith */ 206b4319ba4SBarry Smith 2074aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr); 20885c21eb1SStefano Zampini ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 20985c21eb1SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 21085c21eb1SStefano Zampini if (!issbaij) { 2114aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 2124aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 21385c21eb1SStefano Zampini } else { 21485c21eb1SStefano Zampini Mat newmat; 21585c21eb1SStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr); 21685c21eb1SStefano Zampini ierr = MatGetSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 21785c21eb1SStefano Zampini ierr = MatGetSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 21885c21eb1SStefano Zampini ierr = MatDestroy(&newmat);CHKERRQ(ierr); 21985c21eb1SStefano Zampini } 220b4319ba4SBarry Smith /* 221b4319ba4SBarry Smith Creating work vectors and arrays 222b4319ba4SBarry Smith */ 223892d8026SStefano Zampini ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 224b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 225b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 226b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 227b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 228df187020SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr); 229b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 230b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 231b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 2322a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 233854ce69bSBarry Smith ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr); 234b4319ba4SBarry Smith 235b4319ba4SBarry Smith /* Creating the scatter contexts */ 23623ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 237b4319ba4SBarry Smith ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 23823ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 239b4319ba4SBarry Smith 240831a100dSStefano Zampini /* Creating scaling "matrix" D */ 2410298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); 242831a100dSStefano Zampini if (!pcis->D) { 243b965d45eSStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 244b965d45eSStefano Zampini if (!pcis->use_stiffness_scaling) { 245b965d45eSStefano Zampini ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 246ba1573a8SStefano Zampini } else { 247b965d45eSStefano Zampini ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr); 248b965d45eSStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 249b965d45eSStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 250ba1573a8SStefano Zampini } 251b965d45eSStefano Zampini } 252b965d45eSStefano Zampini ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 2532a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ 254831a100dSStefano Zampini ierr = VecSet(counter,0.0);CHKERRQ(ierr); 255ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 256ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 257ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 258ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 259ba1573a8SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 260892d8026SStefano Zampini ierr = VecDestroy(&counter);CHKERRQ(ierr); 261b4319ba4SBarry Smith 262b4319ba4SBarry Smith /* See historical note 01, at the bottom of this file. */ 263b4319ba4SBarry Smith 264b4319ba4SBarry Smith /* 265b4319ba4SBarry Smith Creating the KSP contexts for the local Dirichlet and Neumann problems. 266b4319ba4SBarry Smith */ 267fcd409f5SStefano Zampini if (pcis->computesolvers) { 268b4319ba4SBarry Smith PC pc_ctx; 269b4319ba4SBarry Smith /* Dirichlet */ 270b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 2711cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 27223ee1639SBarry Smith ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr); 273b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 274b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 275b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 276b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 277b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 278b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 279b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 280b4319ba4SBarry Smith /* Neumann */ 281b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 2821cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); 28323ee1639SBarry Smith ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr); 284b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 285b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 286b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 287b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 288b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 289b4319ba4SBarry Smith { 290ace3abfcSBarry Smith PetscBool damp_fixed = PETSC_FALSE, 29190d69ab7SBarry Smith remove_nullspace_fixed = PETSC_FALSE, 29290d69ab7SBarry Smith set_damping_factor_floating = PETSC_FALSE, 29390d69ab7SBarry Smith not_damp_floating = PETSC_FALSE, 29490d69ab7SBarry Smith not_remove_nullspace_floating = PETSC_FALSE; 295b4319ba4SBarry Smith PetscReal fixed_factor, 296b4319ba4SBarry Smith floating_factor; 297b4319ba4SBarry Smith 2987adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 2992fa5cd67SKarl Rupp if (!damp_fixed) fixed_factor = 0.0; 3000298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr); 301b4319ba4SBarry Smith 3020298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr); 303b4319ba4SBarry Smith 3047adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", 305b4319ba4SBarry Smith &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 3062fa5cd67SKarl Rupp if (!set_damping_factor_floating) floating_factor = 0.0; 3070298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr); 3082fa5cd67SKarl Rupp if (!set_damping_factor_floating) floating_factor = 1.e-12; 309b4319ba4SBarry Smith 3100298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,NULL);CHKERRQ(ierr); 311b4319ba4SBarry Smith 3120298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,NULL);CHKERRQ(ierr); 313b4319ba4SBarry Smith 314b4319ba4SBarry Smith if (pcis->pure_neumann) { /* floating subdomain */ 315b4319ba4SBarry Smith if (!(not_damp_floating)) { 316c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 317c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 318b4319ba4SBarry Smith } 319b4319ba4SBarry Smith if (!(not_remove_nullspace_floating)) { 320b4319ba4SBarry Smith MatNullSpace nullsp; 3210298fd71SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 322*5fa7ec2dSBarry Smith ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 323d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 324b4319ba4SBarry Smith } 325b4319ba4SBarry Smith } else { /* fixed subdomain */ 326b4319ba4SBarry Smith if (damp_fixed) { 327c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 328c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 329b4319ba4SBarry Smith } 330b4319ba4SBarry Smith if (remove_nullspace_fixed) { 331b4319ba4SBarry Smith MatNullSpace nullsp; 3320298fd71SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 333*5fa7ec2dSBarry Smith ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 334d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 335b4319ba4SBarry Smith } 336b4319ba4SBarry Smith } 337b4319ba4SBarry Smith } 338b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 339b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 340b4319ba4SBarry Smith } 341b4319ba4SBarry Smith 342b4319ba4SBarry Smith PetscFunctionReturn(0); 343b4319ba4SBarry Smith } 344b4319ba4SBarry Smith 345b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 346b4319ba4SBarry Smith /* 347b4319ba4SBarry Smith PCISDestroy - 348b4319ba4SBarry Smith */ 349b4319ba4SBarry Smith #undef __FUNCT__ 350b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy" 3517087cfbeSBarry Smith PetscErrorCode PCISDestroy(PC pc) 352b4319ba4SBarry Smith { 353b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 354dfbe8321SBarry Smith PetscErrorCode ierr; 355b4319ba4SBarry Smith 356b4319ba4SBarry Smith PetscFunctionBegin; 357fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 358fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 359fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 360fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 361fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 362fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 363fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 364fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 365fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 366fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 367fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 368fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 369fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 370fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 371fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 372fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 373df187020SStefano Zampini ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr); 374fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 375fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 376fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 377fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 378fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 379fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 380fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 38105b42c5fSBarry Smith ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 3827dbfca69SStefano Zampini if (pcis->n_neigh > -1) { 3838d12c799SStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 3847dbfca69SStefano Zampini } 3858d12c799SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 386bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr); 387bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr); 388bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr); 389b4319ba4SBarry Smith PetscFunctionReturn(0); 390b4319ba4SBarry Smith } 391b4319ba4SBarry Smith 392b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 393b4319ba4SBarry Smith /* 394b4319ba4SBarry Smith PCISCreate - 395b4319ba4SBarry Smith */ 396b4319ba4SBarry Smith #undef __FUNCT__ 397b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate" 3987087cfbeSBarry Smith PetscErrorCode PCISCreate(PC pc) 399b4319ba4SBarry Smith { 400b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 401831a100dSStefano Zampini PetscErrorCode ierr; 402b4319ba4SBarry Smith 403b4319ba4SBarry Smith PetscFunctionBegin; 404b4319ba4SBarry Smith pcis->is_B_local = 0; 405b4319ba4SBarry Smith pcis->is_I_local = 0; 406b4319ba4SBarry Smith pcis->is_B_global = 0; 407b4319ba4SBarry Smith pcis->is_I_global = 0; 408b4319ba4SBarry Smith pcis->A_II = 0; 409b4319ba4SBarry Smith pcis->A_IB = 0; 410b4319ba4SBarry Smith pcis->A_BI = 0; 411b4319ba4SBarry Smith pcis->A_BB = 0; 412b4319ba4SBarry Smith pcis->D = 0; 413b4319ba4SBarry Smith pcis->ksp_N = 0; 414b4319ba4SBarry Smith pcis->ksp_D = 0; 415b4319ba4SBarry Smith pcis->vec1_N = 0; 416b4319ba4SBarry Smith pcis->vec2_N = 0; 417b4319ba4SBarry Smith pcis->vec1_D = 0; 418b4319ba4SBarry Smith pcis->vec2_D = 0; 419b4319ba4SBarry Smith pcis->vec3_D = 0; 420b4319ba4SBarry Smith pcis->vec1_B = 0; 421b4319ba4SBarry Smith pcis->vec2_B = 0; 422b4319ba4SBarry Smith pcis->vec3_B = 0; 423b4319ba4SBarry Smith pcis->vec1_global = 0; 424b4319ba4SBarry Smith pcis->work_N = 0; 425b4319ba4SBarry Smith pcis->global_to_D = 0; 426b4319ba4SBarry Smith pcis->N_to_B = 0; 427b4319ba4SBarry Smith pcis->global_to_B = 0; 428fcd409f5SStefano Zampini pcis->computesolvers = PETSC_TRUE; 4298d12c799SStefano Zampini pcis->mapping = 0; 4307dbfca69SStefano Zampini pcis->n_neigh = -1; 4312fa5cd67SKarl Rupp 432831a100dSStefano Zampini pcis->scaling_factor = 1.0; 433831a100dSStefano Zampini /* composing functions */ 434bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr); 435bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 436bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 437b4319ba4SBarry Smith PetscFunctionReturn(0); 438b4319ba4SBarry Smith } 439b4319ba4SBarry Smith 440b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 441b4319ba4SBarry Smith /* 442b4319ba4SBarry Smith PCISApplySchur - 443b4319ba4SBarry Smith 444b4319ba4SBarry Smith Input parameters: 445b4319ba4SBarry Smith . pc - preconditioner context 446b4319ba4SBarry Smith . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 447b4319ba4SBarry Smith 448b4319ba4SBarry Smith Output parameters: 449b4319ba4SBarry Smith . vec1_B - result of Schur complement applied to chunk 450b4319ba4SBarry Smith . vec2_B - garbage (used as work space), or null (and v is used as workspace) 451b4319ba4SBarry Smith . vec1_D - garbage (used as work space) 452b4319ba4SBarry Smith . vec2_D - garbage (used as work space) 453b4319ba4SBarry Smith 454b4319ba4SBarry Smith */ 455b4319ba4SBarry Smith #undef __FUNCT__ 456831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur" 4577087cfbeSBarry Smith PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 458b4319ba4SBarry Smith { 459dfbe8321SBarry Smith PetscErrorCode ierr; 460b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 461b4319ba4SBarry Smith 462b4319ba4SBarry Smith PetscFunctionBegin; 4632fa5cd67SKarl Rupp if (!vec2_B) vec2_B = v; 464b4319ba4SBarry Smith 465b4319ba4SBarry Smith ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 466b4319ba4SBarry Smith ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 46723ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 468b4319ba4SBarry Smith ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 469efb30889SBarry Smith ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 470b4319ba4SBarry Smith PetscFunctionReturn(0); 471b4319ba4SBarry Smith } 472b4319ba4SBarry Smith 473b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 474b4319ba4SBarry Smith /* 475b4319ba4SBarry Smith PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 476b4319ba4SBarry Smith including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 477b4319ba4SBarry Smith mode. 478b4319ba4SBarry Smith 479b4319ba4SBarry Smith Input parameters: 480b4319ba4SBarry Smith . pc - preconditioner context 481b4319ba4SBarry Smith . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 482b4319ba4SBarry Smith . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 483b4319ba4SBarry Smith 484b4319ba4SBarry Smith Output parameter: 485b4319ba4SBarry Smith . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 486b4319ba4SBarry Smith . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 487b4319ba4SBarry Smith 488b4319ba4SBarry Smith Notes: 489b4319ba4SBarry Smith The entries in the array that do not correspond to interface nodes remain unaltered. 490b4319ba4SBarry Smith */ 491b4319ba4SBarry Smith #undef __FUNCT__ 492b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB" 4937087cfbeSBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 494b4319ba4SBarry Smith { 4955d0c19d7SBarry Smith PetscInt i; 4965d0c19d7SBarry Smith const PetscInt *idex; 4976849ba73SBarry Smith PetscErrorCode ierr; 498b4319ba4SBarry Smith PetscScalar *array_B; 499b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 500b4319ba4SBarry Smith 501b4319ba4SBarry Smith PetscFunctionBegin; 502b4319ba4SBarry Smith ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 503b4319ba4SBarry Smith ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 504b4319ba4SBarry Smith 505b4319ba4SBarry Smith if (smode == SCATTER_FORWARD) { 506b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 5072fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]]; 508b4319ba4SBarry Smith } else { /* ADD_VALUES */ 5092fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]]; 510b4319ba4SBarry Smith } 511b4319ba4SBarry Smith } else { /* SCATTER_REVERSE */ 512b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 5132fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i]; 514b4319ba4SBarry Smith } else { /* ADD_VALUES */ 5152fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i]; 516b4319ba4SBarry Smith } 517b4319ba4SBarry Smith } 518b4319ba4SBarry Smith ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 519b4319ba4SBarry Smith ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 520b4319ba4SBarry Smith PetscFunctionReturn(0); 521b4319ba4SBarry Smith } 522b4319ba4SBarry Smith 523b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 524b4319ba4SBarry Smith /* 525b4319ba4SBarry Smith PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 526b4319ba4SBarry Smith More precisely, solves the problem: 527b4319ba4SBarry Smith [ A_II A_IB ] [ . ] [ 0 ] 528b4319ba4SBarry Smith [ ] [ ] = [ ] 529b4319ba4SBarry Smith [ A_BI A_BB ] [ x ] [ b ] 530b4319ba4SBarry Smith 531b4319ba4SBarry Smith Input parameters: 532b4319ba4SBarry Smith . pc - preconditioner context 533b4319ba4SBarry Smith . b - vector of local interface nodes (including ghosts) 534b4319ba4SBarry Smith 535b4319ba4SBarry Smith Output parameters: 536b4319ba4SBarry Smith . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 537b4319ba4SBarry Smith complement to b 538b4319ba4SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 539b4319ba4SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 540b4319ba4SBarry Smith 541b4319ba4SBarry Smith */ 542b4319ba4SBarry Smith #undef __FUNCT__ 543b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur" 5447087cfbeSBarry Smith PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 545b4319ba4SBarry Smith { 546dfbe8321SBarry Smith PetscErrorCode ierr; 547b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 548b4319ba4SBarry Smith 549b4319ba4SBarry Smith PetscFunctionBegin; 550b4319ba4SBarry Smith /* 551b4319ba4SBarry Smith Neumann solvers. 552b4319ba4SBarry Smith Applying the inverse of the local Schur complement, i.e, solving a Neumann 553b4319ba4SBarry Smith Problem with zero at the interior nodes of the RHS and extracting the interface 554b4319ba4SBarry Smith part of the solution. inverse Schur complement is applied to b and the result 555b4319ba4SBarry Smith is stored in x. 556b4319ba4SBarry Smith */ 557b4319ba4SBarry Smith /* Setting the RHS vec1_N */ 558efb30889SBarry Smith ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 559ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 560ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 561b4319ba4SBarry Smith /* Checking for consistency of the RHS */ 562b4319ba4SBarry Smith { 563ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 5640298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr); 565b4319ba4SBarry Smith if (flg) { 566b4319ba4SBarry Smith PetscScalar average; 5673050cee2SBarry Smith PetscViewer viewer; 568ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 5693050cee2SBarry Smith 570b4319ba4SBarry Smith ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 571b4319ba4SBarry Smith average = average / ((PetscReal)pcis->n); 5727b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 573b4319ba4SBarry Smith if (pcis->pure_neumann) { 5749f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 575b4319ba4SBarry Smith } else { 5769f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 577b4319ba4SBarry Smith } 5783050cee2SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5797b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 580b4319ba4SBarry Smith } 581b4319ba4SBarry Smith } 582b4319ba4SBarry Smith /* Solving the system for vec2_N */ 58323ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 584b4319ba4SBarry Smith /* Extracting the local interface vector out of the solution */ 585ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 586ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 587b4319ba4SBarry Smith PetscFunctionReturn(0); 588b4319ba4SBarry Smith } 589b4319ba4SBarry Smith 590b4319ba4SBarry Smith 591b4319ba4SBarry Smith 592b4319ba4SBarry Smith 593b4319ba4SBarry Smith 594b4319ba4SBarry Smith 595b4319ba4SBarry Smith 596b4319ba4SBarry Smith 597b4319ba4SBarry Smith 598