1 2 #include <../src/ksp/pc/impls/is/pcis.h> /*I "petscpc.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "PCISSetUseStiffnessScaling_IS" 6 static PetscErrorCode PCISSetUseStiffnessScaling_IS(PC pc, PetscBool use) 7 { 8 PC_IS *pcis = (PC_IS*)pc->data; 9 10 PetscFunctionBegin; 11 pcis->use_stiffness_scaling = use; 12 PetscFunctionReturn(0); 13 } 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "PCISSetUseStiffnessScaling" 17 /*@ 18 PCISSetUseStiffnessScaling - Tells PCIS to construct partition of unity using 19 local matrices' diagonal. 20 21 Not collective 22 23 Input Parameters: 24 + pc - the preconditioning context 25 - use - whether or not pcis use matrix diagonal to build partition of unity. 26 27 Level: intermediate 28 29 Notes: 30 31 .seealso: PCBDDC 32 @*/ 33 PetscErrorCode PCISSetUseStiffnessScaling(PC pc, PetscBool use) 34 { 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 39 ierr = PetscTryMethod(pc,"PCISSetUseStiffnessScaling_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "PCISSetSubdomainDiagonalScaling_IS" 45 static PetscErrorCode PCISSetSubdomainDiagonalScaling_IS(PC pc, Vec scaling_factors) 46 { 47 PetscErrorCode ierr; 48 PC_IS *pcis = (PC_IS*)pc->data; 49 50 PetscFunctionBegin; 51 ierr = PetscObjectReference((PetscObject)scaling_factors);CHKERRQ(ierr); 52 ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 53 pcis->D = scaling_factors; 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "PCISSetSubdomainDiagonalScaling" 59 /*@ 60 PCISSetSubdomainDiagonalScaling - Set diagonal scaling for PCIS. 61 62 Not collective 63 64 Input Parameters: 65 + pc - the preconditioning context 66 - scaling_factors - scaling factors for the subdomain 67 68 Level: intermediate 69 70 Notes: 71 Intended to use with jumping coefficients cases. 72 73 .seealso: PCBDDC 74 @*/ 75 PetscErrorCode PCISSetSubdomainDiagonalScaling(PC pc, Vec scaling_factors) 76 { 77 PetscErrorCode ierr; 78 79 PetscFunctionBegin; 80 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 81 ierr = PetscTryMethod(pc,"PCISSetSubdomainDiagonalScaling_C",(PC,Vec),(pc,scaling_factors));CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "PCISSetSubdomainScalingFactor_IS" 87 static PetscErrorCode PCISSetSubdomainScalingFactor_IS(PC pc, PetscScalar scal) 88 { 89 PC_IS *pcis = (PC_IS*)pc->data; 90 91 PetscFunctionBegin; 92 pcis->scaling_factor = scal; 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "PCISSetSubdomainScalingFactor" 98 /*@ 99 PCISSetSubdomainScalingFactor - Set scaling factor for PCIS. 100 101 Not collective 102 103 Input Parameters: 104 + pc - the preconditioning context 105 - scal - scaling factor for the subdomain 106 107 Level: intermediate 108 109 Notes: 110 Intended to use with jumping coefficients cases. 111 112 .seealso: PCBDDC 113 @*/ 114 PetscErrorCode PCISSetSubdomainScalingFactor(PC pc, PetscScalar scal) 115 { 116 PetscErrorCode ierr; 117 118 PetscFunctionBegin; 119 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 120 ierr = PetscTryMethod(pc,"PCISSetSubdomainScalingFactor_C",(PC,PetscScalar),(pc,scal));CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 125 /* -------------------------------------------------------------------------- */ 126 /* 127 PCISSetUp - 128 */ 129 #undef __FUNCT__ 130 #define __FUNCT__ "PCISSetUp" 131 PetscErrorCode PCISSetUp(PC pc, PetscBool computesolvers) 132 { 133 PC_IS *pcis = (PC_IS*)(pc->data); 134 Mat_IS *matis; 135 MatReuse reuse; 136 PetscErrorCode ierr; 137 PetscBool flg,issbaij; 138 Vec counter; 139 140 PetscFunctionBegin; 141 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr); 142 if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 143 matis = (Mat_IS*)pc->pmat->data; 144 145 /* first time creation, get info on substructuring */ 146 if (!pc->setupcalled) { 147 PetscInt n_I; 148 PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 149 PetscInt *array; 150 PetscInt i,j; 151 152 /* get info on mapping */ 153 ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr); 154 ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 155 pcis->mapping = pc->pmat->rmap->mapping; 156 ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr); 157 ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 158 159 /* Identifying interior and interface nodes, in local numbering */ 160 ierr = PetscMalloc1(pcis->n,&array);CHKERRQ(ierr); 161 ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr); 162 for (i=0;i<pcis->n_neigh;i++) 163 for (j=0;j<pcis->n_shared[i];j++) 164 array[pcis->shared[i][j]] += 1; 165 166 /* Creating local and global index sets for interior and inteface nodes. */ 167 ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr); 168 ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr); 169 for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 170 if (!array[i]) { 171 idx_I_local[n_I] = i; 172 n_I++; 173 } else { 174 idx_B_local[pcis->n_B] = i; 175 pcis->n_B++; 176 } 177 } 178 /* Getting the global numbering */ 179 idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 180 idx_I_global = idx_B_local + pcis->n_B; 181 ierr = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 182 ierr = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr); 183 184 /* Creating the index sets */ 185 ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); 186 ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); 187 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); 188 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); 189 190 /* Freeing memory */ 191 ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 192 ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 193 ierr = PetscFree(array);CHKERRQ(ierr); 194 195 /* Creating work vectors and arrays */ 196 ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 197 ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 198 ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 199 ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 200 ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 201 ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr); 202 ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 203 ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 204 ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 205 ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 206 ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr); 207 /* scaling vector */ 208 if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */ 209 ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 210 ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 211 } 212 213 /* Creating the scatter contexts */ 214 ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr); 215 ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 216 ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 217 ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 218 219 /* map from boundary to local */ 220 ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr); 221 } 222 223 /* 224 Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 225 is such that interior nodes come first than the interface ones, we have 226 227 [ A_II | A_IB ] 228 A = [------+------] 229 [ A_BI | A_BB ] 230 */ 231 reuse = MAT_INITIAL_MATRIX; 232 if (pcis->reusesubmatrices && pc->setupcalled) { 233 if (pc->flag == SAME_NONZERO_PATTERN) { 234 reuse = MAT_REUSE_MATRIX; 235 } else { 236 reuse = MAT_INITIAL_MATRIX; 237 } 238 } 239 if (reuse == MAT_INITIAL_MATRIX) { 240 ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 241 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 242 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 243 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 244 } 245 246 ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr); 247 ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr); 248 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 249 if (!issbaij) { 250 ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 251 ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 252 } else { 253 Mat newmat; 254 ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr); 255 ierr = MatGetSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 256 ierr = MatGetSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 257 ierr = MatDestroy(&newmat);CHKERRQ(ierr); 258 } 259 260 /* Creating scaling vector D */ 261 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); 262 if (pcis->use_stiffness_scaling) { 263 ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr); 264 } 265 ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ 266 ierr = VecSet(counter,0.0);CHKERRQ(ierr); 267 ierr = VecScatterBegin(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 268 ierr = VecScatterEnd(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 269 ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 270 ierr = VecScatterEnd(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 271 ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 272 ierr = VecDestroy(&counter);CHKERRQ(ierr); 273 274 /* See historical note 01, at the bottom of this file. */ 275 276 /* Creating the KSP contexts for the local Dirichlet and Neumann problems */ 277 if (computesolvers) { 278 PC pc_ctx; 279 280 pcis->pure_neumann = matis->pure_neumann; 281 /* Dirichlet */ 282 ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 283 ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr); 284 ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 285 ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr); 286 ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 287 ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 288 ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 289 ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 290 ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 291 /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 292 ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 293 /* Neumann */ 294 ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 295 ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr); 296 ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); 297 ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr); 298 ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 299 ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 300 ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 301 ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 302 ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 303 { 304 PetscBool damp_fixed = PETSC_FALSE, 305 remove_nullspace_fixed = PETSC_FALSE, 306 set_damping_factor_floating = PETSC_FALSE, 307 not_damp_floating = PETSC_FALSE, 308 not_remove_nullspace_floating = PETSC_FALSE; 309 PetscReal fixed_factor, 310 floating_factor; 311 312 ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 313 if (!damp_fixed) fixed_factor = 0.0; 314 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr); 315 316 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr); 317 318 ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", 319 &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 320 if (!set_damping_factor_floating) floating_factor = 0.0; 321 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr); 322 if (!set_damping_factor_floating) floating_factor = 1.e-12; 323 324 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,NULL);CHKERRQ(ierr); 325 326 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,NULL);CHKERRQ(ierr); 327 328 if (pcis->pure_neumann) { /* floating subdomain */ 329 if (!(not_damp_floating)) { 330 ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 331 ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 332 } 333 if (!(not_remove_nullspace_floating)) { 334 MatNullSpace nullsp; 335 ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 336 ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 337 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 338 } 339 } else { /* fixed subdomain */ 340 if (damp_fixed) { 341 ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 342 ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 343 } 344 if (remove_nullspace_fixed) { 345 MatNullSpace nullsp; 346 ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 347 ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); 348 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 349 } 350 } 351 } 352 /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 353 ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 354 } 355 356 PetscFunctionReturn(0); 357 } 358 359 /* -------------------------------------------------------------------------- */ 360 /* 361 PCISDestroy - 362 */ 363 #undef __FUNCT__ 364 #define __FUNCT__ "PCISDestroy" 365 PetscErrorCode PCISDestroy(PC pc) 366 { 367 PC_IS *pcis = (PC_IS*)(pc->data); 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 372 ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 373 ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 374 ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 375 ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 376 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 377 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 378 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 379 ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 380 ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 381 ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 382 ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 383 ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 384 ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 385 ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 386 ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 387 ierr = VecDestroy(&pcis->vec4_D);CHKERRQ(ierr); 388 ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 389 ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 390 ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 391 ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 392 ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 393 ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 394 ierr = VecScatterDestroy(&pcis->N_to_D);CHKERRQ(ierr); 395 ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 396 ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 397 if (pcis->n_neigh > -1) { 398 ierr = ISLocalToGlobalMappingRestoreInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 399 } 400 ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); 401 ierr = ISLocalToGlobalMappingDestroy(&pcis->BtoNmap);CHKERRQ(ierr); 402 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr); 403 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr); 404 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr); 405 PetscFunctionReturn(0); 406 } 407 408 /* -------------------------------------------------------------------------- */ 409 /* 410 PCISCreate - 411 */ 412 #undef __FUNCT__ 413 #define __FUNCT__ "PCISCreate" 414 PetscErrorCode PCISCreate(PC pc) 415 { 416 PC_IS *pcis = (PC_IS*)(pc->data); 417 PetscErrorCode ierr; 418 419 PetscFunctionBegin; 420 pcis->is_B_local = 0; 421 pcis->is_I_local = 0; 422 pcis->is_B_global = 0; 423 pcis->is_I_global = 0; 424 pcis->A_II = 0; 425 pcis->A_IB = 0; 426 pcis->A_BI = 0; 427 pcis->A_BB = 0; 428 pcis->D = 0; 429 pcis->ksp_N = 0; 430 pcis->ksp_D = 0; 431 pcis->vec1_N = 0; 432 pcis->vec2_N = 0; 433 pcis->vec1_D = 0; 434 pcis->vec2_D = 0; 435 pcis->vec3_D = 0; 436 pcis->vec1_B = 0; 437 pcis->vec2_B = 0; 438 pcis->vec3_B = 0; 439 pcis->vec1_global = 0; 440 pcis->work_N = 0; 441 pcis->global_to_D = 0; 442 pcis->N_to_B = 0; 443 pcis->N_to_D = 0; 444 pcis->global_to_B = 0; 445 pcis->mapping = 0; 446 pcis->BtoNmap = 0; 447 pcis->n_neigh = -1; 448 pcis->scaling_factor = 1.0; 449 pcis->reusesubmatrices = PETSC_TRUE; 450 /* composing functions */ 451 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr); 452 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 453 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 454 PetscFunctionReturn(0); 455 } 456 457 /* -------------------------------------------------------------------------- */ 458 /* 459 PCISApplySchur - 460 461 Input parameters: 462 . pc - preconditioner context 463 . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 464 465 Output parameters: 466 . vec1_B - result of Schur complement applied to chunk 467 . vec2_B - garbage (used as work space), or null (and v is used as workspace) 468 . vec1_D - garbage (used as work space) 469 . vec2_D - garbage (used as work space) 470 471 */ 472 #undef __FUNCT__ 473 #define __FUNCT__ "PCISApplySchur" 474 PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 475 { 476 PetscErrorCode ierr; 477 PC_IS *pcis = (PC_IS*)(pc->data); 478 479 PetscFunctionBegin; 480 if (!vec2_B) vec2_B = v; 481 482 ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 483 ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 484 ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 485 ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 486 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 487 PetscFunctionReturn(0); 488 } 489 490 /* -------------------------------------------------------------------------- */ 491 /* 492 PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 493 including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 494 mode. 495 496 Input parameters: 497 . pc - preconditioner context 498 . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 499 . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 500 501 Output parameter: 502 . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 503 . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 504 505 Notes: 506 The entries in the array that do not correspond to interface nodes remain unaltered. 507 */ 508 #undef __FUNCT__ 509 #define __FUNCT__ "PCISScatterArrayNToVecB" 510 PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 511 { 512 PetscInt i; 513 const PetscInt *idex; 514 PetscErrorCode ierr; 515 PetscScalar *array_B; 516 PC_IS *pcis = (PC_IS*)(pc->data); 517 518 PetscFunctionBegin; 519 ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 520 ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 521 522 if (smode == SCATTER_FORWARD) { 523 if (imode == INSERT_VALUES) { 524 for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]]; 525 } else { /* ADD_VALUES */ 526 for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]]; 527 } 528 } else { /* SCATTER_REVERSE */ 529 if (imode == INSERT_VALUES) { 530 for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i]; 531 } else { /* ADD_VALUES */ 532 for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i]; 533 } 534 } 535 ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 536 ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 537 PetscFunctionReturn(0); 538 } 539 540 /* -------------------------------------------------------------------------- */ 541 /* 542 PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 543 More precisely, solves the problem: 544 [ A_II A_IB ] [ . ] [ 0 ] 545 [ ] [ ] = [ ] 546 [ A_BI A_BB ] [ x ] [ b ] 547 548 Input parameters: 549 . pc - preconditioner context 550 . b - vector of local interface nodes (including ghosts) 551 552 Output parameters: 553 . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 554 complement to b 555 . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 556 . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 557 558 */ 559 #undef __FUNCT__ 560 #define __FUNCT__ "PCISApplyInvSchur" 561 PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 562 { 563 PetscErrorCode ierr; 564 PC_IS *pcis = (PC_IS*)(pc->data); 565 566 PetscFunctionBegin; 567 /* 568 Neumann solvers. 569 Applying the inverse of the local Schur complement, i.e, solving a Neumann 570 Problem with zero at the interior nodes of the RHS and extracting the interface 571 part of the solution. inverse Schur complement is applied to b and the result 572 is stored in x. 573 */ 574 /* Setting the RHS vec1_N */ 575 ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 576 ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 577 ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 578 /* Checking for consistency of the RHS */ 579 { 580 PetscBool flg = PETSC_FALSE; 581 ierr = PetscOptionsGetBool(NULL,NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr); 582 if (flg) { 583 PetscScalar average; 584 PetscViewer viewer; 585 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 586 587 ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 588 average = average / ((PetscReal)pcis->n); 589 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 590 if (pcis->pure_neumann) { 591 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 592 } else { 593 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 594 } 595 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 596 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 597 } 598 } 599 /* Solving the system for vec2_N */ 600 ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 601 /* Extracting the local interface vector out of the solution */ 602 ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 603 ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606 607 608 609 610 611 612 613 614 615