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