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