123ce1328SBarry Smith 2ccd284c7SBarry Smith #include "../src/ksp/pc/impls/is/pcis.h" /*I "petscpc.h" I*/ 3831a100dSStefano Zampini 4831a100dSStefano Zampini EXTERN_C_BEGIN 5831a100dSStefano Zampini #undef __FUNCT__ 6b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling_IS" 7b965d45eSStefano Zampini static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use) 8b965d45eSStefano Zampini { 9b965d45eSStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 10b965d45eSStefano Zampini 11b965d45eSStefano Zampini PetscFunctionBegin; 12b965d45eSStefano Zampini pcis->use_stiffness_scaling = use; 13b965d45eSStefano Zampini PetscFunctionReturn(0); 14b965d45eSStefano Zampini } 15b965d45eSStefano Zampini EXTERN_C_END 16b965d45eSStefano Zampini 17b965d45eSStefano Zampini #undef __FUNCT__ 18b965d45eSStefano Zampini #define __FUNCT__ "PCISSetUseStiffnessScaling" 19b965d45eSStefano Zampini /*@ 20b965d45eSStefano Zampini PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using 21b965d45eSStefano Zampini local matrices' diagonal. 22b965d45eSStefano Zampini 23b965d45eSStefano Zampini Not collective 24b965d45eSStefano Zampini 25b965d45eSStefano Zampini Input Parameters: 26b965d45eSStefano Zampini + pc - the preconditioning context 27b965d45eSStefano Zampini - use - whether or not pcis use matrix diagonal to build partition of unity. 28b965d45eSStefano Zampini 29b965d45eSStefano Zampini Level: intermediate 30b965d45eSStefano Zampini 31b965d45eSStefano Zampini Notes: 32b965d45eSStefano Zampini 33b965d45eSStefano Zampini .seealso: PCBDDC 34b965d45eSStefano Zampini @*/ 35b965d45eSStefano Zampini PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use) 36b965d45eSStefano Zampini { 37b965d45eSStefano Zampini PetscErrorCode ierr; 38b965d45eSStefano Zampini 39b965d45eSStefano Zampini PetscFunctionBegin; 40b965d45eSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 41b965d45eSStefano Zampini ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr); 42b965d45eSStefano Zampini PetscFunctionReturn(0); 43b965d45eSStefano Zampini } 44b965d45eSStefano Zampini 45b965d45eSStefano Zampini EXTERN_C_BEGIN 46b965d45eSStefano Zampini #undef __FUNCT__ 47ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling_IS" 48ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors) 49ba1573a8SStefano Zampini { 50ba1573a8SStefano Zampini PetscErrorCode ierr; 51ba1573a8SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 52ba1573a8SStefano Zampini 53ba1573a8SStefano Zampini PetscFunctionBegin; 54ba1573a8SStefano Zampini ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 55ba1573a8SStefano Zampini ierr = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr); 56ba1573a8SStefano Zampini pcis->D = scaling_factors; 57ba1573a8SStefano Zampini PetscFunctionReturn(0); 58ba1573a8SStefano Zampini } 59ba1573a8SStefano Zampini EXTERN_C_END 60ba1573a8SStefano Zampini 61ba1573a8SStefano Zampini #undef __FUNCT__ 62ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling" 63ba1573a8SStefano Zampini /*@ 64ba1573a8SStefano Zampini PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS. 65ba1573a8SStefano Zampini 66ba1573a8SStefano Zampini Not collective 67ba1573a8SStefano Zampini 68ba1573a8SStefano Zampini Input Parameters: 69ba1573a8SStefano Zampini + pc - the preconditioning context 70ba1573a8SStefano Zampini - scaling_factors - scaling factors for the subdomain 71ba1573a8SStefano Zampini 72ba1573a8SStefano Zampini Level: intermediate 73ba1573a8SStefano Zampini 74ba1573a8SStefano Zampini Notes: 75ba1573a8SStefano Zampini Intended to use with jumping coefficients cases. 76ba1573a8SStefano Zampini 77ba1573a8SStefano Zampini .seealso: PCBDDC 78ba1573a8SStefano Zampini @*/ 79ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors) 80ba1573a8SStefano Zampini { 81ba1573a8SStefano Zampini PetscErrorCode ierr; 82ba1573a8SStefano Zampini 83ba1573a8SStefano Zampini PetscFunctionBegin; 84ba1573a8SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 85ba1573a8SStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr); 86ba1573a8SStefano Zampini PetscFunctionReturn(0); 87ba1573a8SStefano Zampini } 88ba1573a8SStefano Zampini 89ba1573a8SStefano Zampini EXTERN_C_BEGIN 90ba1573a8SStefano Zampini #undef __FUNCT__ 91831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS" 92831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal) 93831a100dSStefano Zampini { 94831a100dSStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 95831a100dSStefano Zampini 96831a100dSStefano Zampini PetscFunctionBegin; 97831a100dSStefano Zampini pcis->scaling_factor = scal; 98831a100dSStefano Zampini PetscFunctionReturn(0); 99831a100dSStefano Zampini } 100831a100dSStefano Zampini EXTERN_C_END 101831a100dSStefano Zampini 102831a100dSStefano Zampini #undef __FUNCT__ 103831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor" 104831a100dSStefano Zampini /*@ 105831a100dSStefano Zampini PCISSetSubdomainScalingFactor - Set scaling factor for PCIS. 106831a100dSStefano Zampini 107831a100dSStefano Zampini Not collective 108831a100dSStefano Zampini 109831a100dSStefano Zampini Input Parameters: 110831a100dSStefano Zampini + pc - the preconditioning context 111831a100dSStefano Zampini - scal - scaling factor for the subdomain 112831a100dSStefano Zampini 113831a100dSStefano Zampini Level: intermediate 114831a100dSStefano Zampini 115831a100dSStefano Zampini Notes: 116831a100dSStefano Zampini Intended to use with jumping coefficients cases. 117831a100dSStefano Zampini 118831a100dSStefano Zampini .seealso: PCBDDC 119831a100dSStefano Zampini @*/ 120831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal) 121831a100dSStefano Zampini { 122831a100dSStefano Zampini PetscErrorCode ierr; 123831a100dSStefano Zampini 124831a100dSStefano Zampini PetscFunctionBegin; 125831a100dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 126831a100dSStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr); 127831a100dSStefano Zampini PetscFunctionReturn(0); 128831a100dSStefano Zampini } 129831a100dSStefano Zampini 130b4319ba4SBarry Smith 131b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 132b4319ba4SBarry Smith /* 133b4319ba4SBarry Smith PCISSetUp - 134b4319ba4SBarry Smith */ 135b4319ba4SBarry Smith #undef __FUNCT__ 136b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp" 1377087cfbeSBarry Smith PetscErrorCode PCISSetUp(PC pc) 138b4319ba4SBarry Smith { 139b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 140b4319ba4SBarry Smith Mat_IS *matis = (Mat_IS*)pc->mat->data; 14113f74950SBarry Smith PetscInt i; 1426849ba73SBarry Smith PetscErrorCode ierr; 143ace3abfcSBarry Smith PetscBool flg; 144831a100dSStefano Zampini Vec counter; 145b4319ba4SBarry Smith 146b4319ba4SBarry Smith PetscFunctionBegin; 147251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr); 148*ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 149b4319ba4SBarry Smith 150b4319ba4SBarry Smith pcis->pure_neumann = matis->pure_neumann; 151b4319ba4SBarry Smith 152b4319ba4SBarry Smith /* 153b4319ba4SBarry Smith Creating the local vector vec1_N, containing the inverse of the number 154b4319ba4SBarry Smith of subdomains to which each local node (either owned or ghost) 155b4319ba4SBarry Smith pertains. To accomplish that, we scatter local vectors of 1's to 156b4319ba4SBarry Smith a global vector (adding the values); scatter the result back to 157b4319ba4SBarry Smith local vectors and finally invert the result. 158b4319ba4SBarry Smith */ 159b4319ba4SBarry Smith ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 16023ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ 161efb30889SBarry Smith ierr = VecSet(counter,0.0);CHKERRQ(ierr); 162efb30889SBarry Smith ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr); 163ca9f406cSSatish Balay ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 164ca9f406cSSatish Balay ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 165ca9f406cSSatish Balay ierr = VecScatterBegin(matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 166ca9f406cSSatish Balay ierr = VecScatterEnd (matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 167b4319ba4SBarry Smith /* 168b4319ba4SBarry Smith Creating local and global index sets for interior and 169b4319ba4SBarry Smith inteface nodes. Notice that interior nodes have D[i]==1.0. 170b4319ba4SBarry Smith */ 171b4319ba4SBarry Smith { 17213f74950SBarry Smith PetscInt n_I; 17313f74950SBarry Smith PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 174b4319ba4SBarry Smith PetscScalar *array; 175b4319ba4SBarry Smith /* Identifying interior and interface nodes, in local numbering */ 176b4319ba4SBarry Smith ierr = VecGetSize(pcis->vec1_N,&pcis->n);CHKERRQ(ierr); 177b4319ba4SBarry Smith ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 17813f74950SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);CHKERRQ(ierr); 17913f74950SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);CHKERRQ(ierr); 180b4319ba4SBarry Smith for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 1812fa5cd67SKarl Rupp if (array[i] == 1.0) { 1822fa5cd67SKarl Rupp idx_I_local[n_I] = i; 1832fa5cd67SKarl Rupp n_I++; 1842fa5cd67SKarl Rupp } else { 1852fa5cd67SKarl Rupp idx_B_local[pcis->n_B] = i; 1862fa5cd67SKarl Rupp pcis->n_B++; 1872fa5cd67SKarl Rupp } 188b4319ba4SBarry Smith } 189b4319ba4SBarry Smith /* Getting the global numbering */ 190b4319ba4SBarry Smith idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 191b4319ba4SBarry Smith idx_I_global = idx_B_local + pcis->n_B; 192b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 193b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I, idx_I_local,idx_I_global);CHKERRQ(ierr); 1942fa5cd67SKarl Rupp 195b4319ba4SBarry Smith /* Creating the index sets. */ 19670b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); 19770b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); 19870b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); 19970b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); 2002fa5cd67SKarl Rupp 201b4319ba4SBarry Smith /* Freeing memory and restoring arrays */ 202b4319ba4SBarry Smith ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 203b4319ba4SBarry Smith ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 204b4319ba4SBarry Smith ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 205b4319ba4SBarry Smith } 206b4319ba4SBarry Smith 207b4319ba4SBarry Smith /* 208b4319ba4SBarry Smith Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 209b4319ba4SBarry Smith is such that interior nodes come first than the interface ones, we have 210b4319ba4SBarry Smith 211b4319ba4SBarry Smith [ | ] 212b4319ba4SBarry Smith [ A_II | A_IB ] 213b4319ba4SBarry Smith A = [ | ] 214b4319ba4SBarry Smith [-----------+------] 215b4319ba4SBarry Smith [ A_BI | A_BB ] 216b4319ba4SBarry Smith */ 217b4319ba4SBarry Smith 2184aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr); 2194aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 2204aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 2214aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 222b4319ba4SBarry Smith 223b4319ba4SBarry Smith /* 224b4319ba4SBarry Smith Creating work vectors and arrays 225b4319ba4SBarry Smith */ 226b4319ba4SBarry Smith /* pcis->vec1_N has already been created */ 227b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 228b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 229b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 230b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 231b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 232b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 233b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 23423ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 235b4319ba4SBarry Smith ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr); 236b4319ba4SBarry Smith 237b4319ba4SBarry Smith /* Creating the scatter contexts */ 23823ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 239b4319ba4SBarry Smith ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 24023ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 241b4319ba4SBarry Smith 242831a100dSStefano Zampini /* Creating scaling "matrix" D */ 2430298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); 244831a100dSStefano Zampini if (!pcis->D) { 245b965d45eSStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 246b965d45eSStefano Zampini if (!pcis->use_stiffness_scaling) { 247b965d45eSStefano Zampini ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 248ba1573a8SStefano Zampini } else { 249b965d45eSStefano Zampini ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr); 250b965d45eSStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 251b965d45eSStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 252ba1573a8SStefano Zampini } 253b965d45eSStefano Zampini } 254b965d45eSStefano Zampini ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 255831a100dSStefano Zampini ierr = VecSet(counter,0.0);CHKERRQ(ierr); 256ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 257ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 258ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 259ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 260ba1573a8SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);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 */ 267b4319ba4SBarry Smith { 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); 272b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);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); 283b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);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); 322d8fd42c4SBarry Smith ierr = KSPSetNullSpace(pcis->ksp_N,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); 333d8fd42c4SBarry Smith ierr = KSPSetNullSpace(pcis->ksp_N,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 3427c334f02SBarry Smith ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 3432fa5cd67SKarl Rupp 344b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE; 3452fa5cd67SKarl Rupp 346831a100dSStefano Zampini ierr = VecDestroy(&counter);CHKERRQ(ierr); 347b4319ba4SBarry Smith PetscFunctionReturn(0); 348b4319ba4SBarry Smith } 349b4319ba4SBarry Smith 350b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 351b4319ba4SBarry Smith /* 352b4319ba4SBarry Smith PCISDestroy - 353b4319ba4SBarry Smith */ 354b4319ba4SBarry Smith #undef __FUNCT__ 355b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy" 3567087cfbeSBarry Smith PetscErrorCode PCISDestroy(PC pc) 357b4319ba4SBarry Smith { 358b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 359dfbe8321SBarry Smith PetscErrorCode ierr; 360b4319ba4SBarry Smith 361b4319ba4SBarry Smith PetscFunctionBegin; 362fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 363fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 364fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 365fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 366fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 367fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 368fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 369fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 370fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 371fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 372fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 373fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 374fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 375fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 376fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 377fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 378fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 379fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 380fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 381fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 382fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 383fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 384fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 38505b42c5fSBarry Smith ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 386b4319ba4SBarry Smith if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) { 387b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 388b4319ba4SBarry Smith } 3890298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetUseStiffnessScaling_C","",NULL);CHKERRQ(ierr); 3900298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainScalingFactor_C","",NULL);CHKERRQ(ierr); 3910298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C","",NULL);CHKERRQ(ierr); 392b4319ba4SBarry Smith PetscFunctionReturn(0); 393b4319ba4SBarry Smith } 394b4319ba4SBarry Smith 395b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 396b4319ba4SBarry Smith /* 397b4319ba4SBarry Smith PCISCreate - 398b4319ba4SBarry Smith */ 399b4319ba4SBarry Smith #undef __FUNCT__ 400b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate" 4017087cfbeSBarry Smith PetscErrorCode PCISCreate(PC pc) 402b4319ba4SBarry Smith { 403b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 404831a100dSStefano Zampini PetscErrorCode ierr; 405b4319ba4SBarry Smith 406b4319ba4SBarry Smith PetscFunctionBegin; 407b4319ba4SBarry Smith pcis->is_B_local = 0; 408b4319ba4SBarry Smith pcis->is_I_local = 0; 409b4319ba4SBarry Smith pcis->is_B_global = 0; 410b4319ba4SBarry Smith pcis->is_I_global = 0; 411b4319ba4SBarry Smith pcis->A_II = 0; 412b4319ba4SBarry Smith pcis->A_IB = 0; 413b4319ba4SBarry Smith pcis->A_BI = 0; 414b4319ba4SBarry Smith pcis->A_BB = 0; 415b4319ba4SBarry Smith pcis->D = 0; 416b4319ba4SBarry Smith pcis->ksp_N = 0; 417b4319ba4SBarry Smith pcis->ksp_D = 0; 418b4319ba4SBarry Smith pcis->vec1_N = 0; 419b4319ba4SBarry Smith pcis->vec2_N = 0; 420b4319ba4SBarry Smith pcis->vec1_D = 0; 421b4319ba4SBarry Smith pcis->vec2_D = 0; 422b4319ba4SBarry Smith pcis->vec3_D = 0; 423b4319ba4SBarry Smith pcis->vec1_B = 0; 424b4319ba4SBarry Smith pcis->vec2_B = 0; 425b4319ba4SBarry Smith pcis->vec3_B = 0; 426b4319ba4SBarry Smith pcis->vec1_global = 0; 427b4319ba4SBarry Smith pcis->work_N = 0; 428b4319ba4SBarry Smith pcis->global_to_D = 0; 429b4319ba4SBarry Smith pcis->N_to_B = 0; 430b4319ba4SBarry Smith pcis->global_to_B = 0; 4312fa5cd67SKarl Rupp 432b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE; 4332fa5cd67SKarl Rupp 434831a100dSStefano Zampini pcis->scaling_factor = 1.0; 435831a100dSStefano Zampini /* composing functions */ 436b965d45eSStefano Zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetUseStiffnessScaling_C","PCISSetUseStiffnessScaling_IS", 437b965d45eSStefano Zampini PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr); 438831a100dSStefano Zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainScalingFactor_C","PCISSetSubdomainScalingFactor_IS", 439831a100dSStefano Zampini PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 440ba1573a8SStefano Zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C","PCISSetSubdomainDiagonalScaling_IS", 441ba1573a8SStefano Zampini PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 442b4319ba4SBarry Smith PetscFunctionReturn(0); 443b4319ba4SBarry Smith } 444b4319ba4SBarry Smith 445b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 446b4319ba4SBarry Smith /* 447b4319ba4SBarry Smith PCISApplySchur - 448b4319ba4SBarry Smith 449b4319ba4SBarry Smith Input parameters: 450b4319ba4SBarry Smith . pc - preconditioner context 451b4319ba4SBarry Smith . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 452b4319ba4SBarry Smith 453b4319ba4SBarry Smith Output parameters: 454b4319ba4SBarry Smith . vec1_B - result of Schur complement applied to chunk 455b4319ba4SBarry Smith . vec2_B - garbage (used as work space), or null (and v is used as workspace) 456b4319ba4SBarry Smith . vec1_D - garbage (used as work space) 457b4319ba4SBarry Smith . vec2_D - garbage (used as work space) 458b4319ba4SBarry Smith 459b4319ba4SBarry Smith */ 460b4319ba4SBarry Smith #undef __FUNCT__ 461831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur" 4627087cfbeSBarry Smith PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 463b4319ba4SBarry Smith { 464dfbe8321SBarry Smith PetscErrorCode ierr; 465b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 466b4319ba4SBarry Smith 467b4319ba4SBarry Smith PetscFunctionBegin; 4682fa5cd67SKarl Rupp if (!vec2_B) vec2_B = v; 469b4319ba4SBarry Smith 470b4319ba4SBarry Smith ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 471b4319ba4SBarry Smith ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 47223ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 473b4319ba4SBarry Smith ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 474efb30889SBarry Smith ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 475b4319ba4SBarry Smith PetscFunctionReturn(0); 476b4319ba4SBarry Smith } 477b4319ba4SBarry Smith 478b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 479b4319ba4SBarry Smith /* 480b4319ba4SBarry Smith PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 481b4319ba4SBarry Smith including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 482b4319ba4SBarry Smith mode. 483b4319ba4SBarry Smith 484b4319ba4SBarry Smith Input parameters: 485b4319ba4SBarry Smith . pc - preconditioner context 486b4319ba4SBarry Smith . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 487b4319ba4SBarry Smith . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 488b4319ba4SBarry Smith 489b4319ba4SBarry Smith Output parameter: 490b4319ba4SBarry Smith . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 491b4319ba4SBarry Smith . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 492b4319ba4SBarry Smith 493b4319ba4SBarry Smith Notes: 494b4319ba4SBarry Smith The entries in the array that do not correspond to interface nodes remain unaltered. 495b4319ba4SBarry Smith */ 496b4319ba4SBarry Smith #undef __FUNCT__ 497b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB" 4987087cfbeSBarry Smith PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 499b4319ba4SBarry Smith { 5005d0c19d7SBarry Smith PetscInt i; 5015d0c19d7SBarry Smith const PetscInt *idex; 5026849ba73SBarry Smith PetscErrorCode ierr; 503b4319ba4SBarry Smith PetscScalar *array_B; 504b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 505b4319ba4SBarry Smith 506b4319ba4SBarry Smith PetscFunctionBegin; 507b4319ba4SBarry Smith ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 508b4319ba4SBarry Smith ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 509b4319ba4SBarry Smith 510b4319ba4SBarry Smith if (smode == SCATTER_FORWARD) { 511b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 5122fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]]; 513b4319ba4SBarry Smith } else { /* ADD_VALUES */ 5142fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]]; 515b4319ba4SBarry Smith } 516b4319ba4SBarry Smith } else { /* SCATTER_REVERSE */ 517b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 5182fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i]; 519b4319ba4SBarry Smith } else { /* ADD_VALUES */ 5202fa5cd67SKarl Rupp for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i]; 521b4319ba4SBarry Smith } 522b4319ba4SBarry Smith } 523b4319ba4SBarry Smith ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 524b4319ba4SBarry Smith ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 525b4319ba4SBarry Smith PetscFunctionReturn(0); 526b4319ba4SBarry Smith } 527b4319ba4SBarry Smith 528b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 529b4319ba4SBarry Smith /* 530b4319ba4SBarry Smith PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 531b4319ba4SBarry Smith More precisely, solves the problem: 532b4319ba4SBarry Smith [ A_II A_IB ] [ . ] [ 0 ] 533b4319ba4SBarry Smith [ ] [ ] = [ ] 534b4319ba4SBarry Smith [ A_BI A_BB ] [ x ] [ b ] 535b4319ba4SBarry Smith 536b4319ba4SBarry Smith Input parameters: 537b4319ba4SBarry Smith . pc - preconditioner context 538b4319ba4SBarry Smith . b - vector of local interface nodes (including ghosts) 539b4319ba4SBarry Smith 540b4319ba4SBarry Smith Output parameters: 541b4319ba4SBarry Smith . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 542b4319ba4SBarry Smith complement to b 543b4319ba4SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 544b4319ba4SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 545b4319ba4SBarry Smith 546b4319ba4SBarry Smith */ 547b4319ba4SBarry Smith #undef __FUNCT__ 548b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur" 5497087cfbeSBarry Smith PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 550b4319ba4SBarry Smith { 551dfbe8321SBarry Smith PetscErrorCode ierr; 552b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 553b4319ba4SBarry Smith 554b4319ba4SBarry Smith PetscFunctionBegin; 555b4319ba4SBarry Smith /* 556b4319ba4SBarry Smith Neumann solvers. 557b4319ba4SBarry Smith Applying the inverse of the local Schur complement, i.e, solving a Neumann 558b4319ba4SBarry Smith Problem with zero at the interior nodes of the RHS and extracting the interface 559b4319ba4SBarry Smith part of the solution. inverse Schur complement is applied to b and the result 560b4319ba4SBarry Smith is stored in x. 561b4319ba4SBarry Smith */ 562b4319ba4SBarry Smith /* Setting the RHS vec1_N */ 563efb30889SBarry Smith ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 564ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 565ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 566b4319ba4SBarry Smith /* Checking for consistency of the RHS */ 567b4319ba4SBarry Smith { 568ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 5690298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr); 570b4319ba4SBarry Smith if (flg) { 571b4319ba4SBarry Smith PetscScalar average; 5723050cee2SBarry Smith PetscViewer viewer; 573*ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 5743050cee2SBarry Smith 575b4319ba4SBarry Smith ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 576b4319ba4SBarry Smith average = average / ((PetscReal)pcis->n); 5777b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 578b4319ba4SBarry Smith if (pcis->pure_neumann) { 5799f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 580b4319ba4SBarry Smith } else { 5819f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 582b4319ba4SBarry Smith } 5833050cee2SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5847b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 585b4319ba4SBarry Smith } 586b4319ba4SBarry Smith } 587b4319ba4SBarry Smith /* Solving the system for vec2_N */ 58823ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 589b4319ba4SBarry Smith /* Extracting the local interface vector out of the solution */ 590ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 591ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 592b4319ba4SBarry Smith PetscFunctionReturn(0); 593b4319ba4SBarry Smith } 594b4319ba4SBarry Smith 595b4319ba4SBarry Smith 596b4319ba4SBarry Smith 597b4319ba4SBarry Smith 598b4319ba4SBarry Smith 599b4319ba4SBarry Smith 600b4319ba4SBarry Smith 601b4319ba4SBarry Smith 602b4319ba4SBarry Smith 603