1*b4319ba4SBarry Smith /*$Id: is.c,v 1.9 2001/08/07 03:03:41 balay Exp $*/ 2*b4319ba4SBarry Smith #include "src/ksp/pc/impls/is/pcis.h" 3*b4319ba4SBarry Smith 4*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 5*b4319ba4SBarry Smith /* 6*b4319ba4SBarry Smith PCISSetUp - 7*b4319ba4SBarry Smith */ 8*b4319ba4SBarry Smith #undef __FUNCT__ 9*b4319ba4SBarry Smith #define __FUNCT__ "PCISSetUp" 10*b4319ba4SBarry Smith int PCISSetUp(PC pc) 11*b4319ba4SBarry Smith { 12*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 13*b4319ba4SBarry Smith Mat_IS *matis = (Mat_IS*)pc->mat->data; 14*b4319ba4SBarry Smith int i, ierr; 15*b4319ba4SBarry Smith PetscTruth flg; 16*b4319ba4SBarry Smith 17*b4319ba4SBarry Smith PetscFunctionBegin; 18*b4319ba4SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr); 19*b4319ba4SBarry Smith if (!flg){ 20*b4319ba4SBarry Smith SETERRQ(1,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 21*b4319ba4SBarry Smith } 22*b4319ba4SBarry Smith 23*b4319ba4SBarry Smith pcis->pure_neumann = matis->pure_neumann; 24*b4319ba4SBarry Smith 25*b4319ba4SBarry Smith /* 26*b4319ba4SBarry Smith Creating the local vector vec1_N, containing the inverse of the number 27*b4319ba4SBarry Smith of subdomains to which each local node (either owned or ghost) 28*b4319ba4SBarry Smith pertains. To accomplish that, we scatter local vectors of 1's to 29*b4319ba4SBarry Smith a global vector (adding the values); scatter the result back to 30*b4319ba4SBarry Smith local vectors and finally invert the result. 31*b4319ba4SBarry Smith */ 32*b4319ba4SBarry Smith { 33*b4319ba4SBarry Smith Vec counter; 34*b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 35*b4319ba4SBarry Smith ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 36*b4319ba4SBarry Smith ierr = VecDuplicate(pc->vec,&counter);CHKERRQ(ierr); /* temporary auxiliar vector */ 37*b4319ba4SBarry Smith ierr = VecSet(&zero,counter);CHKERRQ(ierr); 38*b4319ba4SBarry Smith ierr = VecSet(&one,pcis->vec1_N);CHKERRQ(ierr); 39*b4319ba4SBarry Smith ierr = VecScatterBegin(pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);CHKERRQ(ierr); 40*b4319ba4SBarry Smith ierr = VecScatterEnd (pcis->vec1_N,counter,ADD_VALUES,SCATTER_REVERSE,matis->ctx);CHKERRQ(ierr); 41*b4319ba4SBarry Smith ierr = VecScatterBegin(counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);CHKERRQ(ierr); 42*b4319ba4SBarry Smith ierr = VecScatterEnd (counter,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD,matis->ctx);CHKERRQ(ierr); 43*b4319ba4SBarry Smith ierr = VecDestroy(counter);CHKERRQ(ierr); 44*b4319ba4SBarry Smith } 45*b4319ba4SBarry Smith /* 46*b4319ba4SBarry Smith Creating local and global index sets for interior and 47*b4319ba4SBarry Smith inteface nodes. Notice that interior nodes have D[i]==1.0. 48*b4319ba4SBarry Smith */ 49*b4319ba4SBarry Smith { 50*b4319ba4SBarry Smith int n_I; 51*b4319ba4SBarry Smith int *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 52*b4319ba4SBarry Smith PetscScalar *array; 53*b4319ba4SBarry Smith /* Identifying interior and interface nodes, in local numbering */ 54*b4319ba4SBarry Smith ierr = VecGetSize(pcis->vec1_N,&pcis->n);CHKERRQ(ierr); 55*b4319ba4SBarry Smith ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 56*b4319ba4SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(int),&idx_I_local);CHKERRQ(ierr); 57*b4319ba4SBarry Smith ierr = PetscMalloc(pcis->n*sizeof(int),&idx_B_local);CHKERRQ(ierr); 58*b4319ba4SBarry Smith for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 59*b4319ba4SBarry Smith if (array[i] == 1.0) { idx_I_local[n_I] = i; n_I++; } 60*b4319ba4SBarry Smith else { idx_B_local[pcis->n_B] = i; pcis->n_B++; } 61*b4319ba4SBarry Smith } 62*b4319ba4SBarry Smith /* Getting the global numbering */ 63*b4319ba4SBarry Smith idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 64*b4319ba4SBarry Smith idx_I_global = idx_B_local + pcis->n_B; 65*b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 66*b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I, idx_I_local,idx_I_global);CHKERRQ(ierr); 67*b4319ba4SBarry Smith /* Creating the index sets. */ 68*b4319ba4SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_local, &pcis->is_B_local);CHKERRQ(ierr); 69*b4319ba4SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,pcis->n_B,idx_B_global,&pcis->is_B_global);CHKERRQ(ierr); 70*b4319ba4SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_local, &pcis->is_I_local);CHKERRQ(ierr); 71*b4319ba4SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF,n_I ,idx_I_global,&pcis->is_I_global);CHKERRQ(ierr); 72*b4319ba4SBarry Smith /* Freeing memory and restoring arrays */ 73*b4319ba4SBarry Smith ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 74*b4319ba4SBarry Smith ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 75*b4319ba4SBarry Smith ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 76*b4319ba4SBarry Smith } 77*b4319ba4SBarry Smith 78*b4319ba4SBarry Smith /* 79*b4319ba4SBarry Smith Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 80*b4319ba4SBarry Smith is such that interior nodes come first than the interface ones, we have 81*b4319ba4SBarry Smith 82*b4319ba4SBarry Smith [ | ] 83*b4319ba4SBarry Smith [ A_II | A_IB ] 84*b4319ba4SBarry Smith A = [ | ] 85*b4319ba4SBarry Smith [-----------+------] 86*b4319ba4SBarry Smith [ A_BI | A_BB ] 87*b4319ba4SBarry Smith */ 88*b4319ba4SBarry Smith 89*b4319ba4SBarry Smith ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr); 90*b4319ba4SBarry Smith ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 91*b4319ba4SBarry Smith ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 92*b4319ba4SBarry Smith ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,PETSC_DECIDE,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 93*b4319ba4SBarry Smith 94*b4319ba4SBarry Smith /* 95*b4319ba4SBarry Smith Creating work vectors and arrays 96*b4319ba4SBarry Smith */ 97*b4319ba4SBarry Smith /* pcis->vec1_N has already been created */ 98*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 99*b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 100*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 101*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 102*b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 103*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 104*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 105*b4319ba4SBarry Smith { 106*b4319ba4SBarry Smith Vec global; 107*b4319ba4SBarry Smith ierr = PCGetVector(pc,&global);CHKERRQ(ierr); 108*b4319ba4SBarry Smith ierr = VecDuplicate(global,&pcis->vec1_global);CHKERRQ(ierr); 109*b4319ba4SBarry Smith } 110*b4319ba4SBarry Smith ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr); 111*b4319ba4SBarry Smith 112*b4319ba4SBarry Smith /* Creating the scatter contexts */ 113*b4319ba4SBarry Smith ierr = VecScatterCreate(pc->vec,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 114*b4319ba4SBarry Smith ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 115*b4319ba4SBarry Smith ierr = VecScatterCreate(pc->vec,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 116*b4319ba4SBarry Smith 117*b4319ba4SBarry Smith /* Creating scaling "matrix" D, from information in vec1_N */ 118*b4319ba4SBarry Smith ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 119*b4319ba4SBarry Smith ierr = VecScatterBegin(pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr); 120*b4319ba4SBarry Smith ierr = VecScatterEnd (pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr); 121*b4319ba4SBarry Smith ierr = VecReciprocal(pcis->D);CHKERRQ(ierr); 122*b4319ba4SBarry Smith 123*b4319ba4SBarry Smith /* See historical note 01, at the bottom of this file. */ 124*b4319ba4SBarry Smith 125*b4319ba4SBarry Smith /* 126*b4319ba4SBarry Smith Creating the KSP contexts for the local Dirichlet and Neumann problems. 127*b4319ba4SBarry Smith */ 128*b4319ba4SBarry Smith { 129*b4319ba4SBarry Smith PC pc_ctx; 130*b4319ba4SBarry Smith /* Dirichlet */ 131*b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 132*b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr); 133*b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 134*b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 135*b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 136*b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 137*b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 138*b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 139*b4319ba4SBarry Smith ierr = KSPSetRhs(pcis->ksp_D,pcis->vec1_D);CHKERRQ(ierr); 140*b4319ba4SBarry Smith ierr = KSPSetSolution(pcis->ksp_D,pcis->vec2_D);CHKERRQ(ierr); 141*b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 142*b4319ba4SBarry Smith /* Neumann */ 143*b4319ba4SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 144*b4319ba4SBarry Smith ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr); 145*b4319ba4SBarry Smith ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 146*b4319ba4SBarry Smith ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 147*b4319ba4SBarry Smith ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 148*b4319ba4SBarry Smith ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 149*b4319ba4SBarry Smith ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 150*b4319ba4SBarry Smith { 151*b4319ba4SBarry Smith PetscTruth damp_fixed, 152*b4319ba4SBarry Smith remove_nullspace_fixed, 153*b4319ba4SBarry Smith set_damping_factor_floating, 154*b4319ba4SBarry Smith not_damp_floating, 155*b4319ba4SBarry Smith not_remove_nullspace_floating; 156*b4319ba4SBarry Smith PetscReal fixed_factor, 157*b4319ba4SBarry Smith floating_factor; 158*b4319ba4SBarry Smith 159*b4319ba4SBarry Smith ierr = PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 160*b4319ba4SBarry Smith if (!damp_fixed) { fixed_factor = 0.0; } 161*b4319ba4SBarry Smith ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_damp_fixed",&damp_fixed);CHKERRQ(ierr); 162*b4319ba4SBarry Smith 163*b4319ba4SBarry Smith ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed);CHKERRQ(ierr); 164*b4319ba4SBarry Smith 165*b4319ba4SBarry Smith ierr = PetscOptionsGetReal(pc_ctx->prefix,"-pc_is_set_damping_factor_floating", 166*b4319ba4SBarry Smith &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 167*b4319ba4SBarry Smith if (!set_damping_factor_floating) { floating_factor = 0.0; } 168*b4319ba4SBarry Smith ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating);CHKERRQ(ierr); 169*b4319ba4SBarry Smith if (!set_damping_factor_floating) { floating_factor = 1.e-12; } 170*b4319ba4SBarry Smith 171*b4319ba4SBarry Smith ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_damp_floating",¬_damp_floating);CHKERRQ(ierr); 172*b4319ba4SBarry Smith 173*b4319ba4SBarry Smith ierr = PetscOptionsHasName(pc_ctx->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating);CHKERRQ(ierr); 174*b4319ba4SBarry Smith 175*b4319ba4SBarry Smith if (pcis->pure_neumann) { /* floating subdomain */ 176*b4319ba4SBarry Smith if (!(not_damp_floating)) { 177*b4319ba4SBarry Smith ierr = PCLUSetDamping (pc_ctx,floating_factor);CHKERRQ(ierr); 178*b4319ba4SBarry Smith ierr = PCILUSetDamping(pc_ctx,floating_factor);CHKERRQ(ierr); 179*b4319ba4SBarry Smith } 180*b4319ba4SBarry Smith if (!(not_remove_nullspace_floating)){ 181*b4319ba4SBarry Smith MatNullSpace nullsp; 182*b4319ba4SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,1,0,PETSC_NULL,&nullsp);CHKERRQ(ierr); 183*b4319ba4SBarry Smith ierr = PCNullSpaceAttach(pc_ctx,nullsp);CHKERRQ(ierr); 184*b4319ba4SBarry Smith ierr = MatNullSpaceDestroy(nullsp);CHKERRQ(ierr); 185*b4319ba4SBarry Smith } 186*b4319ba4SBarry Smith } else { /* fixed subdomain */ 187*b4319ba4SBarry Smith if (damp_fixed) { 188*b4319ba4SBarry Smith ierr = PCLUSetDamping (pc_ctx,fixed_factor);CHKERRQ(ierr); 189*b4319ba4SBarry Smith ierr = PCILUSetDamping(pc_ctx,fixed_factor);CHKERRQ(ierr); 190*b4319ba4SBarry Smith } 191*b4319ba4SBarry Smith if (remove_nullspace_fixed) { 192*b4319ba4SBarry Smith MatNullSpace nullsp; 193*b4319ba4SBarry Smith ierr = MatNullSpaceCreate(PETSC_COMM_SELF,1,0,PETSC_NULL,&nullsp);CHKERRQ(ierr); 194*b4319ba4SBarry Smith ierr = PCNullSpaceAttach(pc_ctx,nullsp);CHKERRQ(ierr); 195*b4319ba4SBarry Smith ierr = MatNullSpaceDestroy(nullsp);CHKERRQ(ierr); 196*b4319ba4SBarry Smith } 197*b4319ba4SBarry Smith } 198*b4319ba4SBarry Smith } 199*b4319ba4SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 200*b4319ba4SBarry Smith ierr = KSPSetRhs(pcis->ksp_N,pcis->vec1_N);CHKERRQ(ierr); 201*b4319ba4SBarry Smith ierr = KSPSetSolution(pcis->ksp_N,pcis->vec2_N);CHKERRQ(ierr); 202*b4319ba4SBarry Smith ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 203*b4319ba4SBarry Smith } 204*b4319ba4SBarry Smith 205*b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh), 206*b4319ba4SBarry Smith &(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 207*b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE; 208*b4319ba4SBarry Smith 209*b4319ba4SBarry Smith PetscFunctionReturn(0); 210*b4319ba4SBarry Smith } 211*b4319ba4SBarry Smith 212*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 213*b4319ba4SBarry Smith /* 214*b4319ba4SBarry Smith PCISDestroy - 215*b4319ba4SBarry Smith */ 216*b4319ba4SBarry Smith #undef __FUNCT__ 217*b4319ba4SBarry Smith #define __FUNCT__ "PCISDestroy" 218*b4319ba4SBarry Smith int PCISDestroy(PC pc) 219*b4319ba4SBarry Smith { 220*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 221*b4319ba4SBarry Smith int ierr; 222*b4319ba4SBarry Smith 223*b4319ba4SBarry Smith PetscFunctionBegin; 224*b4319ba4SBarry Smith 225*b4319ba4SBarry Smith if (pcis->is_B_local) {ierr = ISDestroy(pcis->is_B_local);CHKERRQ(ierr);} 226*b4319ba4SBarry Smith if (pcis->is_I_local) {ierr = ISDestroy(pcis->is_I_local);CHKERRQ(ierr);} 227*b4319ba4SBarry Smith if (pcis->is_B_global) {ierr = ISDestroy(pcis->is_B_global);CHKERRQ(ierr);} 228*b4319ba4SBarry Smith if (pcis->is_I_global) {ierr = ISDestroy(pcis->is_I_global);CHKERRQ(ierr);} 229*b4319ba4SBarry Smith if (pcis->A_II) {ierr = MatDestroy(pcis->A_II);CHKERRQ(ierr);} 230*b4319ba4SBarry Smith if (pcis->A_IB) {ierr = MatDestroy(pcis->A_IB);CHKERRQ(ierr);} 231*b4319ba4SBarry Smith if (pcis->A_BI) {ierr = MatDestroy(pcis->A_BI);CHKERRQ(ierr);} 232*b4319ba4SBarry Smith if (pcis->A_BB) {ierr = MatDestroy(pcis->A_BB);CHKERRQ(ierr);} 233*b4319ba4SBarry Smith if (pcis->D) {ierr = VecDestroy(pcis->D);CHKERRQ(ierr);} 234*b4319ba4SBarry Smith if (pcis->ksp_N) {ierr = KSPDestroy(pcis->ksp_N);CHKERRQ(ierr);} 235*b4319ba4SBarry Smith if (pcis->ksp_D) {ierr = KSPDestroy(pcis->ksp_D);CHKERRQ(ierr);} 236*b4319ba4SBarry Smith if (pcis->vec1_N) {ierr = VecDestroy(pcis->vec1_N);CHKERRQ(ierr);} 237*b4319ba4SBarry Smith if (pcis->vec2_N) {ierr = VecDestroy(pcis->vec2_N);CHKERRQ(ierr);} 238*b4319ba4SBarry Smith if (pcis->vec1_D) {ierr = VecDestroy(pcis->vec1_D);CHKERRQ(ierr);} 239*b4319ba4SBarry Smith if (pcis->vec2_D) {ierr = VecDestroy(pcis->vec2_D);CHKERRQ(ierr);} 240*b4319ba4SBarry Smith if (pcis->vec3_D) {ierr = VecDestroy(pcis->vec3_D);CHKERRQ(ierr);} 241*b4319ba4SBarry Smith if (pcis->vec1_B) {ierr = VecDestroy(pcis->vec1_B);CHKERRQ(ierr);} 242*b4319ba4SBarry Smith if (pcis->vec2_B) {ierr = VecDestroy(pcis->vec2_B);CHKERRQ(ierr);} 243*b4319ba4SBarry Smith if (pcis->vec3_B) {ierr = VecDestroy(pcis->vec3_B);CHKERRQ(ierr);} 244*b4319ba4SBarry Smith if (pcis->vec1_global) {ierr = VecDestroy(pcis->vec1_global);CHKERRQ(ierr);} 245*b4319ba4SBarry Smith if (pcis->work_N) {ierr = PetscFree(pcis->work_N);CHKERRQ(ierr);} 246*b4319ba4SBarry Smith if (pcis->global_to_D) {ierr = VecScatterDestroy(pcis->global_to_D);CHKERRQ(ierr);} 247*b4319ba4SBarry Smith if (pcis->N_to_B) {ierr = VecScatterDestroy(pcis->N_to_B);CHKERRQ(ierr);} 248*b4319ba4SBarry Smith if (pcis->global_to_B) {ierr = VecScatterDestroy(pcis->global_to_B);CHKERRQ(ierr);} 249*b4319ba4SBarry Smith if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) { 250*b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 251*b4319ba4SBarry Smith } 252*b4319ba4SBarry Smith 253*b4319ba4SBarry Smith PetscFunctionReturn(0); 254*b4319ba4SBarry Smith } 255*b4319ba4SBarry Smith 256*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 257*b4319ba4SBarry Smith /* 258*b4319ba4SBarry Smith PCISCreate - 259*b4319ba4SBarry Smith */ 260*b4319ba4SBarry Smith #undef __FUNCT__ 261*b4319ba4SBarry Smith #define __FUNCT__ "PCISCreate" 262*b4319ba4SBarry Smith int PCISCreate(PC pc) 263*b4319ba4SBarry Smith { 264*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 265*b4319ba4SBarry Smith 266*b4319ba4SBarry Smith PetscFunctionBegin; 267*b4319ba4SBarry Smith 268*b4319ba4SBarry Smith pcis->is_B_local = 0; 269*b4319ba4SBarry Smith pcis->is_I_local = 0; 270*b4319ba4SBarry Smith pcis->is_B_global = 0; 271*b4319ba4SBarry Smith pcis->is_I_global = 0; 272*b4319ba4SBarry Smith pcis->A_II = 0; 273*b4319ba4SBarry Smith pcis->A_IB = 0; 274*b4319ba4SBarry Smith pcis->A_BI = 0; 275*b4319ba4SBarry Smith pcis->A_BB = 0; 276*b4319ba4SBarry Smith pcis->D = 0; 277*b4319ba4SBarry Smith pcis->ksp_N = 0; 278*b4319ba4SBarry Smith pcis->ksp_D = 0; 279*b4319ba4SBarry Smith pcis->vec1_N = 0; 280*b4319ba4SBarry Smith pcis->vec2_N = 0; 281*b4319ba4SBarry Smith pcis->vec1_D = 0; 282*b4319ba4SBarry Smith pcis->vec2_D = 0; 283*b4319ba4SBarry Smith pcis->vec3_D = 0; 284*b4319ba4SBarry Smith pcis->vec1_B = 0; 285*b4319ba4SBarry Smith pcis->vec2_B = 0; 286*b4319ba4SBarry Smith pcis->vec3_B = 0; 287*b4319ba4SBarry Smith pcis->vec1_global = 0; 288*b4319ba4SBarry Smith pcis->work_N = 0; 289*b4319ba4SBarry Smith pcis->global_to_D = 0; 290*b4319ba4SBarry Smith pcis->N_to_B = 0; 291*b4319ba4SBarry Smith pcis->global_to_B = 0; 292*b4319ba4SBarry Smith pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE; 293*b4319ba4SBarry Smith 294*b4319ba4SBarry Smith PetscFunctionReturn(0); 295*b4319ba4SBarry Smith } 296*b4319ba4SBarry Smith 297*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 298*b4319ba4SBarry Smith /* 299*b4319ba4SBarry Smith PCISApplySchur - 300*b4319ba4SBarry Smith 301*b4319ba4SBarry Smith Input parameters: 302*b4319ba4SBarry Smith . pc - preconditioner context 303*b4319ba4SBarry Smith . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 304*b4319ba4SBarry Smith 305*b4319ba4SBarry Smith Output parameters: 306*b4319ba4SBarry Smith . vec1_B - result of Schur complement applied to chunk 307*b4319ba4SBarry Smith . vec2_B - garbage (used as work space), or null (and v is used as workspace) 308*b4319ba4SBarry Smith . vec1_D - garbage (used as work space) 309*b4319ba4SBarry Smith . vec2_D - garbage (used as work space) 310*b4319ba4SBarry Smith 311*b4319ba4SBarry Smith */ 312*b4319ba4SBarry Smith #undef __FUNCT__ 313*b4319ba4SBarry Smith #define __FUNCT__ "PCIterSuApplySchur" 314*b4319ba4SBarry Smith int PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 315*b4319ba4SBarry Smith { 316*b4319ba4SBarry Smith int ierr; 317*b4319ba4SBarry Smith PetscScalar m_one = -1.0; 318*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 319*b4319ba4SBarry Smith 320*b4319ba4SBarry Smith PetscFunctionBegin; 321*b4319ba4SBarry Smith 322*b4319ba4SBarry Smith if (vec2_B == (Vec)0) { vec2_B = v; } 323*b4319ba4SBarry Smith 324*b4319ba4SBarry Smith ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 325*b4319ba4SBarry Smith ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 326*b4319ba4SBarry Smith ierr = KSPSetRhs(pcis->ksp_D,vec1_D);CHKERRQ(ierr); 327*b4319ba4SBarry Smith ierr = KSPSetSolution(pcis->ksp_D,vec2_D);CHKERRQ(ierr); 328*b4319ba4SBarry Smith ierr = KSPSolve(pcis->ksp_D);CHKERRQ(ierr); 329*b4319ba4SBarry Smith ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 330*b4319ba4SBarry Smith ierr = VecAXPY(&m_one,vec2_B,vec1_B);CHKERRQ(ierr); 331*b4319ba4SBarry Smith 332*b4319ba4SBarry Smith PetscFunctionReturn(0); 333*b4319ba4SBarry Smith } 334*b4319ba4SBarry Smith 335*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 336*b4319ba4SBarry Smith /* 337*b4319ba4SBarry Smith PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 338*b4319ba4SBarry Smith including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 339*b4319ba4SBarry Smith mode. 340*b4319ba4SBarry Smith 341*b4319ba4SBarry Smith Input parameters: 342*b4319ba4SBarry Smith . pc - preconditioner context 343*b4319ba4SBarry Smith . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 344*b4319ba4SBarry Smith . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 345*b4319ba4SBarry Smith 346*b4319ba4SBarry Smith Output parameter: 347*b4319ba4SBarry Smith . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 348*b4319ba4SBarry Smith . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 349*b4319ba4SBarry Smith 350*b4319ba4SBarry Smith Notes: 351*b4319ba4SBarry Smith The entries in the array that do not correspond to interface nodes remain unaltered. 352*b4319ba4SBarry Smith */ 353*b4319ba4SBarry Smith #undef __FUNCT__ 354*b4319ba4SBarry Smith #define __FUNCT__ "PCISScatterArrayNToVecB" 355*b4319ba4SBarry Smith int PCISScatterArrayNToVecB (PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 356*b4319ba4SBarry Smith { 357*b4319ba4SBarry Smith int i, ierr, *idex; 358*b4319ba4SBarry Smith PetscScalar *array_B; 359*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 360*b4319ba4SBarry Smith 361*b4319ba4SBarry Smith PetscFunctionBegin; 362*b4319ba4SBarry Smith 363*b4319ba4SBarry Smith ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 364*b4319ba4SBarry Smith ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 365*b4319ba4SBarry Smith 366*b4319ba4SBarry Smith if (smode == SCATTER_FORWARD) { 367*b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 368*b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_B[i] = array_N[idex[i]]; } 369*b4319ba4SBarry Smith } else { /* ADD_VALUES */ 370*b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_B[i] += array_N[idex[i]]; } 371*b4319ba4SBarry Smith } 372*b4319ba4SBarry Smith } else { /* SCATTER_REVERSE */ 373*b4319ba4SBarry Smith if (imode == INSERT_VALUES) { 374*b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] = array_B[i]; } 375*b4319ba4SBarry Smith } else { /* ADD_VALUES */ 376*b4319ba4SBarry Smith for (i=0; i<pcis->n_B; i++) { array_N[idex[i]] += array_B[i]; } 377*b4319ba4SBarry Smith } 378*b4319ba4SBarry Smith } 379*b4319ba4SBarry Smith 380*b4319ba4SBarry Smith ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 381*b4319ba4SBarry Smith ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 382*b4319ba4SBarry Smith 383*b4319ba4SBarry Smith PetscFunctionReturn(0); 384*b4319ba4SBarry Smith } 385*b4319ba4SBarry Smith 386*b4319ba4SBarry Smith /* -------------------------------------------------------------------------- */ 387*b4319ba4SBarry Smith /* 388*b4319ba4SBarry Smith PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 389*b4319ba4SBarry Smith More precisely, solves the problem: 390*b4319ba4SBarry Smith [ A_II A_IB ] [ . ] [ 0 ] 391*b4319ba4SBarry Smith [ ] [ ] = [ ] 392*b4319ba4SBarry Smith [ A_BI A_BB ] [ x ] [ b ] 393*b4319ba4SBarry Smith 394*b4319ba4SBarry Smith Input parameters: 395*b4319ba4SBarry Smith . pc - preconditioner context 396*b4319ba4SBarry Smith . b - vector of local interface nodes (including ghosts) 397*b4319ba4SBarry Smith 398*b4319ba4SBarry Smith Output parameters: 399*b4319ba4SBarry Smith . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 400*b4319ba4SBarry Smith complement to b 401*b4319ba4SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 402*b4319ba4SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 403*b4319ba4SBarry Smith 404*b4319ba4SBarry Smith */ 405*b4319ba4SBarry Smith #undef __FUNCT__ 406*b4319ba4SBarry Smith #define __FUNCT__ "PCISApplyInvSchur" 407*b4319ba4SBarry Smith int PCISApplyInvSchur (PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 408*b4319ba4SBarry Smith { 409*b4319ba4SBarry Smith int ierr; 410*b4319ba4SBarry Smith PC_IS *pcis = (PC_IS*)(pc->data); 411*b4319ba4SBarry Smith PetscScalar zero = 0.0; 412*b4319ba4SBarry Smith 413*b4319ba4SBarry Smith PetscFunctionBegin; 414*b4319ba4SBarry Smith 415*b4319ba4SBarry Smith /* 416*b4319ba4SBarry Smith Neumann solvers. 417*b4319ba4SBarry Smith Applying the inverse of the local Schur complement, i.e, solving a Neumann 418*b4319ba4SBarry Smith Problem with zero at the interior nodes of the RHS and extracting the interface 419*b4319ba4SBarry Smith part of the solution. inverse Schur complement is applied to b and the result 420*b4319ba4SBarry Smith is stored in x. 421*b4319ba4SBarry Smith */ 422*b4319ba4SBarry Smith /* Setting the RHS vec1_N */ 423*b4319ba4SBarry Smith ierr = VecSet(&zero,vec1_N);CHKERRQ(ierr); 424*b4319ba4SBarry Smith ierr = VecScatterBegin(b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);CHKERRQ(ierr); 425*b4319ba4SBarry Smith ierr = VecScatterEnd (b,vec1_N,INSERT_VALUES,SCATTER_REVERSE,pcis->N_to_B);CHKERRQ(ierr); 426*b4319ba4SBarry Smith /* Checking for consistency of the RHS */ 427*b4319ba4SBarry Smith { 428*b4319ba4SBarry Smith PetscTruth flg; 429*b4319ba4SBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-pc_is_check_consistency",&flg);CHKERRQ(ierr); 430*b4319ba4SBarry Smith if (flg) { 431*b4319ba4SBarry Smith PetscScalar average; 432*b4319ba4SBarry Smith ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 433*b4319ba4SBarry Smith average = average / ((PetscReal)pcis->n); 434*b4319ba4SBarry Smith if (pcis->pure_neumann) { 435*b4319ba4SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is floating. Average = % 1.14e\n", 436*b4319ba4SBarry Smith PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 437*b4319ba4SBarry Smith } else { 438*b4319ba4SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_(pc->comm),"Subdomain %04d is fixed. Average = % 1.14e\n", 439*b4319ba4SBarry Smith PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 440*b4319ba4SBarry Smith } 441*b4319ba4SBarry Smith PetscViewerFlush(PETSC_VIEWER_STDOUT_(pc->comm)); 442*b4319ba4SBarry Smith } 443*b4319ba4SBarry Smith } 444*b4319ba4SBarry Smith /* Solving the system for vec2_N */ 445*b4319ba4SBarry Smith ierr = KSPSetRhs(pcis->ksp_N,vec1_N);CHKERRQ(ierr); 446*b4319ba4SBarry Smith ierr = KSPSetSolution(pcis->ksp_N,vec2_N);CHKERRQ(ierr); 447*b4319ba4SBarry Smith ierr = KSPSolve(pcis->ksp_N);CHKERRQ(ierr); 448*b4319ba4SBarry Smith /* Extracting the local interface vector out of the solution */ 449*b4319ba4SBarry Smith ierr = VecScatterBegin(vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr); 450*b4319ba4SBarry Smith ierr = VecScatterEnd (vec2_N,x,INSERT_VALUES,SCATTER_FORWARD,pcis->N_to_B);CHKERRQ(ierr); 451*b4319ba4SBarry Smith 452*b4319ba4SBarry Smith PetscFunctionReturn(0); 453*b4319ba4SBarry Smith } 454*b4319ba4SBarry Smith 455*b4319ba4SBarry Smith 456*b4319ba4SBarry Smith 457*b4319ba4SBarry Smith 458*b4319ba4SBarry Smith 459*b4319ba4SBarry Smith 460*b4319ba4SBarry Smith 461*b4319ba4SBarry Smith 462*b4319ba4SBarry Smith 463