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