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