1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <petscblaslapack.h> 4 5 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 6 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*); 7 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC,Vec,Vec); 8 static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC,Vec,Vec); 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "PCBDDCMumpsCorrectionSolve_Private" 12 static PetscErrorCode PCBDDCMumpsCorrectionSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose) 13 { 14 PCBDDCReuseMumps ctx; 15 #if defined(PETSC_HAVE_MUMPS) 16 PetscInt ival; 17 #endif 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr); 22 #if defined(PETSC_HAVE_MUMPS) 23 ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr); 24 ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr); 25 #endif 26 if (transpose) { 27 ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr); 28 } else { 29 ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr); 30 } 31 #if defined(PETSC_HAVE_MUMPS) 32 ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr); 33 #endif 34 PetscFunctionReturn(0); 35 } 36 37 #undef __FUNCT__ 38 #define __FUNCT__ "PCBDDCMumpsCorrectionSolve" 39 static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol) 40 { 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr); 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose" 50 static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol) 51 { 52 PetscErrorCode ierr; 53 54 PetscFunctionBegin; 55 ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr); 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "PCBDDCReuseMumpsReset" 61 static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse) 62 { 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 ierr = MatDestroy(&reuse->F);CHKERRQ(ierr); 67 ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr); 68 ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr); 69 ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr); 70 ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr); 71 ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr); 72 ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr); 73 ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr); 74 ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr); 75 ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "PCBDDCMumpsInteriorSolve_Private" 81 static PetscErrorCode PCBDDCMumpsInteriorSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose) 82 { 83 PCBDDCReuseMumps ctx; 84 PetscScalar *array,*array_mumps; 85 #if defined(PETSC_HAVE_MUMPS) 86 PetscInt ival; 87 #endif 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr); 92 #if defined(PETSC_HAVE_MUMPS) 93 ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr); 94 ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr); 95 #endif 96 /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */ 97 ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 98 ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr); 99 ierr = PetscMemcpy(array_mumps,array,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr); 100 ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr); 101 ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 102 103 if (transpose) { 104 ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 105 } else { 106 ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 107 } 108 109 /* get back data to caller worskpace */ 110 ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr); 111 ierr = VecGetArray(sol,&array);CHKERRQ(ierr); 112 ierr = PetscMemcpy(array,array_mumps,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr); 113 ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr); 114 ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr); 115 #if defined(PETSC_HAVE_MUMPS) 116 ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr); 117 #endif 118 PetscFunctionReturn(0); 119 } 120 121 #undef __FUNCT__ 122 #define __FUNCT__ "PCBDDCMumpsInteriorSolve" 123 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol) 124 { 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr); 129 PetscFunctionReturn(0); 130 } 131 132 #undef __FUNCT__ 133 #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose" 134 static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol) 135 { 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr); 140 PetscFunctionReturn(0); 141 } 142 143 #undef __FUNCT__ 144 #define __FUNCT__ "PCBDDCComputeExplicitSchur" 145 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 146 { 147 Mat B, C, D, Bd, Cd, AinvBd; 148 KSP ksp; 149 PC pc; 150 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 151 PetscReal fill = 2.0; 152 PetscInt n_I; 153 PetscMPIInt size; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 158 if (size != 1) { 159 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 160 } 161 if (reuse == MAT_REUSE_MATRIX) { 162 PetscBool Sdense; 163 164 ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr); 165 if (!Sdense) { 166 SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense"); 167 } 168 } 169 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 170 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 171 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 172 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 173 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 174 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 175 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 176 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 177 ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr); 178 if (n_I) { 179 if (!Bdense) { 180 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 181 } else { 182 Bd = B; 183 } 184 185 if (isLU || isILU || isCHOL) { 186 Mat fact; 187 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 188 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 189 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 190 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 191 } else { 192 PetscBool ex = PETSC_TRUE; 193 194 if (ex) { 195 Mat Ainvd; 196 197 ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 198 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 199 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 200 } else { 201 Vec sol,rhs; 202 PetscScalar *arrayrhs,*arraysol; 203 PetscInt i,nrhs,n; 204 205 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 206 ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr); 207 ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr); 208 ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr); 209 ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr); 210 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 211 for (i=0;i<nrhs;i++) { 212 ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr); 213 ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr); 214 ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr); 215 ierr = VecResetArray(rhs);CHKERRQ(ierr); 216 ierr = VecResetArray(sol);CHKERRQ(ierr); 217 } 218 ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr); 219 ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr); 220 } 221 } 222 if (!Bdense & !issym) { 223 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 224 } 225 226 if (!issym) { 227 if (!Cdense) { 228 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 229 } else { 230 Cd = C; 231 } 232 ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 233 if (!Cdense) { 234 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 235 } 236 } else { 237 ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 238 if (!Bdense) { 239 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 240 } 241 } 242 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 243 } 244 245 if (D) { 246 Mat Dd; 247 PetscBool Ddense; 248 249 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 250 if (!Ddense) { 251 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 252 } else { 253 Dd = D; 254 } 255 if (n_I) { 256 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 257 } else { 258 if (reuse == MAT_INITIAL_MATRIX) { 259 ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr); 260 } else { 261 ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 262 } 263 } 264 if (!Ddense) { 265 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 266 } 267 } else { 268 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "PCBDDCSubSchursSetUp" 275 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers) 276 { 277 Mat F,A_II,A_IB,A_BI,A_BB,AE_II; 278 Mat S_all; 279 Mat global_schur_subsets,work_mat,*submats; 280 ISLocalToGlobalMapping l2gmap_subsets; 281 IS is_I,is_I_layer; 282 IS all_subsets,all_subsets_mult,all_subsets_n; 283 PetscInt *nnz,*all_local_idx_N; 284 PetscInt *auxnum1,*auxnum2; 285 PetscInt i,subset_size,max_subset_size; 286 PetscInt extra,local_size,global_size; 287 PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 288 PetscScalar *Bwork; 289 PetscSubcomm subcomm; 290 PetscMPIInt color,rank; 291 MPI_Comm comm_n; 292 PetscErrorCode ierr; 293 294 PetscFunctionBegin; 295 /* update info in sub_schurs */ 296 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 297 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 298 sub_schurs->is_hermitian = PETSC_FALSE; 299 sub_schurs->is_posdef = PETSC_FALSE; 300 if (Ain) { 301 PetscBool isseqaij; 302 /* determine if we are dealing with hermitian positive definite problems */ 303 #if !defined(PETSC_USE_COMPLEX) 304 if (Ain->symmetric_set) { 305 sub_schurs->is_hermitian = Ain->symmetric; 306 } 307 #else 308 if (Ain->hermitian_set) { 309 sub_schurs->is_hermitian = Ain->hermitian; 310 } 311 #endif 312 if (Ain->spd_set) { 313 sub_schurs->is_posdef = Ain->spd; 314 } 315 316 /* check */ 317 ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 318 if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) { 319 PetscInt lsize; 320 321 ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr); 322 if (lsize) { 323 PetscScalar val; 324 PetscReal norm; 325 Vec vec1,vec2,vec3; 326 327 ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr); 328 ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr); 329 ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr); 330 ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr); 331 #if !defined(PETSC_USE_COMPLEX) 332 ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr); 333 #else 334 ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr); 335 #endif 336 ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr); 337 ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr); 338 if (norm > PetscSqrtReal(PETSC_SMALL)) { 339 sub_schurs->is_hermitian = PETSC_FALSE; 340 } else { 341 sub_schurs->is_hermitian = PETSC_TRUE; 342 } 343 ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr); 344 if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE; 345 ierr = VecDestroy(&vec1);CHKERRQ(ierr); 346 ierr = VecDestroy(&vec2);CHKERRQ(ierr); 347 ierr = VecDestroy(&vec3);CHKERRQ(ierr); 348 } else { 349 sub_schurs->is_hermitian = PETSC_TRUE; 350 sub_schurs->is_posdef = PETSC_TRUE; 351 } 352 } 353 if (isseqaij) { 354 ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr); 355 sub_schurs->A = Ain; 356 } else { 357 ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr); 358 } 359 } 360 if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) { 361 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"General matrix pencils are not currently supported (%d,%d)",sub_schurs->is_hermitian,sub_schurs->is_posdef); 362 } 363 364 ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr); 365 sub_schurs->S = Sin; 366 if (sub_schurs->use_mumps) { 367 sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A); 368 } 369 370 /* preliminary checks */ 371 if (!sub_schurs->use_mumps && compute_Stilda) { 372 SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS"); 373 } 374 375 /* restrict work on active processes */ 376 color = 0; 377 if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */ 378 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr); 379 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr); 380 ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr); 381 ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr); 382 ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr); 383 ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr); 384 if (!sub_schurs->n_subs) { 385 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 386 PetscFunctionReturn(0); 387 } 388 389 /* get Schur complement matrices */ 390 if (!sub_schurs->use_mumps) { 391 Mat tA_IB,tA_BI,tA_BB; 392 PetscBool isseqsbaij; 393 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr); 394 ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 395 if (isseqsbaij) { 396 ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr); 397 ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr); 398 ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr); 399 } else { 400 ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr); 401 A_BB = tA_BB; 402 ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr); 403 A_IB = tA_IB; 404 ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr); 405 A_BI = tA_BI; 406 } 407 } else { 408 A_II = NULL; 409 A_IB = NULL; 410 A_BI = NULL; 411 A_BB = NULL; 412 } 413 S_all = NULL; 414 415 /* determine interior problems */ 416 ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 417 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 418 PetscBT touched; 419 const PetscInt* idx_B; 420 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 421 422 if (xadj == NULL || adjncy == NULL) { 423 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 424 } 425 /* get sizes */ 426 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 427 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 428 429 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 430 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 431 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 432 433 /* all boundary dofs must be skipped when adding layers */ 434 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 435 for (j=0;j<n_B;j++) { 436 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 437 } 438 ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 439 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 440 441 /* add prescribed number of layers of dofs */ 442 n_local_dofs = n_B; 443 n_prev_added = n_B; 444 for (layer=0;layer<nlayers;layer++) { 445 PetscInt n_added; 446 if (n_local_dofs == n_I+n_B) break; 447 if (n_local_dofs > n_I+n_B) { 448 SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error querying layer %d. Out of bound access (%d > %d)",layer,n_local_dofs,n_I+n_B); 449 } 450 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 451 n_prev_added = n_added; 452 n_local_dofs += n_added; 453 if (!n_added) break; 454 } 455 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 456 457 /* IS for I layer dofs in original numbering */ 458 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);CHKERRQ(ierr); 459 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 460 ierr = ISSort(is_I_layer);CHKERRQ(ierr); 461 /* IS for I layer dofs in I numbering */ 462 if (!sub_schurs->use_mumps) { 463 ISLocalToGlobalMapping ItoNmap; 464 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 465 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr); 466 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 467 468 /* II block */ 469 ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 470 } 471 } else { 472 PetscInt n_I; 473 474 /* IS for I dofs in original numbering */ 475 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 476 is_I_layer = sub_schurs->is_I; 477 478 /* IS for I dofs in I numbering (strided 1) */ 479 if (!sub_schurs->use_mumps) { 480 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 481 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 482 483 /* II block is the same */ 484 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 485 AE_II = A_II; 486 } 487 } 488 489 /* Get info on subset sizes and sum of all subsets sizes */ 490 max_subset_size = 0; 491 local_size = 0; 492 for (i=0;i<sub_schurs->n_subs;i++) { 493 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 494 max_subset_size = PetscMax(subset_size,max_subset_size); 495 local_size += subset_size; 496 } 497 498 /* Work arrays for local indices */ 499 extra = 0; 500 if (sub_schurs->use_mumps) { 501 ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr); 502 } 503 ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 504 if (extra) { 505 const PetscInt *idxs; 506 ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr); 507 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 508 ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr); 509 } 510 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 511 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr); 512 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr); 513 514 /* Get local indices in local numbering */ 515 local_size = 0; 516 for (i=0;i<sub_schurs->n_subs;i++) { 517 PetscInt j; 518 const PetscInt *idxs; 519 520 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 521 ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 522 /* start (smallest in global ordering) and multiplicity */ 523 auxnum1[i] = idxs[0]; 524 auxnum2[i] = subset_size; 525 /* subset indices in local numbering */ 526 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 527 ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 528 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 529 local_size += subset_size; 530 } 531 532 /* allocate extra workspace needed only for GETRI */ 533 Bwork = NULL; 534 pivots = NULL; 535 if (!sub_schurs->is_hermitian || !sub_schurs->is_posdef) { 536 PetscScalar lwork; 537 538 B_lwork = -1; 539 ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 540 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 541 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr)); 542 ierr = PetscFPTrapPop();CHKERRQ(ierr); 543 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 544 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 545 ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr); 546 } 547 548 /* prepare parallel matrices for summing up properly schurs on subsets */ 549 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr); 550 ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr); 551 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 552 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr); 553 ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr); 554 ierr = ISDestroy(&all_subsets);CHKERRQ(ierr); 555 ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr); 556 ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr); 557 if (i != local_size) { 558 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size); 559 } 560 ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr); 561 ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 562 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 563 ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 564 ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 565 ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 566 567 /* subset indices in local boundary numbering */ 568 if (!sub_schurs->is_Ej_all) { 569 PetscInt *all_local_idx_B; 570 571 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 572 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 573 if (subset_size != local_size) { 574 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 575 } 576 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 577 } 578 579 /* Local matrix of all local Schur on subsets (transposed) */ 580 if (!sub_schurs->S_Ej_all) { 581 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 582 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 583 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 584 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 585 } 586 587 /* Compute Schur complements explicitly */ 588 F = NULL; 589 if (!sub_schurs->use_mumps) { 590 Mat S_Ej_expl; 591 PetscScalar *work; 592 PetscInt j,*dummy_idx; 593 PetscBool Sdense; 594 595 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 596 local_size = 0; 597 for (i=0;i<sub_schurs->n_subs;i++) { 598 IS is_subset_B; 599 Mat AE_EE,AE_IE,AE_EI,S_Ej; 600 601 /* subsets in original and boundary numbering */ 602 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr); 603 /* EE block */ 604 ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr); 605 /* IE block */ 606 ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr); 607 /* EI block */ 608 if (sub_schurs->is_hermitian) { 609 ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 610 } else { 611 ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr); 612 } 613 ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr); 614 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr); 615 ierr = MatDestroy(&AE_EE);CHKERRQ(ierr); 616 ierr = MatDestroy(&AE_IE);CHKERRQ(ierr); 617 ierr = MatDestroy(&AE_EI);CHKERRQ(ierr); 618 if (AE_II == A_II) { /* we can reuse the same ksp */ 619 KSP ksp; 620 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 621 ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr); 622 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 623 KSP origksp,schurksp; 624 PC origpc,schurpc; 625 KSPType ksp_type; 626 PetscInt n_internal; 627 PetscBool ispcnone; 628 629 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 630 ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr); 631 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 632 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 633 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 634 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 635 ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr); 636 if (!ispcnone) { 637 PCType pc_type; 638 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 639 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 640 } else { 641 ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr); 642 } 643 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 644 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 645 MatSolverPackage solver=NULL; 646 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 647 if (solver) { 648 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 649 } 650 } 651 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 652 } 653 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 654 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 655 ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 656 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 657 if (Sdense) { 658 for (j=0;j<subset_size;j++) { 659 dummy_idx[j]=local_size+j; 660 } 661 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 662 } else { 663 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 664 } 665 ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 666 ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr); 667 local_size += subset_size; 668 } 669 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 670 /* free */ 671 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 672 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 673 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 674 } else { 675 Mat A; 676 IS is_A_all; 677 PetscScalar *work,*S_data; 678 PetscInt *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2; 679 PetscBool mumps_S; 680 681 /* get working mat */ 682 ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr); 683 if (!sub_schurs->is_dir) { 684 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 685 size_schur = local_size; 686 } else { 687 IS list[2]; 688 689 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr); 690 list[1] = sub_schurs->is_dir; 691 ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr); 692 ierr = ISDestroy(&list[0]);CHKERRQ(ierr); 693 ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr); 694 size_schur += local_size; 695 } 696 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 697 size_active_schur = local_size; /* size active schurs does not count any dirichlet dof on the interface */ 698 ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 699 ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr); 700 ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr); 701 ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 702 ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 703 704 if (n_I) { 705 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 706 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 707 } else { 708 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 709 } 710 /* subsets ordered last */ 711 ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr); 712 for (i=0;i<size_schur;i++) { 713 idxs_schur[i] = n_I+i+1; 714 } 715 #if defined(PETSC_HAVE_MUMPS) 716 ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr); 717 #endif 718 ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 719 720 /* factorization step */ 721 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 722 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 723 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 724 ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr); 725 #endif 726 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 727 } else { 728 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 729 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 730 ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr); 731 #endif 732 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 733 } 734 735 /* get explicit Schur Complement computed during numeric factorization */ 736 #if defined(PETSC_HAVE_MUMPS) 737 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 738 #endif 739 740 /* we can reuse the solvers if we are not using the economic version */ 741 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr); 742 reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all)); 743 mumps_S = PETSC_TRUE; 744 } else { /* we can't use MUMPS when size_schur == size_of_the_problem */ 745 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 746 reuse_solvers = PETSC_FALSE; 747 mumps_S = PETSC_FALSE; 748 } 749 750 if (reuse_solvers) { 751 Mat A_II; 752 Vec vec1_B; 753 PCBDDCReuseMumps msolv_ctx; 754 755 if (sub_schurs->reuse_mumps) { 756 ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr); 757 } else { 758 ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr); 759 } 760 msolv_ctx = sub_schurs->reuse_mumps; 761 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 762 ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr); 763 ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 764 msolv_ctx->F = F; 765 ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr); 766 767 /* interior solver */ 768 ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr); 769 ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr); 770 ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr); 771 ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr); 772 ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr); 773 ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr); 774 775 /* correction solver */ 776 ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr); 777 ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr); 778 ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr); 779 ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr); 780 ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr); 781 ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr); 782 783 /* scatter and vecs for Schur complement solver */ 784 ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr); 785 ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr); 786 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr); 787 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 788 ierr = VecDestroy(&vec1_B);CHKERRQ(ierr); 789 ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr); 790 msolv_ctx->is_R = is_A_all; 791 } 792 ierr = MatDestroy(&A);CHKERRQ(ierr); 793 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 794 795 /* Work arrays */ 796 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 797 798 /* matrices for adaptive selection */ 799 if (compute_Stilda) { 800 if (!sub_schurs->sum_S_Ej_tilda_all) { 801 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 802 ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr); 803 ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 804 ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 805 } 806 if (!sub_schurs->sum_S_Ej_inv_all) { 807 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 808 ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr); 809 ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr); 810 ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr); 811 } 812 } 813 814 /* S_Ej_all */ 815 cum = cum2 = 0; 816 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 817 for (i=0;i<sub_schurs->n_subs;i++) { 818 PetscInt j; 819 820 /* get S_E */ 821 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 822 if (mumps_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */ 823 PetscInt k; 824 for (k=0;k<subset_size;k++) { 825 for (j=k;j<subset_size;j++) { 826 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 827 work[j*subset_size+k] = S_data[cum2+k*size_schur+j]; 828 } 829 } 830 } else { /* copy to workspace */ 831 PetscInt k; 832 for (k=0;k<subset_size;k++) { 833 for (j=0;j<subset_size;j++) { 834 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 835 } 836 } 837 } 838 /* insert S_E values */ 839 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 840 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 841 842 /* if adaptivity is requested, invert S_E block */ 843 if (compute_Stilda) { 844 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 845 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 846 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */ 847 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr)); 848 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 849 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr)); 850 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 851 } else { 852 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr)); 853 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 854 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 855 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 856 } 857 ierr = PetscFPTrapPop();CHKERRQ(ierr); 858 ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 859 } 860 cum += subset_size; 861 cum2 += subset_size*(size_schur + 1); 862 } 863 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 864 865 #if defined(PETSC_HAVE_MUMPS) 866 if (mumps_S) { 867 ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr); 868 } 869 #endif 870 871 if (compute_Stilda) { 872 if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */ 873 PetscInt j; 874 for (j=0;j<size_schur;j++) dummy_idx[j] = j; 875 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 876 } else { 877 if (mumps_S) { /* use MatMumps calls to invert S */ 878 #if defined(PETSC_HAVE_MUMPS) 879 ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr); 880 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 881 #endif 882 } else { /* we need to invert explicitly since we are not using MUMPS for S */ 883 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 884 ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr); 885 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 886 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */ 887 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr)); 888 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 889 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr)); 890 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 891 } else { 892 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr)); 893 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 894 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 895 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 896 } 897 ierr = PetscFPTrapPop();CHKERRQ(ierr); 898 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 899 } 900 /* S_Ej_tilda_all */ 901 cum = cum2 = 0; 902 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 903 for (i=0;i<sub_schurs->n_subs;i++) { 904 PetscInt j; 905 906 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 907 /* get (St^-1)_E */ 908 if (sub_schurs->is_hermitian) { /* Here I don't need to expand to upper triangular (column oriented) */ 909 PetscInt k; 910 for (k=0;k<subset_size;k++) { 911 for (j=k;j<subset_size;j++) { 912 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 913 } 914 } 915 } else { /* copy to workspace */ 916 PetscInt k; 917 for (k=0;k<subset_size;k++) { 918 for (j=0;j<subset_size;j++) { 919 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 920 } 921 } 922 } 923 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 924 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 925 cum += subset_size; 926 cum2 += subset_size*(size_schur + 1); 927 } 928 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 929 #if defined(PETSC_HAVE_MUMPS) 930 if (mumps_S) { 931 ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr); 932 } 933 #endif 934 } 935 } 936 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 937 } 938 ierr = PetscFree(nnz);CHKERRQ(ierr); 939 ierr = MatDestroy(&F);CHKERRQ(ierr); 940 ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr); 941 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 942 ierr = MatDestroy(&A_BB);CHKERRQ(ierr); 943 ierr = MatDestroy(&A_IB);CHKERRQ(ierr); 944 ierr = MatDestroy(&A_BI);CHKERRQ(ierr); 945 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 946 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 947 if (compute_Stilda) { 948 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 949 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 950 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 951 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 952 } 953 954 /* Global matrix of all assembled Schur on subsets */ 955 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 956 ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 957 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 958 959 /* Get local part of (\sum_j S_Ej) */ 960 if (!sub_schurs->sum_S_Ej_all) { 961 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr); 962 sub_schurs->sum_S_Ej_all = submats[0]; 963 } else { 964 ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr); 965 submats[0] = sub_schurs->sum_S_Ej_all; 966 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 967 } 968 969 /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */ 970 if (faster_deluxe) { 971 Mat tmpmat; 972 PetscScalar *array; 973 PetscInt cum; 974 975 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 976 cum = 0; 977 for (i=0;i<sub_schurs->n_subs;i++) { 978 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 979 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 980 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 981 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 982 PetscInt j,k; 983 984 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 985 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 986 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 987 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 988 for (j=0;j<B_N;j++) { 989 for (k=j+1;k<B_N;k++) { 990 array[k*B_N+j+cum] = array[j*B_N+k+cum]; 991 } 992 } 993 } else { 994 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 995 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 996 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 997 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 998 } 999 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1000 cum += subset_size*subset_size; 1001 } 1002 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 1003 ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr); 1004 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1005 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1006 sub_schurs->S_Ej_all = tmpmat; 1007 } 1008 1009 /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */ 1010 if (compute_Stilda) { 1011 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1012 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1013 submats[0] = sub_schurs->sum_S_Ej_tilda_all; 1014 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1015 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1016 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1017 submats[0] = sub_schurs->sum_S_Ej_inv_all; 1018 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1019 } 1020 1021 /* free workspace */ 1022 ierr = PetscFree(submats);CHKERRQ(ierr); 1023 ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 1024 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 1025 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 1026 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 1027 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 1028 PetscFunctionReturn(0); 1029 } 1030 1031 #undef __FUNCT__ 1032 #define __FUNCT__ "PCBDDCSubSchursInit" 1033 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap) 1034 { 1035 IS *faces,*edges,*all_cc,vertices; 1036 PetscInt i,n_faces,n_edges,n_all_cc; 1037 PetscBool is_sorted; 1038 PetscErrorCode ierr; 1039 1040 PetscFunctionBegin; 1041 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 1042 if (!is_sorted) { 1043 SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 1044 } 1045 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 1046 if (!is_sorted) { 1047 SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 1048 } 1049 1050 /* reset any previous data */ 1051 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 1052 1053 /* get index sets for faces and edges (already sorted by global ordering) */ 1054 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1055 n_all_cc = n_faces+n_edges; 1056 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 1057 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 1058 for (i=0;i<n_faces;i++) { 1059 all_cc[i] = faces[i]; 1060 } 1061 for (i=0;i<n_edges;i++) { 1062 all_cc[n_faces+i] = edges[i]; 1063 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 1064 } 1065 ierr = PetscFree(faces);CHKERRQ(ierr); 1066 ierr = PetscFree(edges);CHKERRQ(ierr); 1067 sub_schurs->is_dir = NULL; 1068 ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr); 1069 1070 /* Determine if MUMPS can be used */ 1071 sub_schurs->use_mumps = PETSC_FALSE; 1072 #if defined(PETSC_HAVE_MUMPS) 1073 sub_schurs->use_mumps = PETSC_TRUE; 1074 #endif 1075 1076 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 1077 sub_schurs->is_I = is_I; 1078 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 1079 sub_schurs->is_B = is_B; 1080 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 1081 sub_schurs->l2gmap = graph->l2gmap; 1082 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 1083 sub_schurs->BtoNmap = BtoNmap; 1084 sub_schurs->n_subs = n_all_cc; 1085 sub_schurs->is_subs = all_cc; 1086 if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */ 1087 for (i=0;i<sub_schurs->n_subs;i++) { 1088 ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 1089 } 1090 } 1091 sub_schurs->is_vertices = vertices; 1092 sub_schurs->S_Ej_all = NULL; 1093 sub_schurs->sum_S_Ej_all = NULL; 1094 sub_schurs->sum_S_Ej_inv_all = NULL; 1095 sub_schurs->sum_S_Ej_tilda_all = NULL; 1096 sub_schurs->is_Ej_all = NULL; 1097 PetscFunctionReturn(0); 1098 } 1099 1100 #undef __FUNCT__ 1101 #define __FUNCT__ "PCBDDCSubSchursCreate" 1102 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 1103 { 1104 PCBDDCSubSchurs schurs_ctx; 1105 PetscErrorCode ierr; 1106 1107 PetscFunctionBegin; 1108 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 1109 schurs_ctx->n_subs = 0; 1110 *sub_schurs = schurs_ctx; 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "PCBDDCSubSchursDestroy" 1116 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 1117 { 1118 PetscErrorCode ierr; 1119 1120 PetscFunctionBegin; 1121 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 1122 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 1123 PetscFunctionReturn(0); 1124 } 1125 1126 #undef __FUNCT__ 1127 #define __FUNCT__ "PCBDDCSubSchursReset" 1128 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 1129 { 1130 PetscInt i; 1131 PetscErrorCode ierr; 1132 1133 PetscFunctionBegin; 1134 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 1135 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 1136 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 1137 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 1138 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 1139 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 1140 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1141 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1142 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1143 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1144 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 1145 ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr); 1146 ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr); 1147 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 1148 for (i=0;i<sub_schurs->n_subs;i++) { 1149 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 1150 } 1151 if (sub_schurs->n_subs) { 1152 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 1153 } 1154 if (sub_schurs->reuse_mumps) { 1155 ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr); 1156 } 1157 ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr); 1158 sub_schurs->n_subs = 0; 1159 PetscFunctionReturn(0); 1160 } 1161 1162 #undef __FUNCT__ 1163 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 1164 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 1165 { 1166 PetscInt i,j,n; 1167 PetscErrorCode ierr; 1168 1169 PetscFunctionBegin; 1170 n = 0; 1171 for (i=-n_prev;i<0;i++) { 1172 PetscInt start_dof = queue_tip[i]; 1173 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 1174 PetscInt dof = adjncy[j]; 1175 if (!PetscBTLookup(touched,dof)) { 1176 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 1177 queue_tip[n] = dof; 1178 n++; 1179 } 1180 } 1181 } 1182 *n_added = n; 1183 PetscFunctionReturn(0); 1184 } 1185