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 = VecDestroy(&pcis->D);CHKERRQ(ierr); 52 ierr = PetscObjectReference((PetscObject)scaling_factors);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) 132 { 133 PC_IS *pcis = (PC_IS*)(pc->data); 134 Mat_IS *matis = (Mat_IS*)pc->mat->data; 135 PetscErrorCode ierr; 136 PetscBool flg; 137 Vec counter; 138 139 PetscFunctionBegin; 140 ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&flg);CHKERRQ(ierr); 141 if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); 142 143 pcis->pure_neumann = matis->pure_neumann; 144 145 /* get info on mapping */ 146 ierr = ISLocalToGlobalMappingGetInfo(((Mat_IS*)(pc->mat->data))->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 147 pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_TRUE; 148 149 /* Creating local and global index sets for interior and inteface nodes. */ 150 { 151 PetscInt n_I; 152 PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; 153 PetscInt *array; 154 PetscInt i,j; 155 156 /* Identifying interior and interface nodes, in local numbering */ 157 ierr = VecGetSize(matis->x,&pcis->n);CHKERRQ(ierr); 158 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array);CHKERRQ(ierr); 159 ierr = PetscMemzero(array,pcis->n*sizeof(PetscInt));CHKERRQ(ierr); 160 for (i=0;i<pcis->n_neigh;i++) 161 for (j=0;j<pcis->n_shared[i];j++) 162 array[pcis->shared[i][j]] += 1; 163 164 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_I_local);CHKERRQ(ierr); 165 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&idx_B_local);CHKERRQ(ierr); 166 for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { 167 if (!array[i]) { 168 idx_I_local[n_I] = i; 169 n_I++; 170 } else { 171 idx_B_local[pcis->n_B] = i; 172 pcis->n_B++; 173 } 174 } 175 /* Getting the global numbering */ 176 idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ 177 idx_I_global = idx_B_local + pcis->n_B; 178 ierr = ISLocalToGlobalMappingApply(matis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); 179 ierr = ISLocalToGlobalMappingApply(matis->mapping,n_I, idx_I_local,idx_I_global);CHKERRQ(ierr); 180 181 /* Creating the index sets. */ 182 ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); 183 ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); 184 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); 185 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); 186 187 /* Freeing memory and restoring arrays */ 188 ierr = PetscFree(idx_B_local);CHKERRQ(ierr); 189 ierr = PetscFree(idx_I_local);CHKERRQ(ierr); 190 ierr = PetscFree(array);CHKERRQ(ierr); 191 } 192 193 /* 194 Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering 195 is such that interior nodes come first than the interface ones, we have 196 197 [ | ] 198 [ A_II | A_IB ] 199 A = [ | ] 200 [-----------+------] 201 [ A_BI | A_BB ] 202 */ 203 204 ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_II);CHKERRQ(ierr); 205 ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 206 ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 207 ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 208 209 /* 210 Creating work vectors and arrays 211 */ 212 ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); 213 ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); 214 ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); 215 ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); 216 ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); 217 ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); 218 ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); 219 ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); 220 ierr = MatGetVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); 221 ierr = PetscMalloc((pcis->n)*sizeof(PetscScalar),&pcis->work_N);CHKERRQ(ierr); 222 223 /* Creating the scatter contexts */ 224 ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); 225 ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); 226 ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); 227 228 /* Creating scaling "matrix" D */ 229 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); 230 if (!pcis->D) { 231 ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); 232 if (!pcis->use_stiffness_scaling) { 233 ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); 234 } else { 235 ierr = MatGetDiagonal(matis->A,pcis->vec1_N);CHKERRQ(ierr); 236 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 237 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 238 } 239 } 240 ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 241 ierr = MatGetVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ 242 ierr = VecSet(counter,0.0);CHKERRQ(ierr); 243 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 244 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 245 ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 246 ierr = VecScatterEnd (pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 247 ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 248 ierr = VecDestroy(&counter);CHKERRQ(ierr); 249 250 /* See historical note 01, at the bottom of this file. */ 251 252 /* 253 Creating the KSP contexts for the local Dirichlet and Neumann problems. 254 */ 255 { 256 PC pc_ctx; 257 /* Dirichlet */ 258 ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); 259 ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 260 ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II,SAME_PRECONDITIONER);CHKERRQ(ierr); 261 ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); 262 ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); 263 ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 264 ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); 265 ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); 266 /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 267 ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); 268 /* Neumann */ 269 ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); 270 ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); 271 ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A,SAME_PRECONDITIONER);CHKERRQ(ierr); 272 ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); 273 ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); 274 ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); 275 ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); 276 ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); 277 { 278 PetscBool damp_fixed = PETSC_FALSE, 279 remove_nullspace_fixed = PETSC_FALSE, 280 set_damping_factor_floating = PETSC_FALSE, 281 not_damp_floating = PETSC_FALSE, 282 not_remove_nullspace_floating = PETSC_FALSE; 283 PetscReal fixed_factor, 284 floating_factor; 285 286 ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); 287 if (!damp_fixed) fixed_factor = 0.0; 288 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr); 289 290 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr); 291 292 ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", 293 &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); 294 if (!set_damping_factor_floating) floating_factor = 0.0; 295 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr); 296 if (!set_damping_factor_floating) floating_factor = 1.e-12; 297 298 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,NULL);CHKERRQ(ierr); 299 300 ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,NULL);CHKERRQ(ierr); 301 302 if (pcis->pure_neumann) { /* floating subdomain */ 303 if (!(not_damp_floating)) { 304 ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 305 ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 306 } 307 if (!(not_remove_nullspace_floating)) { 308 MatNullSpace nullsp; 309 ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 310 ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr); 311 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 312 } 313 } else { /* fixed subdomain */ 314 if (damp_fixed) { 315 ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); 316 ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); 317 } 318 if (remove_nullspace_fixed) { 319 MatNullSpace nullsp; 320 ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); 321 ierr = KSPSetNullSpace(pcis->ksp_N,nullsp);CHKERRQ(ierr); 322 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 323 } 324 } 325 } 326 /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 327 ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); 328 } 329 330 PetscFunctionReturn(0); 331 } 332 333 /* -------------------------------------------------------------------------- */ 334 /* 335 PCISDestroy - 336 */ 337 #undef __FUNCT__ 338 #define __FUNCT__ "PCISDestroy" 339 PetscErrorCode PCISDestroy(PC pc) 340 { 341 PC_IS *pcis = (PC_IS*)(pc->data); 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 ierr = ISDestroy(&pcis->is_B_local);CHKERRQ(ierr); 346 ierr = ISDestroy(&pcis->is_I_local);CHKERRQ(ierr); 347 ierr = ISDestroy(&pcis->is_B_global);CHKERRQ(ierr); 348 ierr = ISDestroy(&pcis->is_I_global);CHKERRQ(ierr); 349 ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 350 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 351 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 352 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 353 ierr = VecDestroy(&pcis->D);CHKERRQ(ierr); 354 ierr = KSPDestroy(&pcis->ksp_N);CHKERRQ(ierr); 355 ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 356 ierr = VecDestroy(&pcis->vec1_N);CHKERRQ(ierr); 357 ierr = VecDestroy(&pcis->vec2_N);CHKERRQ(ierr); 358 ierr = VecDestroy(&pcis->vec1_D);CHKERRQ(ierr); 359 ierr = VecDestroy(&pcis->vec2_D);CHKERRQ(ierr); 360 ierr = VecDestroy(&pcis->vec3_D);CHKERRQ(ierr); 361 ierr = VecDestroy(&pcis->vec1_B);CHKERRQ(ierr); 362 ierr = VecDestroy(&pcis->vec2_B);CHKERRQ(ierr); 363 ierr = VecDestroy(&pcis->vec3_B);CHKERRQ(ierr); 364 ierr = VecDestroy(&pcis->vec1_global);CHKERRQ(ierr); 365 ierr = VecScatterDestroy(&pcis->global_to_D);CHKERRQ(ierr); 366 ierr = VecScatterDestroy(&pcis->N_to_B);CHKERRQ(ierr); 367 ierr = VecScatterDestroy(&pcis->global_to_B);CHKERRQ(ierr); 368 ierr = PetscFree(pcis->work_N);CHKERRQ(ierr); 369 if (pcis->ISLocalToGlobalMappingGetInfoWasCalled) { 370 ierr = ISLocalToGlobalMappingRestoreInfo((ISLocalToGlobalMapping)0,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); 371 } 372 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",NULL);CHKERRQ(ierr); 373 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",NULL);CHKERRQ(ierr); 374 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",NULL);CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 /* -------------------------------------------------------------------------- */ 379 /* 380 PCISCreate - 381 */ 382 #undef __FUNCT__ 383 #define __FUNCT__ "PCISCreate" 384 PetscErrorCode PCISCreate(PC pc) 385 { 386 PC_IS *pcis = (PC_IS*)(pc->data); 387 PetscErrorCode ierr; 388 389 PetscFunctionBegin; 390 pcis->is_B_local = 0; 391 pcis->is_I_local = 0; 392 pcis->is_B_global = 0; 393 pcis->is_I_global = 0; 394 pcis->A_II = 0; 395 pcis->A_IB = 0; 396 pcis->A_BI = 0; 397 pcis->A_BB = 0; 398 pcis->D = 0; 399 pcis->ksp_N = 0; 400 pcis->ksp_D = 0; 401 pcis->vec1_N = 0; 402 pcis->vec2_N = 0; 403 pcis->vec1_D = 0; 404 pcis->vec2_D = 0; 405 pcis->vec3_D = 0; 406 pcis->vec1_B = 0; 407 pcis->vec2_B = 0; 408 pcis->vec3_B = 0; 409 pcis->vec1_global = 0; 410 pcis->work_N = 0; 411 pcis->global_to_D = 0; 412 pcis->N_to_B = 0; 413 pcis->global_to_B = 0; 414 415 pcis->ISLocalToGlobalMappingGetInfoWasCalled = PETSC_FALSE; 416 417 pcis->scaling_factor = 1.0; 418 /* composing functions */ 419 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetUseStiffnessScaling_C",PCISSetUseStiffnessScaling_IS);CHKERRQ(ierr); 420 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainScalingFactor_C",PCISSetSubdomainScalingFactor_IS);CHKERRQ(ierr); 421 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCISSetSubdomainDiagonalScaling_C",PCISSetSubdomainDiagonalScaling_IS);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 /* -------------------------------------------------------------------------- */ 426 /* 427 PCISApplySchur - 428 429 Input parameters: 430 . pc - preconditioner context 431 . v - vector to which the Schur complement is to be applied (it is NOT modified inside this function, UNLESS vec2_B is null) 432 433 Output parameters: 434 . vec1_B - result of Schur complement applied to chunk 435 . vec2_B - garbage (used as work space), or null (and v is used as workspace) 436 . vec1_D - garbage (used as work space) 437 . vec2_D - garbage (used as work space) 438 439 */ 440 #undef __FUNCT__ 441 #define __FUNCT__ "PCISApplySchur" 442 PetscErrorCode PCISApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 443 { 444 PetscErrorCode ierr; 445 PC_IS *pcis = (PC_IS*)(pc->data); 446 447 PetscFunctionBegin; 448 if (!vec2_B) vec2_B = v; 449 450 ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 451 ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 452 ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 453 ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 454 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 458 /* -------------------------------------------------------------------------- */ 459 /* 460 PCISScatterArrayNToVecB - Scatters interface node values from a big array (of all local nodes, interior or interface, 461 including ghosts) into an interface vector, when in SCATTER_FORWARD mode, or vice-versa, when in SCATTER_REVERSE 462 mode. 463 464 Input parameters: 465 . pc - preconditioner context 466 . array_N - [when in SCATTER_FORWARD mode] Array to be scattered into the vector 467 . v_B - [when in SCATTER_REVERSE mode] Vector to be scattered into the array 468 469 Output parameter: 470 . array_N - [when in SCATTER_REVERSE mode] Array to receive the scattered vector 471 . v_B - [when in SCATTER_FORWARD mode] Vector to receive the scattered array 472 473 Notes: 474 The entries in the array that do not correspond to interface nodes remain unaltered. 475 */ 476 #undef __FUNCT__ 477 #define __FUNCT__ "PCISScatterArrayNToVecB" 478 PetscErrorCode PCISScatterArrayNToVecB(PetscScalar *array_N, Vec v_B, InsertMode imode, ScatterMode smode, PC pc) 479 { 480 PetscInt i; 481 const PetscInt *idex; 482 PetscErrorCode ierr; 483 PetscScalar *array_B; 484 PC_IS *pcis = (PC_IS*)(pc->data); 485 486 PetscFunctionBegin; 487 ierr = VecGetArray(v_B,&array_B);CHKERRQ(ierr); 488 ierr = ISGetIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 489 490 if (smode == SCATTER_FORWARD) { 491 if (imode == INSERT_VALUES) { 492 for (i=0; i<pcis->n_B; i++) array_B[i] = array_N[idex[i]]; 493 } else { /* ADD_VALUES */ 494 for (i=0; i<pcis->n_B; i++) array_B[i] += array_N[idex[i]]; 495 } 496 } else { /* SCATTER_REVERSE */ 497 if (imode == INSERT_VALUES) { 498 for (i=0; i<pcis->n_B; i++) array_N[idex[i]] = array_B[i]; 499 } else { /* ADD_VALUES */ 500 for (i=0; i<pcis->n_B; i++) array_N[idex[i]] += array_B[i]; 501 } 502 } 503 ierr = ISRestoreIndices(pcis->is_B_local,&idex);CHKERRQ(ierr); 504 ierr = VecRestoreArray(v_B,&array_B);CHKERRQ(ierr); 505 PetscFunctionReturn(0); 506 } 507 508 /* -------------------------------------------------------------------------- */ 509 /* 510 PCISApplyInvSchur - Solves the Neumann problem related to applying the inverse of the Schur complement. 511 More precisely, solves the problem: 512 [ A_II A_IB ] [ . ] [ 0 ] 513 [ ] [ ] = [ ] 514 [ A_BI A_BB ] [ x ] [ b ] 515 516 Input parameters: 517 . pc - preconditioner context 518 . b - vector of local interface nodes (including ghosts) 519 520 Output parameters: 521 . x - vector of local interface nodes (including ghosts); returns the application of the inverse of the Schur 522 complement to b 523 . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 524 . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 525 526 */ 527 #undef __FUNCT__ 528 #define __FUNCT__ "PCISApplyInvSchur" 529 PetscErrorCode PCISApplyInvSchur(PC pc, Vec b, Vec x, Vec vec1_N, Vec vec2_N) 530 { 531 PetscErrorCode ierr; 532 PC_IS *pcis = (PC_IS*)(pc->data); 533 534 PetscFunctionBegin; 535 /* 536 Neumann solvers. 537 Applying the inverse of the local Schur complement, i.e, solving a Neumann 538 Problem with zero at the interior nodes of the RHS and extracting the interface 539 part of the solution. inverse Schur complement is applied to b and the result 540 is stored in x. 541 */ 542 /* Setting the RHS vec1_N */ 543 ierr = VecSet(vec1_N,0.0);CHKERRQ(ierr); 544 ierr = VecScatterBegin(pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 545 ierr = VecScatterEnd (pcis->N_to_B,b,vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 546 /* Checking for consistency of the RHS */ 547 { 548 PetscBool flg = PETSC_FALSE; 549 ierr = PetscOptionsGetBool(NULL,"-pc_is_check_consistency",&flg,NULL);CHKERRQ(ierr); 550 if (flg) { 551 PetscScalar average; 552 PetscViewer viewer; 553 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 554 555 ierr = VecSum(vec1_N,&average);CHKERRQ(ierr); 556 average = average / ((PetscReal)pcis->n); 557 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 558 if (pcis->pure_neumann) { 559 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is floating. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 560 } else { 561 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Subdomain %04d is fixed. Average = % 1.14e\n",PetscGlobalRank,PetscAbsScalar(average));CHKERRQ(ierr); 562 } 563 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 564 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 565 } 566 } 567 /* Solving the system for vec2_N */ 568 ierr = KSPSolve(pcis->ksp_N,vec1_N,vec2_N);CHKERRQ(ierr); 569 /* Extracting the local interface vector out of the solution */ 570 ierr = VecScatterBegin(pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 571 ierr = VecScatterEnd (pcis->N_to_B,vec2_N,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 572 PetscFunctionReturn(0); 573 } 574 575 576 577 578 579 580 581 582 583