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__ 6ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling_IS" 7ba1573a8SStefano Zampini static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors) 8ba1573a8SStefano Zampini { 9ba1573a8SStefano Zampini PetscErrorCode ierr; 10ba1573a8SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 11ba1573a8SStefano Zampini 12ba1573a8SStefano Zampini PetscFunctionBegin; 13ba1573a8SStefano Zampini ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 14ba1573a8SStefano Zampini ierr = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr); 15ba1573a8SStefano Zampini pcis->D = scaling_factors; 16ba1573a8SStefano Zampini PetscFunctionReturn(0); 17ba1573a8SStefano Zampini } 18ba1573a8SStefano Zampini EXTERN_C_END 19ba1573a8SStefano Zampini 20ba1573a8SStefano Zampini #undef __FUNCT__ 21ba1573a8SStefano Zampini #define __FUNCT__ "PCISSetSubdomainDiagonalScaling" 22ba1573a8SStefano Zampini /*@ 23ba1573a8SStefano Zampini PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS. 24ba1573a8SStefano Zampini 25ba1573a8SStefano Zampini Not collective 26ba1573a8SStefano Zampini 27ba1573a8SStefano Zampini Input Parameters: 28ba1573a8SStefano Zampini + pc - the preconditioning context 29ba1573a8SStefano Zampini - scaling_factors - scaling factors for the subdomain 30ba1573a8SStefano Zampini 31ba1573a8SStefano Zampini Level: intermediate 32ba1573a8SStefano Zampini 33ba1573a8SStefano Zampini Notes: 34ba1573a8SStefano Zampini Intended to use with jumping coefficients cases. 35ba1573a8SStefano Zampini 36ba1573a8SStefano Zampini .seealso: PCBDDC 37ba1573a8SStefano Zampini @*/ 38ba1573a8SStefano Zampini PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors) 39ba1573a8SStefano Zampini { 40ba1573a8SStefano Zampini PetscErrorCode ierr; 41ba1573a8SStefano Zampini 42ba1573a8SStefano Zampini PetscFunctionBegin; 43ba1573a8SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 44ba1573a8SStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr); 45ba1573a8SStefano Zampini PetscFunctionReturn(0); 46ba1573a8SStefano Zampini } 47ba1573a8SStefano Zampini 48ba1573a8SStefano Zampini EXTERN_C_BEGIN 49ba1573a8SStefano Zampini #undef __FUNCT__ 50831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS" 51831a100dSStefano Zampini static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal) 52831a100dSStefano Zampini { 53831a100dSStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 54831a100dSStefano Zampini 55831a100dSStefano Zampini PetscFunctionBegin; 56831a100dSStefano Zampini pcis->scaling_factor = scal; 57831a100dSStefano Zampini PetscFunctionReturn(0); 58831a100dSStefano Zampini } 59831a100dSStefano Zampini EXTERN_C_END 60831a100dSStefano Zampini 61831a100dSStefano Zampini #undef __FUNCT__ 62831a100dSStefano Zampini #define __FUNCT__ "PCISSetSubdomainScalingFactor" 63831a100dSStefano Zampini /*@ 64831a100dSStefano Zampini PCISSetSubdomainScalingFactor - Set scaling factor for PCIS. 65831a100dSStefano Zampini 66831a100dSStefano Zampini Not collective 67831a100dSStefano Zampini 68831a100dSStefano Zampini Input Parameters: 69831a100dSStefano Zampini + pc - the preconditioning context 70831a100dSStefano Zampini - scal - scaling factor for the subdomain 71831a100dSStefano Zampini 72831a100dSStefano Zampini Level: intermediate 73831a100dSStefano Zampini 74831a100dSStefano Zampini Notes: 75831a100dSStefano Zampini Intended to use with jumping coefficients cases. 76831a100dSStefano Zampini 77831a100dSStefano Zampini .seealso: PCBDDC 78831a100dSStefano Zampini @*/ 79831a100dSStefano Zampini PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal) 80831a100dSStefano Zampini { 81831a100dSStefano Zampini PetscErrorCode ierr; 82831a100dSStefano Zampini 83831a100dSStefano Zampini PetscFunctionBegin; 84831a100dSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 85831a100dSStefano Zampini ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr); 86831a100dSStefano Zampini PetscFunctionReturn(0); 87831a100dSStefano Zampini } 88831a100dSStefano Zampini 89b4319ba4SBarry Smith 90b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 91b4319ba4SBarry Smith /* 92b4319ba4SBarry Smith PCISSetUp - 93b4319ba4SBarry Smith */ 94b4319ba4SBarry Smith #undef __FUNCT__ 95b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp" 967087cfbeSBarry Smith PetscErrorCode PCISSetUp(PC pc) 97b4319ba4SBarry Smith { 98b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 99b4319ba4SBarry Smith Mat_IS *matis = (Mat_IS*)pc->mat->data; 10013f74950SBarry Smith PetscInt i; 1016849ba73SBarry Smith PetscErrorCode ierr; 102ace3abfcSBarry Smith PetscBool flg; 103831a100dSStefano Zampini Vec counter; 104b4319ba4SBarry Smith 105b4319ba4SBarry Smith PetscFunctionBegin; 106251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr); 107e7e72b3dSBarry Smith if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 108b4319ba4SBarry Smith 109b4319ba4SBarry Smith pcis->pure_neumann = matis->pure_neumann; 110b4319ba4SBarry Smith 111b4319ba4SBarry Smith /* 112b4319ba4SBarry Smith Creating the local vector vec1_N, containing the inverse of the number 113b4319ba4SBarry Smith of subdomains to which each local node (either owned or ghost) 114b4319ba4SBarry Smith pertains. To accomplish that, we scatter local vectors of 1's to 115b4319ba4SBarry Smith a global vector (adding the values); scatter the result back to 116b4319ba4SBarry Smith local vectors and finally invert the result. 117b4319ba4SBarry Smith */ 118b4319ba4SBarry Smith ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 11923ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ 120efb30889SBarry Smith ierr = VecSet(counter,0.0);CHKERRQ(ierr); 121efb30889SBarry Smith ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr); 122ca9f406cSSatish Balay ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 123ca9f406cSSatish Balay ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 124ca9f406cSSatish Balay ierr = VecScatterBegin(matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125ca9f406cSSatish Balay ierr = VecScatterEnd (matis->ctx,counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126b4319ba4SBarry Smith /* 127b4319ba4SBarry Smith Creating local and global index sets for interior and 128b4319ba4SBarry Smith inteface nodes. Notice that interior nodes have D[i]==1.0. 129b4319ba4SBarry Smith */ 130b4319ba4SBarry Smith { 13113f74950SBarry Smith PetscInt n_I; 13213f74950SBarry Smith PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 133b4319ba4SBarry Smith PetscScalar *array; 134b4319ba4SBarry Smith /* Identifying interior and interface nodes, in local numbering */ 135b4319ba4SBarry Smith ierr = VecGetSize(pcis->vec1_N,&pcis->n);CHKERRQ(ierr); 136b4319ba4SBarry Smith ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 13713f74950SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);CHKERRQ(ierr); 13813f74950SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);CHKERRQ(ierr); 139b4319ba4SBarry Smith for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 140b4319ba4SBarry Smith if (array[i] == 1.0) { idx_I_local[n_I] = i; n_I++; } 141b4319ba4SBarry Smith else { idx_B_local[pcis->n_B] = i; pcis->n_B++; } 142b4319ba4SBarry Smith } 143b4319ba4SBarry Smith /* Getting the global numbering */ 144b4319ba4SBarry Smith idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 145b4319ba4SBarry Smith idx_I_global = idx_B_local + pcis->n_B; 146b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 147b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I, idx_I_local,idx_I_global);CHKERRQ(ierr); 148b4319ba4SBarry Smith /* Creating the index sets. */ 14970b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); 15070b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); 15170b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); 15270b3c8c7SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); 153b4319ba4SBarry Smith /* Freeing memory and restoring arrays */ 154b4319ba4SBarry Smith ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 155b4319ba4SBarry Smith ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 156b4319ba4SBarry Smith ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 157b4319ba4SBarry Smith } 158b4319ba4SBarry Smith 159b4319ba4SBarry Smith /* 160b4319ba4SBarry Smith Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 161b4319ba4SBarry Smith is such that interior nodes come first than the interface ones, we have 162b4319ba4SBarry Smith 163b4319ba4SBarry Smith [ | ] 164b4319ba4SBarry Smith [ A_II | A_IB ] 165b4319ba4SBarry Smith A = [ | ] 166b4319ba4SBarry Smith [-----------+------] 167b4319ba4SBarry Smith [ A_BI | A_BB ] 168b4319ba4SBarry Smith */ 169b4319ba4SBarry Smith 1704aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr); 1714aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 1724aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 1734aa3045dSJed Brown ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 174b4319ba4SBarry Smith 175b4319ba4SBarry Smith /* 176b4319ba4SBarry Smith Creating work vectors and arrays 177b4319ba4SBarry Smith */ 178b4319ba4SBarry Smith /* pcis->vec1_N has already been created */ 179b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 180b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 181b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 182b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 183b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 184b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 185b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 18623ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 187b4319ba4SBarry Smith ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr); 188b4319ba4SBarry Smith 189b4319ba4SBarry Smith /* Creating the scatter contexts */ 19023ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 191b4319ba4SBarry Smith ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 19223ce1328SBarry Smith ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 193b4319ba4SBarry Smith 194831a100dSStefano Zampini /* Creating scaling "matrix" D */ 195831a100dSStefano Zampini if( !pcis->D ) { 196ba1573a8SStefano Zampini ierr = VecSet(pcis->vec1_B,pcis->scaling_factor);CHKERRQ(ierr); 197ba1573a8SStefano Zampini } else { 198ba1573a8SStefano Zampini ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 199ba1573a8SStefano Zampini } 200831a100dSStefano Zampini ierr = VecSet(counter,0.0);CHKERRQ(ierr); 201ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 202ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 203ba1573a8SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 204ba1573a8SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 205ba1573a8SStefano Zampini if( !pcis->D ) { 206b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 207ba1573a8SStefano Zampini ierr = VecCopy(pcis->vec1_B,pcis->D);CHKERRQ(ierr); 208b4319ba4SBarry Smith ierr = VecReciprocal(pcis->D);CHKERRQ(ierr); 209831a100dSStefano Zampini ierr = VecScale(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 210ba1573a8SStefano Zampini } else { 211ba1573a8SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 212831a100dSStefano Zampini } 213b4319ba4SBarry Smith 214b4319ba4SBarry Smith /* See historical note 01, at the bottom of this file. */ 215b4319ba4SBarry Smith 216b4319ba4SBarry Smith /* 217b4319ba4SBarry Smith Creating the KSP contexts for the local Dirichlet and Neumann problems. 218b4319ba4SBarry Smith */ 219b4319ba4SBarry Smith { 220b4319ba4SBarry Smith PC pc_ctx; 221b4319ba4SBarry Smith /* Dirichlet */ 222b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 2231cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 224b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr); 225b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 226b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 227b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 228b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 229b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 230b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 231b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 232b4319ba4SBarry Smith /* Neumann */ 233b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 2341cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); 235b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr); 236b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 237b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 238b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 239b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 240b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 241b4319ba4SBarry Smith { 242ace3abfcSBarry Smith PetscBool damp_fixed = PETSC_FALSE, 24390d69ab7SBarry Smith remove_nullspace_fixed = PETSC_FALSE, 24490d69ab7SBarry Smith set_damping_factor_floating = PETSC_FALSE, 24590d69ab7SBarry Smith not_damp_floating = PETSC_FALSE, 24690d69ab7SBarry Smith not_remove_nullspace_floating = PETSC_FALSE; 247b4319ba4SBarry Smith PetscReal fixed_factor, 248b4319ba4SBarry Smith floating_factor; 249b4319ba4SBarry Smith 2507adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 251b4319ba4SBarry Smith if (!damp_fixed) { fixed_factor = 0.0; } 252acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,PETSC_NULL);CHKERRQ(ierr); 253b4319ba4SBarry Smith 254acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,PETSC_NULL);CHKERRQ(ierr); 255b4319ba4SBarry Smith 2567adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", 257b4319ba4SBarry Smith &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 258b4319ba4SBarry Smith if (!set_damping_factor_floating) { floating_factor = 0.0; } 259acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,PETSC_NULL);CHKERRQ(ierr); 260b4319ba4SBarry Smith if (!set_damping_factor_floating) { floating_factor = 1.e-12; } 261b4319ba4SBarry Smith 262acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,PETSC_NULL);CHKERRQ(ierr); 263b4319ba4SBarry Smith 264acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,PETSC_NULL);CHKERRQ(ierr); 265b4319ba4SBarry Smith 266b4319ba4SBarry Smith if (pcis->pure_neumann) { /* floating subdomain */ 267b4319ba4SBarry Smith if (!(not_damp_floating)) { 268c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 269c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 270b4319ba4SBarry Smith } 271b4319ba4SBarry Smith if (!(not_remove_nullspace_floating)){ 272b4319ba4SBarry Smith MatNullSpace nullsp; 273ee512c37SSatish Balay ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr); 274d8fd42c4SBarry Smith ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr); 275d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 276b4319ba4SBarry Smith } 277b4319ba4SBarry Smith } else { /* fixed subdomain */ 278b4319ba4SBarry Smith if (damp_fixed) { 279c88783f4SHong Zhang ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 280c88783f4SHong Zhang ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 281b4319ba4SBarry Smith } 282b4319ba4SBarry Smith if (remove_nullspace_fixed) { 283b4319ba4SBarry Smith MatNullSpace nullsp; 284ee512c37SSatish Balay ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,PETSC_NULL,&nullsp);CHKERRQ(ierr); 285d8fd42c4SBarry Smith ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr); 286d34fcf5fSBarry Smith ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 287b4319ba4SBarry Smith } 288b4319ba4SBarry Smith } 289b4319ba4SBarry Smith } 290b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 291b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 292b4319ba4SBarry Smith } 293b4319ba4SBarry Smith 2947c334f02SBarry Smith ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 295b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE; 296831a100dSStefano Zampini ierr = VecDestroy(&counter);CHKERRQ(ierr); 297b4319ba4SBarry Smith PetscFunctionReturn(0); 298b4319ba4SBarry Smith } 299b4319ba4SBarry Smith 300b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 301b4319ba4SBarry Smith /* 302b4319ba4SBarry Smith PCISDestroy - 303b4319ba4SBarry Smith */ 304b4319ba4SBarry Smith #undef __FUNCT__ 305b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy" 3067087cfbeSBarry Smith PetscErrorCode PCISDestroy(PC pc) 307b4319ba4SBarry Smith { 308b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 309dfbe8321SBarry Smith PetscErrorCode ierr; 310b4319ba4SBarry Smith 311b4319ba4SBarry Smith PetscFunctionBegin; 312fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 313fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 314fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 315fcfd50ebSBarry Smith ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 316fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 317fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 318fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 319fcfd50ebSBarry Smith ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 320fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 321fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 322fcfd50ebSBarry Smith ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 323fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 324fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 325fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 326fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 327fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 328fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 329fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 330fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 331fcfd50ebSBarry Smith ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 332fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 333fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 334fcfd50ebSBarry Smith ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 33505b42c5fSBarry Smith ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 336b4319ba4SBarry Smith if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) { 337b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 338b4319ba4SBarry Smith } 339*46caae47SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainScalingFactor_C","",PETSC_NULL);CHKERRQ(ierr); 340*46caae47SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C","",PETSC_NULL);CHKERRQ(ierr); 341b4319ba4SBarry Smith PetscFunctionReturn(0); 342b4319ba4SBarry Smith } 343b4319ba4SBarry Smith 344b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 345b4319ba4SBarry Smith /* 346b4319ba4SBarry Smith PCISCreate - 347b4319ba4SBarry Smith */ 348b4319ba4SBarry Smith #undef __FUNCT__ 349b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate" 3507087cfbeSBarry Smith PetscErrorCode PCISCreate(PC pc) 351b4319ba4SBarry Smith { 352b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 353831a100dSStefano Zampini PetscErrorCode ierr; 354b4319ba4SBarry Smith 355b4319ba4SBarry Smith PetscFunctionBegin; 356b4319ba4SBarry Smith pcis->is_B_local = 0; 357b4319ba4SBarry Smith pcis->is_I_local = 0; 358b4319ba4SBarry Smith pcis->is_B_global = 0; 359b4319ba4SBarry Smith pcis->is_I_global = 0; 360b4319ba4SBarry Smith pcis->A_II = 0; 361b4319ba4SBarry Smith pcis->A_IB = 0; 362b4319ba4SBarry Smith pcis->A_BI = 0; 363b4319ba4SBarry Smith pcis->A_BB = 0; 364b4319ba4SBarry Smith pcis->D = 0; 365b4319ba4SBarry Smith pcis->ksp_N = 0; 366b4319ba4SBarry Smith pcis->ksp_D = 0; 367b4319ba4SBarry Smith pcis->vec1_N = 0; 368b4319ba4SBarry Smith pcis->vec2_N = 0; 369b4319ba4SBarry Smith pcis->vec1_D = 0; 370b4319ba4SBarry Smith pcis->vec2_D = 0; 371b4319ba4SBarry Smith pcis->vec3_D = 0; 372b4319ba4SBarry Smith pcis->vec1_B = 0; 373b4319ba4SBarry Smith pcis->vec2_B = 0; 374b4319ba4SBarry Smith pcis->vec3_B = 0; 375b4319ba4SBarry Smith pcis->vec1_global = 0; 376b4319ba4SBarry Smith pcis->work_N = 0; 377b4319ba4SBarry Smith pcis->global_to_D = 0; 378b4319ba4SBarry Smith pcis->N_to_B = 0; 379b4319ba4SBarry Smith pcis->global_to_B = 0; 380b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE; 381831a100dSStefano Zampini pcis->scaling_factor = 1.0; 382831a100dSStefano Zampini /* composing functions */ 383831a100dSStefano Zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainScalingFactor_C","PCISSetSubdomainScalingFactor_IS", 384831a100dSStefano Zampini PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 385ba1573a8SStefano Zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C","PCISSetSubdomainDiagonalScaling_IS", 386ba1573a8SStefano Zampini PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 387b4319ba4SBarry Smith PetscFunctionReturn(0); 388b4319ba4SBarry Smith } 389b4319ba4SBarry Smith 390b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 391b4319ba4SBarry Smith /* 392b4319ba4SBarry Smith PCISApplySchur - 393b4319ba4SBarry Smith 394b4319ba4SBarry Smith Input parameters: 395b4319ba4SBarry Smith . pc - preconditioner context 396b4319ba4SBarry Smith . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 397b4319ba4SBarry Smith 398b4319ba4SBarry Smith Output parameters: 399b4319ba4SBarry Smith . vec1_B - result of Schur complement applied to chunk 400b4319ba4SBarry Smith . vec2_B - garbage (used as work space), or null (and v is used as workspace) 401b4319ba4SBarry Smith . vec1_D - garbage (used as work space) 402b4319ba4SBarry Smith . vec2_D - garbage (used as work space) 403b4319ba4SBarry Smith 404b4319ba4SBarry Smith */ 405b4319ba4SBarry Smith #undef __FUNCT__ 406831a100dSStefano Zampini #define __FUNCT__ "PCISApplySchur" 4077087cfbeSBarry Smith PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 408b4319ba4SBarry Smith { 409dfbe8321SBarry Smith PetscErrorCode ierr; 410b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 411b4319ba4SBarry Smith 412b4319ba4SBarry Smith PetscFunctionBegin; 41313f74950SBarry Smith if (!vec2_B) { vec2_B = v; } 414b4319ba4SBarry Smith 415b4319ba4SBarry Smith ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 416b4319ba4SBarry Smith ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 41723ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 418b4319ba4SBarry Smith ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 419efb30889SBarry Smith ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 420b4319ba4SBarry Smith PetscFunctionReturn(0); 421b4319ba4SBarry Smith } 422b4319ba4SBarry Smith 423b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 424b4319ba4SBarry Smith /* 425b4319ba4SBarry Smith PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 426b4319ba4SBarry Smith including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 427b4319ba4SBarry Smith mode. 428b4319ba4SBarry Smith 429b4319ba4SBarry Smith Input parameters: 430b4319ba4SBarry Smith . pc - preconditioner context 431b4319ba4SBarry Smith . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 432b4319ba4SBarry Smith . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 433b4319ba4SBarry Smith 434b4319ba4SBarry Smith Output parameter: 435b4319ba4SBarry Smith . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 436b4319ba4SBarry Smith . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 437b4319ba4SBarry Smith 438b4319ba4SBarry Smith Notes: 439b4319ba4SBarry Smith The entries in the array that do not correspond to interface nodes remain unaltered. 440b4319ba4SBarry Smith */ 441b4319ba4SBarry Smith #undef __FUNCT__ 442b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB" 4437087cfbeSBarry Smith PetscErrorCode PCISScatterArrayNToVecB (PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 444b4319ba4SBarry Smith { 4455d0c19d7SBarry Smith PetscInt i; 4465d0c19d7SBarry Smith const PetscInt *idex; 4476849ba73SBarry Smith PetscErrorCode ierr; 448b4319ba4SBarry Smith PetscScalar *array_B; 449b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 450b4319ba4SBarry Smith 451b4319ba4SBarry Smith PetscFunctionBegin; 452b4319ba4SBarry Smith ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 453b4319ba4SBarry Smith ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 454b4319ba4SBarry Smith 455b4319ba4SBarry Smith if (smode == SCATTER_FORWARD) { 456b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 457b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_B[i] = array_N[idex[i]]; } 458b4319ba4SBarry Smith } else { /* ADD_VALUES */ 459b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_B[i] += array_N[idex[i]]; } 460b4319ba4SBarry Smith } 461b4319ba4SBarry Smith } else { /* SCATTER_REVERSE */ 462b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 463b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] = array_B[i]; } 464b4319ba4SBarry Smith } else { /* ADD_VALUES */ 465b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] += array_B[i]; } 466b4319ba4SBarry Smith } 467b4319ba4SBarry Smith } 468b4319ba4SBarry Smith ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 469b4319ba4SBarry Smith ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 470b4319ba4SBarry Smith PetscFunctionReturn(0); 471b4319ba4SBarry Smith } 472b4319ba4SBarry Smith 473b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 474b4319ba4SBarry Smith /* 475b4319ba4SBarry Smith PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 476b4319ba4SBarry Smith More precisely, solves the problem: 477b4319ba4SBarry Smith [ A_II A_IB ] [ . ] [ 0 ] 478b4319ba4SBarry Smith [ ] [ ] = [ ] 479b4319ba4SBarry Smith [ A_BI A_BB ] [ x ] [ b ] 480b4319ba4SBarry Smith 481b4319ba4SBarry Smith Input parameters: 482b4319ba4SBarry Smith . pc - preconditioner context 483b4319ba4SBarry Smith . b - vector of local interface nodes (including ghosts) 484b4319ba4SBarry Smith 485b4319ba4SBarry Smith Output parameters: 486b4319ba4SBarry Smith . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 487b4319ba4SBarry Smith complement to b 488b4319ba4SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 489b4319ba4SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 490b4319ba4SBarry Smith 491b4319ba4SBarry Smith */ 492b4319ba4SBarry Smith #undef __FUNCT__ 493b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur" 4947087cfbeSBarry Smith PetscErrorCode PCISApplyInvSchur (PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 495b4319ba4SBarry Smith { 496dfbe8321SBarry Smith PetscErrorCode ierr; 497b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 498b4319ba4SBarry Smith 499b4319ba4SBarry Smith PetscFunctionBegin; 500b4319ba4SBarry Smith /* 501b4319ba4SBarry Smith Neumann solvers. 502b4319ba4SBarry Smith Applying the inverse of the local Schur complement, i.e, solving a Neumann 503b4319ba4SBarry Smith Problem with zero at the interior nodes of the RHS and extracting the interface 504b4319ba4SBarry Smith part of the solution. inverse Schur complement is applied to b and the result 505b4319ba4SBarry Smith is stored in x. 506b4319ba4SBarry Smith */ 507b4319ba4SBarry Smith /* Setting the RHS vec1_N */ 508efb30889SBarry Smith ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 509ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 510ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 511b4319ba4SBarry Smith /* Checking for consistency of the RHS */ 512b4319ba4SBarry Smith { 513ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 514acfcf0e5SJed Brown ierr = PetscOptionsGetBool(PETSC_NULL,"-pc_is_check_consistency",&flg,PETSC_NULL);CHKERRQ(ierr); 515b4319ba4SBarry Smith if (flg) { 516b4319ba4SBarry Smith PetscScalar average; 5173050cee2SBarry Smith PetscViewer viewer; 5187adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);CHKERRQ(ierr); 5193050cee2SBarry Smith 520b4319ba4SBarry Smith ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 521b4319ba4SBarry Smith average = average / ((PetscReal)pcis->n); 5227b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 523b4319ba4SBarry Smith if (pcis->pure_neumann) { 5249f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 525b4319ba4SBarry Smith } else { 5269f73f8ecSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 527b4319ba4SBarry Smith } 5283050cee2SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5297b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 530b4319ba4SBarry Smith } 531b4319ba4SBarry Smith } 532b4319ba4SBarry Smith /* Solving the system for vec2_N */ 53323ce1328SBarry Smith ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 534b4319ba4SBarry Smith /* Extracting the local interface vector out of the solution */ 535ca9f406cSSatish Balay ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 536ca9f406cSSatish Balay ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 537b4319ba4SBarry Smith PetscFunctionReturn(0); 538b4319ba4SBarry Smith } 539b4319ba4SBarry Smith 540b4319ba4SBarry Smith 541b4319ba4SBarry Smith 542b4319ba4SBarry Smith 543b4319ba4SBarry Smith 544b4319ba4SBarry Smith 545b4319ba4SBarry Smith 546b4319ba4SBarry Smith 547b4319ba4SBarry Smith 548