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 /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */ 389 390 /* get Schur complement matrices */ 391 if (!sub_schurs->use_mumps) { 392 Mat tA_IB,tA_BI,tA_BB; 393 PetscBool isseqsbaij; 394 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr); 395 ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 396 if (isseqsbaij) { 397 ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr); 398 ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr); 399 ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr); 400 } else { 401 ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr); 402 A_BB = tA_BB; 403 ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr); 404 A_IB = tA_IB; 405 ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr); 406 A_BI = tA_BI; 407 } 408 } else { 409 A_II = NULL; 410 A_IB = NULL; 411 A_BI = NULL; 412 A_BB = NULL; 413 } 414 S_all = NULL; 415 416 /* determine interior problems */ 417 ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 418 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 419 PetscBT touched; 420 const PetscInt* idx_B; 421 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 422 423 if (xadj == NULL || adjncy == NULL) { 424 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 425 } 426 /* get sizes */ 427 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 428 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 429 430 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 431 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 432 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 433 434 /* all boundary dofs must be skipped when adding layers */ 435 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 436 for (j=0;j<n_B;j++) { 437 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 438 } 439 ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 440 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 441 442 /* add prescribed number of layers of dofs */ 443 n_local_dofs = n_B; 444 n_prev_added = n_B; 445 for (layer=0;layer<nlayers;layer++) { 446 PetscInt n_added; 447 if (n_local_dofs == n_I+n_B) break; 448 if (n_local_dofs > n_I+n_B) { 449 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); 450 } 451 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 452 n_prev_added = n_added; 453 n_local_dofs += n_added; 454 if (!n_added) break; 455 } 456 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 457 458 /* IS for I layer dofs in original numbering */ 459 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); 460 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 461 ierr = ISSort(is_I_layer);CHKERRQ(ierr); 462 /* IS for I layer dofs in I numbering */ 463 if (!sub_schurs->use_mumps) { 464 ISLocalToGlobalMapping ItoNmap; 465 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 466 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr); 467 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 468 469 /* II block */ 470 ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 471 } 472 } else { 473 PetscInt n_I; 474 475 /* IS for I dofs in original numbering */ 476 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 477 is_I_layer = sub_schurs->is_I; 478 479 /* IS for I dofs in I numbering (strided 1) */ 480 if (!sub_schurs->use_mumps) { 481 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 482 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 483 484 /* II block is the same */ 485 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 486 AE_II = A_II; 487 } 488 } 489 490 /* Get info on subset sizes and sum of all subsets sizes */ 491 max_subset_size = 0; 492 local_size = 0; 493 for (i=0;i<sub_schurs->n_subs;i++) { 494 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 495 max_subset_size = PetscMax(subset_size,max_subset_size); 496 local_size += subset_size; 497 } 498 499 /* Work arrays for local indices */ 500 extra = 0; 501 if (sub_schurs->use_mumps && is_I_layer) { 502 ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr); 503 } 504 ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 505 if (extra) { 506 const PetscInt *idxs; 507 ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr); 508 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 509 ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr); 510 } 511 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 512 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr); 513 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr); 514 515 /* Get local indices in local numbering */ 516 local_size = 0; 517 for (i=0;i<sub_schurs->n_subs;i++) { 518 PetscInt j; 519 const PetscInt *idxs; 520 521 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 522 ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 523 /* start (smallest in global ordering) and multiplicity */ 524 auxnum1[i] = idxs[0]; 525 auxnum2[i] = subset_size; 526 /* subset indices in local numbering */ 527 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 528 ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 529 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 530 local_size += subset_size; 531 } 532 533 /* allocate extra workspace needed only for GETRI */ 534 Bwork = NULL; 535 pivots = NULL; 536 if (local_size && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) { 537 PetscScalar lwork; 538 539 B_lwork = -1; 540 ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 541 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 542 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr)); 543 ierr = PetscFPTrapPop();CHKERRQ(ierr); 544 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 545 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 546 ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr); 547 } 548 549 /* prepare parallel matrices for summing up properly schurs on subsets */ 550 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr); 551 ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr); 552 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 553 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr); 554 ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr); 555 ierr = ISDestroy(&all_subsets);CHKERRQ(ierr); 556 ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr); 557 ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr); 558 if (i != local_size) { 559 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size); 560 } 561 ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr); 562 ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 563 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 564 ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 565 ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 566 ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 567 568 /* subset indices in local boundary numbering */ 569 if (!sub_schurs->is_Ej_all) { 570 PetscInt *all_local_idx_B; 571 572 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 573 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 574 if (subset_size != local_size) { 575 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 576 } 577 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 578 } 579 580 /* Local matrix of all local Schur on subsets (transposed) */ 581 if (!sub_schurs->S_Ej_all) { 582 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 583 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 584 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 585 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 586 } 587 588 /* Compute Schur complements explicitly */ 589 F = NULL; 590 if (!sub_schurs->use_mumps) { 591 Mat S_Ej_expl; 592 PetscScalar *work; 593 PetscInt j,*dummy_idx; 594 PetscBool Sdense; 595 596 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 597 local_size = 0; 598 for (i=0;i<sub_schurs->n_subs;i++) { 599 IS is_subset_B; 600 Mat AE_EE,AE_IE,AE_EI,S_Ej; 601 602 /* subsets in original and boundary numbering */ 603 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr); 604 /* EE block */ 605 ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr); 606 /* IE block */ 607 ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr); 608 /* EI block */ 609 if (sub_schurs->is_hermitian) { 610 ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 611 } else { 612 ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr); 613 } 614 ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr); 615 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr); 616 ierr = MatDestroy(&AE_EE);CHKERRQ(ierr); 617 ierr = MatDestroy(&AE_IE);CHKERRQ(ierr); 618 ierr = MatDestroy(&AE_EI);CHKERRQ(ierr); 619 if (AE_II == A_II) { /* we can reuse the same ksp */ 620 KSP ksp; 621 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 622 ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr); 623 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 624 KSP origksp,schurksp; 625 PC origpc,schurpc; 626 KSPType ksp_type; 627 PetscInt n_internal; 628 PetscBool ispcnone; 629 630 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 631 ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr); 632 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 633 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 634 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 635 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 636 ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr); 637 if (!ispcnone) { 638 PCType pc_type; 639 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 640 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 641 } else { 642 ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr); 643 } 644 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 645 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 646 MatSolverPackage solver=NULL; 647 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 648 if (solver) { 649 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 650 } 651 } 652 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 653 } 654 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 655 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 656 ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 657 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 658 if (Sdense) { 659 for (j=0;j<subset_size;j++) { 660 dummy_idx[j]=local_size+j; 661 } 662 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 663 } else { 664 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 665 } 666 ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 667 ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr); 668 local_size += subset_size; 669 } 670 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 671 /* free */ 672 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 673 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 674 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 675 } else { 676 Mat A; 677 IS is_A_all; 678 PetscScalar *work,*S_data; 679 PetscInt *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2; 680 PetscBool mumps_S; 681 682 /* get working mat */ 683 n_I = 0; 684 if (is_I_layer) { 685 ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr); 686 } 687 if (!sub_schurs->is_dir) { 688 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 689 size_schur = local_size; 690 } else { 691 IS list[2]; 692 693 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr); 694 list[1] = sub_schurs->is_dir; 695 ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr); 696 ierr = ISDestroy(&list[0]);CHKERRQ(ierr); 697 ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr); 698 size_schur += local_size; 699 } 700 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 701 size_active_schur = local_size; /* size active schurs does not count any dirichlet dof on the interface */ 702 ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 703 ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr); 704 ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr); 705 ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 706 ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 707 708 if (n_I) { 709 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 710 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 711 } else { 712 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 713 } 714 /* subsets ordered last */ 715 ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr); 716 for (i=0;i<size_schur;i++) { 717 idxs_schur[i] = n_I+i+1; 718 } 719 #if defined(PETSC_HAVE_MUMPS) 720 ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr); 721 #endif 722 ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 723 724 /* factorization step */ 725 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 726 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 727 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 728 ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr); 729 #endif 730 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 731 } else { 732 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 733 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 734 ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr); 735 #endif 736 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 737 } 738 739 /* get explicit Schur Complement computed during numeric factorization */ 740 #if defined(PETSC_HAVE_MUMPS) 741 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 742 #endif 743 744 /* we can reuse the solvers if we are not using the economic version */ 745 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr); 746 reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all)); 747 mumps_S = PETSC_TRUE; 748 } else { /* we can't use MUMPS when size_schur == size_of_the_problem */ 749 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 750 reuse_solvers = PETSC_FALSE; 751 mumps_S = PETSC_FALSE; 752 } 753 754 if (reuse_solvers) { 755 Mat A_II; 756 Vec vec1_B; 757 PCBDDCReuseMumps msolv_ctx; 758 759 if (sub_schurs->reuse_mumps) { 760 ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr); 761 } else { 762 ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr); 763 } 764 msolv_ctx = sub_schurs->reuse_mumps; 765 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 766 ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr); 767 ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 768 msolv_ctx->F = F; 769 ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr); 770 771 /* interior solver */ 772 ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr); 773 ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr); 774 ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr); 775 ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr); 776 ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr); 777 ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr); 778 779 /* correction solver */ 780 ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr); 781 ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr); 782 ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr); 783 ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr); 784 ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr); 785 ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr); 786 787 /* scatter and vecs for Schur complement solver */ 788 ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr); 789 ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr); 790 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr); 791 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 792 ierr = VecDestroy(&vec1_B);CHKERRQ(ierr); 793 ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr); 794 msolv_ctx->is_R = is_A_all; 795 } 796 ierr = MatDestroy(&A);CHKERRQ(ierr); 797 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 798 799 /* Work arrays */ 800 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 801 802 /* matrices for adaptive selection */ 803 if (compute_Stilda) { 804 if (!sub_schurs->sum_S_Ej_tilda_all) { 805 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 806 ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr); 807 ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 808 ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 809 } 810 if (!sub_schurs->sum_S_Ej_inv_all) { 811 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 812 ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr); 813 ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr); 814 ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr); 815 } 816 } 817 818 /* S_Ej_all */ 819 cum = cum2 = 0; 820 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 821 for (i=0;i<sub_schurs->n_subs;i++) { 822 PetscInt j; 823 824 /* get S_E */ 825 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 826 if (mumps_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */ 827 PetscInt k; 828 for (k=0;k<subset_size;k++) { 829 for (j=k;j<subset_size;j++) { 830 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 831 work[j*subset_size+k] = S_data[cum2+k*size_schur+j]; 832 } 833 } 834 } else { /* copy to workspace */ 835 PetscInt k; 836 for (k=0;k<subset_size;k++) { 837 for (j=0;j<subset_size;j++) { 838 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 839 } 840 } 841 } 842 /* insert S_E values */ 843 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 844 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 845 846 /* if adaptivity is requested, invert S_E block */ 847 if (compute_Stilda) { 848 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 849 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 850 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */ 851 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr)); 852 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 853 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr)); 854 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 855 } else { 856 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr)); 857 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 858 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 859 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 860 } 861 ierr = PetscFPTrapPop();CHKERRQ(ierr); 862 ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 863 } 864 cum += subset_size; 865 cum2 += subset_size*(size_schur + 1); 866 } 867 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 868 869 #if defined(PETSC_HAVE_MUMPS) 870 if (mumps_S) { 871 ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr); 872 } 873 #endif 874 875 if (compute_Stilda && size_active_schur) { 876 if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */ 877 PetscInt j; 878 for (j=0;j<size_schur;j++) dummy_idx[j] = j; 879 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 880 } else { 881 if (mumps_S) { /* use MatMumps calls to invert S */ 882 #if defined(PETSC_HAVE_MUMPS) 883 ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr); 884 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 885 #endif 886 } else { /* we need to invert explicitly since we are not using MUMPS for S */ 887 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 888 ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr); 889 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 890 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */ 891 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr)); 892 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 893 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr)); 894 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 895 } else { 896 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr)); 897 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 898 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 899 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 900 } 901 ierr = PetscFPTrapPop();CHKERRQ(ierr); 902 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 903 } 904 /* S_Ej_tilda_all */ 905 cum = cum2 = 0; 906 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 907 for (i=0;i<sub_schurs->n_subs;i++) { 908 PetscInt j; 909 910 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 911 /* get (St^-1)_E */ 912 if (sub_schurs->is_hermitian) { /* Here I don't need to expand to upper triangular (column oriented) */ 913 PetscInt k; 914 for (k=0;k<subset_size;k++) { 915 for (j=k;j<subset_size;j++) { 916 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 917 } 918 } 919 } else { /* copy to workspace */ 920 PetscInt k; 921 for (k=0;k<subset_size;k++) { 922 for (j=0;j<subset_size;j++) { 923 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 924 } 925 } 926 } 927 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 928 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 929 cum += subset_size; 930 cum2 += subset_size*(size_schur + 1); 931 } 932 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 933 #if defined(PETSC_HAVE_MUMPS) 934 if (mumps_S) { 935 ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr); 936 } 937 #endif 938 } 939 } 940 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 941 } 942 ierr = PetscFree(nnz);CHKERRQ(ierr); 943 ierr = MatDestroy(&F);CHKERRQ(ierr); 944 ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr); 945 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 946 ierr = MatDestroy(&A_BB);CHKERRQ(ierr); 947 ierr = MatDestroy(&A_IB);CHKERRQ(ierr); 948 ierr = MatDestroy(&A_BI);CHKERRQ(ierr); 949 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 950 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 951 if (compute_Stilda) { 952 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 953 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 954 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 955 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 956 } 957 958 /* Global matrix of all assembled Schur on subsets */ 959 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 960 ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 961 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 962 963 /* Get local part of (\sum_j S_Ej) */ 964 if (!sub_schurs->sum_S_Ej_all) { 965 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr); 966 sub_schurs->sum_S_Ej_all = submats[0]; 967 } else { 968 ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr); 969 submats[0] = sub_schurs->sum_S_Ej_all; 970 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 971 } 972 973 /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */ 974 if (faster_deluxe) { 975 Mat tmpmat; 976 PetscScalar *array; 977 PetscInt cum; 978 979 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 980 cum = 0; 981 for (i=0;i<sub_schurs->n_subs;i++) { 982 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 983 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 984 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 985 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 986 PetscInt j,k; 987 988 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 989 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 990 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 991 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 992 for (j=0;j<B_N;j++) { 993 for (k=j+1;k<B_N;k++) { 994 array[k*B_N+j+cum] = array[j*B_N+k+cum]; 995 } 996 } 997 } else { 998 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 999 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 1000 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1001 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 1002 } 1003 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1004 cum += subset_size*subset_size; 1005 } 1006 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 1007 ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr); 1008 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1009 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1010 sub_schurs->S_Ej_all = tmpmat; 1011 } 1012 1013 /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */ 1014 if (compute_Stilda) { 1015 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1016 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1017 submats[0] = sub_schurs->sum_S_Ej_tilda_all; 1018 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1019 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1020 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1021 submats[0] = sub_schurs->sum_S_Ej_inv_all; 1022 ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1023 } 1024 1025 /* free workspace */ 1026 ierr = PetscFree(submats);CHKERRQ(ierr); 1027 ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 1028 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 1029 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 1030 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 1031 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 1032 PetscFunctionReturn(0); 1033 } 1034 1035 #undef __FUNCT__ 1036 #define __FUNCT__ "PCBDDCSubSchursInit" 1037 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap) 1038 { 1039 IS *faces,*edges,*all_cc,vertices; 1040 PetscInt i,n_faces,n_edges,n_all_cc; 1041 PetscBool is_sorted; 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 1046 if (!is_sorted) { 1047 SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 1048 } 1049 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 1050 if (!is_sorted) { 1051 SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 1052 } 1053 1054 /* reset any previous data */ 1055 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 1056 1057 /* get index sets for faces and edges (already sorted by global ordering) */ 1058 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1059 n_all_cc = n_faces+n_edges; 1060 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 1061 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 1062 for (i=0;i<n_faces;i++) { 1063 all_cc[i] = faces[i]; 1064 } 1065 for (i=0;i<n_edges;i++) { 1066 all_cc[n_faces+i] = edges[i]; 1067 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 1068 } 1069 ierr = PetscFree(faces);CHKERRQ(ierr); 1070 ierr = PetscFree(edges);CHKERRQ(ierr); 1071 sub_schurs->is_dir = NULL; 1072 ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr); 1073 1074 /* Determine if MUMPS can be used */ 1075 sub_schurs->use_mumps = PETSC_FALSE; 1076 #if defined(PETSC_HAVE_MUMPS) 1077 sub_schurs->use_mumps = PETSC_TRUE; 1078 #endif 1079 1080 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 1081 sub_schurs->is_I = is_I; 1082 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 1083 sub_schurs->is_B = is_B; 1084 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 1085 sub_schurs->l2gmap = graph->l2gmap; 1086 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 1087 sub_schurs->BtoNmap = BtoNmap; 1088 sub_schurs->n_subs = n_all_cc; 1089 sub_schurs->is_subs = all_cc; 1090 if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */ 1091 for (i=0;i<sub_schurs->n_subs;i++) { 1092 ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 1093 } 1094 } 1095 sub_schurs->is_vertices = vertices; 1096 sub_schurs->S_Ej_all = NULL; 1097 sub_schurs->sum_S_Ej_all = NULL; 1098 sub_schurs->sum_S_Ej_inv_all = NULL; 1099 sub_schurs->sum_S_Ej_tilda_all = NULL; 1100 sub_schurs->is_Ej_all = NULL; 1101 PetscFunctionReturn(0); 1102 } 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "PCBDDCSubSchursCreate" 1106 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 1107 { 1108 PCBDDCSubSchurs schurs_ctx; 1109 PetscErrorCode ierr; 1110 1111 PetscFunctionBegin; 1112 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 1113 schurs_ctx->n_subs = 0; 1114 *sub_schurs = schurs_ctx; 1115 PetscFunctionReturn(0); 1116 } 1117 1118 #undef __FUNCT__ 1119 #define __FUNCT__ "PCBDDCSubSchursDestroy" 1120 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 1121 { 1122 PetscErrorCode ierr; 1123 1124 PetscFunctionBegin; 1125 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 1126 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 1127 PetscFunctionReturn(0); 1128 } 1129 1130 #undef __FUNCT__ 1131 #define __FUNCT__ "PCBDDCSubSchursReset" 1132 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 1133 { 1134 PetscInt i; 1135 PetscErrorCode ierr; 1136 1137 PetscFunctionBegin; 1138 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 1139 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 1140 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 1141 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 1142 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 1143 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 1144 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1145 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1146 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1147 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1148 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 1149 ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr); 1150 ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr); 1151 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 1152 for (i=0;i<sub_schurs->n_subs;i++) { 1153 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 1154 } 1155 if (sub_schurs->n_subs) { 1156 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 1157 } 1158 if (sub_schurs->reuse_mumps) { 1159 ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr); 1160 } 1161 ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr); 1162 sub_schurs->n_subs = 0; 1163 PetscFunctionReturn(0); 1164 } 1165 1166 #undef __FUNCT__ 1167 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 1168 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 1169 { 1170 PetscInt i,j,n; 1171 PetscErrorCode ierr; 1172 1173 PetscFunctionBegin; 1174 n = 0; 1175 for (i=-n_prev;i<0;i++) { 1176 PetscInt start_dof = queue_tip[i]; 1177 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 1178 PetscInt dof = adjncy[j]; 1179 if (!PetscBTLookup(touched,dof)) { 1180 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 1181 queue_tip[n] = dof; 1182 n++; 1183 } 1184 } 1185 } 1186 *n_added = n; 1187 PetscFunctionReturn(0); 1188 } 1189