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