1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <../src/mat/impls/dense/seq/dense.h> 4 #include <petscblaslapack.h> 5 6 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 7 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*); 8 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec); 9 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec); 10 11 /* if v2 is not present, correction is done in-place */ 12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full) 13 { 14 PetscScalar *array; 15 PetscScalar *array2; 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 if (!ctx->benign_n) PetscFunctionReturn(0); 20 if (sol && full) { 21 PetscInt n_I,size_schur; 22 23 /* get sizes */ 24 ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr); 25 ierr = VecGetSize(v,&n_I);CHKERRQ(ierr); 26 n_I = n_I - size_schur; 27 /* get schur sol from array */ 28 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 29 ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr); 30 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 31 /* apply interior sol correction */ 32 ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr); 33 ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 34 ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr); 35 } 36 if (v2) { 37 PetscInt nl; 38 39 ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 40 ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr); 41 ierr = VecGetArray(v2,&array2);CHKERRQ(ierr); 42 ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr); 43 } else { 44 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 45 array2 = array; 46 } 47 if (!sol) { /* change rhs */ 48 PetscInt n; 49 for (n=0;n<ctx->benign_n;n++) { 50 PetscScalar sum = 0.; 51 const PetscInt *cols; 52 PetscInt nz,i; 53 54 ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr); 55 ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 56 for (i=0;i<nz-1;i++) sum += array[cols[i]]; 57 #if defined(PETSC_USE_COMPLEX) 58 sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz)); 59 #else 60 sum = -sum/nz; 61 #endif 62 for (i=0;i<nz-1;i++) array2[cols[i]] += sum; 63 ctx->benign_save_vals[n] = array2[cols[nz-1]]; 64 array2[cols[nz-1]] = sum; 65 ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 66 } 67 } else { 68 PetscInt n; 69 for (n=0;n<ctx->benign_n;n++) { 70 PetscScalar sum = 0.; 71 const PetscInt *cols; 72 PetscInt nz,i; 73 ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr); 74 ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 75 for (i=0;i<nz-1;i++) sum += array[cols[i]]; 76 #if defined(PETSC_USE_COMPLEX) 77 sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz)); 78 #else 79 sum = -sum/nz; 80 #endif 81 for (i=0;i<nz-1;i++) array2[cols[i]] += sum; 82 array2[cols[nz-1]] = ctx->benign_save_vals[n]; 83 ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 84 } 85 } 86 if (v2) { 87 ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 88 ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr); 89 } else { 90 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 91 } 92 if (!sol && full) { 93 Vec usedv; 94 PetscInt n_I,size_schur; 95 96 /* get sizes */ 97 ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr); 98 ierr = VecGetSize(v,&n_I);CHKERRQ(ierr); 99 n_I = n_I - size_schur; 100 /* compute schur rhs correction */ 101 if (v2) { 102 usedv = v2; 103 } else { 104 usedv = v; 105 } 106 /* apply schur rhs correction */ 107 ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr); 108 ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr); 109 ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr); 110 ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr); 111 ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 112 ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 113 } 114 PetscFunctionReturn(0); 115 } 116 117 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full) 118 { 119 PCBDDCReuseSolvers ctx; 120 PetscBool copy = PETSC_FALSE; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr); 125 if (full) { 126 #if defined(PETSC_HAVE_MUMPS) 127 ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr); 128 #endif 129 #if defined(PETSC_HAVE_MKL_PARDISO) 130 ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr); 131 #endif 132 copy = ctx->has_vertices; 133 } else { /* interior solver */ 134 #if defined(PETSC_HAVE_MUMPS) 135 ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr); 136 #endif 137 #if defined(PETSC_HAVE_MKL_PARDISO) 138 ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr); 139 #endif 140 copy = PETSC_TRUE; 141 } 142 /* copy rhs into factored matrix workspace */ 143 if (copy) { 144 PetscInt n; 145 PetscScalar *array,*array_solver; 146 147 ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr); 148 ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 149 ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr); 150 ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr); 151 ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr); 152 ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 153 154 ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr); 155 if (transpose) { 156 ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 157 } else { 158 ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 159 } 160 ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr); 161 162 /* get back data to caller worskpace */ 163 ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr); 164 ierr = VecGetArray(sol,&array);CHKERRQ(ierr); 165 ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr); 166 ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr); 167 ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr); 168 } else { 169 if (ctx->benign_n) { 170 ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr); 171 if (transpose) { 172 ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr); 173 } else { 174 ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr); 175 } 176 ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr); 177 } else { 178 if (transpose) { 179 ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr); 180 } else { 181 ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr); 182 } 183 } 184 } 185 /* restore defaults */ 186 #if defined(PETSC_HAVE_MUMPS) 187 ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr); 188 #endif 189 #if defined(PETSC_HAVE_MKL_PARDISO) 190 ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr); 191 #endif 192 PetscFunctionReturn(0); 193 } 194 195 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol) 196 { 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol) 205 { 206 PetscErrorCode ierr; 207 208 PetscFunctionBegin; 209 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr); 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 222 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol) 223 { 224 PetscErrorCode ierr; 225 226 PetscFunctionBegin; 227 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse) 232 { 233 PetscInt i; 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 ierr = MatDestroy(&reuse->F);CHKERRQ(ierr); 238 ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr); 239 ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr); 240 ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr); 241 ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr); 242 ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr); 243 ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr); 244 ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr); 245 ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr); 246 ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr); 247 for (i=0;i<reuse->benign_n;i++) { 248 ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr); 249 } 250 ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr); 251 ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr); 252 ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr); 253 ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr); 254 ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr); 255 ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr); 256 PetscFunctionReturn(0); 257 } 258 259 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 260 { 261 Mat B, C, D, Bd, Cd, AinvBd; 262 KSP ksp; 263 PC pc; 264 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 265 PetscReal fill = 2.0; 266 PetscInt n_I; 267 PetscMPIInt size; 268 PetscErrorCode ierr; 269 270 PetscFunctionBegin; 271 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 272 if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 273 if (reuse == MAT_REUSE_MATRIX) { 274 PetscBool Sdense; 275 276 ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr); 277 if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense"); 278 } 279 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 280 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 281 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 282 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 283 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 284 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 285 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 286 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 287 ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr); 288 if (n_I) { 289 if (!Bdense) { 290 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 291 } else { 292 Bd = B; 293 } 294 295 if (isLU || isILU || isCHOL) { 296 Mat fact; 297 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 298 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 299 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 300 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 301 } else { 302 PetscBool ex = PETSC_TRUE; 303 304 if (ex) { 305 Mat Ainvd; 306 307 ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 308 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 309 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 310 } else { 311 Vec sol,rhs; 312 PetscScalar *arrayrhs,*arraysol; 313 PetscInt i,nrhs,n; 314 315 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 316 ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr); 317 ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr); 318 ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr); 319 ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr); 320 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 321 for (i=0;i<nrhs;i++) { 322 ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr); 323 ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr); 324 ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr); 325 ierr = VecResetArray(rhs);CHKERRQ(ierr); 326 ierr = VecResetArray(sol);CHKERRQ(ierr); 327 } 328 ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr); 329 ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr); 330 } 331 } 332 if (!Bdense & !issym) { 333 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 334 } 335 336 if (!issym) { 337 if (!Cdense) { 338 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 339 } else { 340 Cd = C; 341 } 342 ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 343 if (!Cdense) { 344 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 345 } 346 } else { 347 ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 348 if (!Bdense) { 349 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 350 } 351 } 352 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 353 } 354 355 if (D) { 356 Mat Dd; 357 PetscBool Ddense; 358 359 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 360 if (!Ddense) { 361 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 362 } else { 363 Dd = D; 364 } 365 if (n_I) { 366 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 367 } else { 368 if (reuse == MAT_INITIAL_MATRIX) { 369 ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr); 370 } else { 371 ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 372 } 373 } 374 if (!Ddense) { 375 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 376 } 377 } else { 378 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 379 } 380 PetscFunctionReturn(0); 381 } 382 383 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal) 384 { 385 Mat F,A_II,A_IB,A_BI,A_BB,AE_II; 386 Mat S_all; 387 Mat global_schur_subsets,work_mat,*submats; 388 ISLocalToGlobalMapping l2gmap_subsets; 389 IS is_I,is_I_layer; 390 IS all_subsets,all_subsets_mult,all_subsets_n; 391 PetscInt *nnz,*all_local_idx_N; 392 PetscInt *auxnum1,*auxnum2; 393 PetscInt i,subset_size,max_subset_size; 394 PetscInt n_B,extra,local_size,global_size; 395 PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 396 PetscScalar *Bwork; 397 PetscSubcomm subcomm; 398 PetscMPIInt color,rank; 399 MPI_Comm comm_n; 400 PetscBool deluxe = PETSC_TRUE, use_getr = PETSC_FALSE; 401 PetscErrorCode ierr; 402 403 PetscFunctionBegin; 404 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 405 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 406 407 sub_schurs->is_hermitian = PETSC_TRUE; 408 sub_schurs->is_posdef = PETSC_TRUE; 409 if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE; 410 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr); 411 ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr); 412 ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr); 413 ierr = PetscOptionsEnd();CHKERRQ(ierr); 414 415 /* convert matrix if needed */ 416 if (Ain) { 417 PetscBool isseqaij; 418 ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 419 if (isseqaij) { 420 ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr); 421 sub_schurs->A = Ain; 422 } else { 423 ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr); 424 } 425 } 426 427 ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr); 428 sub_schurs->S = Sin; 429 if (sub_schurs->schur_explicit) { 430 sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A); 431 } 432 433 /* preliminary checks */ 434 if (!sub_schurs->schur_explicit && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO"); 435 436 /* restrict work on active processes */ 437 color = 0; 438 if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */ 439 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr); 440 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr); 441 ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr); 442 ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr); 443 ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr); 444 ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr); 445 if (!sub_schurs->n_subs) { 446 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 447 PetscFunctionReturn(0); 448 } 449 /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */ 450 451 /* get Schur complement matrices */ 452 if (!sub_schurs->schur_explicit) { 453 Mat tA_IB,tA_BI,tA_BB; 454 PetscBool isseqsbaij; 455 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr); 456 ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 457 if (isseqsbaij) { 458 ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr); 459 ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr); 460 ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr); 461 } else { 462 ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr); 463 A_BB = tA_BB; 464 ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr); 465 A_IB = tA_IB; 466 ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr); 467 A_BI = tA_BI; 468 } 469 } else { 470 A_II = NULL; 471 A_IB = NULL; 472 A_BI = NULL; 473 A_BB = NULL; 474 } 475 S_all = NULL; 476 477 /* determine interior problems */ 478 ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 479 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 480 PetscBT touched; 481 const PetscInt* idx_B; 482 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 483 484 if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 485 /* get sizes */ 486 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 487 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 488 489 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 490 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 491 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 492 493 /* all boundary dofs must be skipped when adding layers */ 494 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 495 for (j=0;j<n_B;j++) { 496 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 497 } 498 ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 499 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 500 501 /* add prescribed number of layers of dofs */ 502 n_local_dofs = n_B; 503 n_prev_added = n_B; 504 for (layer=0;layer<nlayers;layer++) { 505 PetscInt n_added; 506 if (n_local_dofs == n_I+n_B) break; 507 if (n_local_dofs > n_I+n_B) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error querying layer %D. Out of bound access (%D > %D)",layer,n_local_dofs,n_I+n_B); 508 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 509 n_prev_added = n_added; 510 n_local_dofs += n_added; 511 if (!n_added) break; 512 } 513 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 514 515 /* IS for I layer dofs in original numbering */ 516 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); 517 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 518 ierr = ISSort(is_I_layer);CHKERRQ(ierr); 519 /* IS for I layer dofs in I numbering */ 520 if (!sub_schurs->schur_explicit) { 521 ISLocalToGlobalMapping ItoNmap; 522 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 523 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr); 524 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 525 526 /* II block */ 527 ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 528 } 529 } else { 530 PetscInt n_I; 531 532 /* IS for I dofs in original numbering */ 533 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 534 is_I_layer = sub_schurs->is_I; 535 536 /* IS for I dofs in I numbering (strided 1) */ 537 if (!sub_schurs->schur_explicit) { 538 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 539 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 540 541 /* II block is the same */ 542 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 543 AE_II = A_II; 544 } 545 } 546 547 /* Get info on subset sizes and sum of all subsets sizes */ 548 max_subset_size = 0; 549 local_size = 0; 550 for (i=0;i<sub_schurs->n_subs;i++) { 551 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 552 max_subset_size = PetscMax(subset_size,max_subset_size); 553 local_size += subset_size; 554 } 555 556 /* Work arrays for local indices */ 557 extra = 0; 558 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 559 if (sub_schurs->schur_explicit && is_I_layer) { 560 ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr); 561 } 562 ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr); 563 if (extra) { 564 const PetscInt *idxs; 565 ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr); 566 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 567 ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr); 568 } 569 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 570 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr); 571 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr); 572 573 /* Get local indices in local numbering */ 574 local_size = 0; 575 for (i=0;i<sub_schurs->n_subs;i++) { 576 PetscInt j; 577 const PetscInt *idxs; 578 579 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 580 ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 581 /* start (smallest in global ordering) and multiplicity */ 582 auxnum1[i] = idxs[0]; 583 auxnum2[i] = subset_size; 584 /* subset indices in local numbering */ 585 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 586 ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 587 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 588 local_size += subset_size; 589 } 590 591 /* allocate extra workspace needed only for GETRI */ 592 if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) { 593 PetscScalar lwork,dummyscalar = 0.; 594 PetscBLASInt dummyint = 0; 595 596 use_getr = PETSC_TRUE; 597 B_lwork = -1; 598 ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 599 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 600 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr)); 601 ierr = PetscFPTrapPop();CHKERRQ(ierr); 602 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 603 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 604 ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr); 605 } else { 606 Bwork = NULL; 607 pivots = NULL; 608 } 609 610 /* prepare parallel matrices for summing up properly schurs on subsets */ 611 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr); 612 ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr); 613 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 614 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr); 615 ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr); 616 ierr = ISDestroy(&all_subsets);CHKERRQ(ierr); 617 ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr); 618 ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr); 619 if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size); 620 ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr); 621 ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr); 622 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 623 ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 624 ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 625 ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 626 627 /* subset indices in local boundary numbering */ 628 if (!sub_schurs->is_Ej_all) { 629 PetscInt *all_local_idx_B; 630 631 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 632 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 633 if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D\n",subset_size,local_size); 634 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 635 } 636 637 if (change) { 638 ISLocalToGlobalMapping BtoS; 639 IS change_primal_B; 640 IS change_primal_all; 641 642 if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 643 if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 644 ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr); 645 for (i=0;i<sub_schurs->n_subs;i++) { 646 ISLocalToGlobalMapping NtoS; 647 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr); 648 ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr); 649 ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr); 650 } 651 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr); 652 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr); 653 ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr); 654 ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr); 655 ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr); 656 ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr); 657 for (i=0;i<sub_schurs->n_subs;i++) { 658 Mat change_sub; 659 660 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 661 ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr); 662 ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr); 663 if (!sub_schurs->change_with_qr) { 664 ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr); 665 } else { 666 Mat change_subt; 667 ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr); 668 ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr); 669 ierr = MatDestroy(&change_subt);CHKERRQ(ierr); 670 } 671 ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr); 672 ierr = MatDestroy(&change_sub);CHKERRQ(ierr); 673 ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr); 674 ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr); 675 } 676 ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr); 677 } 678 679 /* Local matrix of all local Schur on subsets (transposed) */ 680 if (!sub_schurs->S_Ej_all) { 681 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 682 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 683 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 684 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 685 } 686 687 /* Compute Schur complements explicitly */ 688 F = NULL; 689 if (!sub_schurs->schur_explicit) { /* this code branch is used when MatFactor with Schur complemnent support is not present; it is not very efficient, unless the economic version of the scaling is required */ 690 Mat S_Ej_expl; 691 PetscScalar *work; 692 PetscInt j,*dummy_idx; 693 PetscBool Sdense; 694 695 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 696 local_size = 0; 697 for (i=0;i<sub_schurs->n_subs;i++) { 698 IS is_subset_B; 699 Mat AE_EE,AE_IE,AE_EI,S_Ej; 700 701 /* subsets in original and boundary numbering */ 702 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr); 703 /* EE block */ 704 ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr); 705 /* IE block */ 706 ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr); 707 /* EI block */ 708 if (sub_schurs->is_hermitian) { 709 ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 710 } else { 711 ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr); 712 } 713 ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr); 714 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr); 715 ierr = MatDestroy(&AE_EE);CHKERRQ(ierr); 716 ierr = MatDestroy(&AE_IE);CHKERRQ(ierr); 717 ierr = MatDestroy(&AE_EI);CHKERRQ(ierr); 718 if (AE_II == A_II) { /* we can reuse the same ksp */ 719 KSP ksp; 720 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 721 ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr); 722 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 723 KSP origksp,schurksp; 724 PC origpc,schurpc; 725 KSPType ksp_type; 726 PetscInt n_internal; 727 PetscBool ispcnone; 728 729 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 730 ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr); 731 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 732 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 733 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 734 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 735 ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr); 736 if (!ispcnone) { 737 PCType pc_type; 738 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 739 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 740 } else { 741 ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr); 742 } 743 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 744 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 745 MatSolverPackage solver=NULL; 746 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 747 if (solver) { 748 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 749 } 750 } 751 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 752 } 753 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 754 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 755 ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 756 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 757 if (Sdense) { 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 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 763 ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 764 ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr); 765 local_size += subset_size; 766 } 767 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 768 /* free */ 769 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 770 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 771 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 772 } else { 773 Mat A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL; 774 Vec Dall = NULL; 775 IS is_A_all,*is_p_r = NULL; 776 PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL; 777 PetscInt n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2; 778 PetscBool economic,solver_S,S_lower_triangular = PETSC_FALSE; 779 PetscBool schur_has_vertices,factor_workaround; 780 char solver[256]; 781 782 /* get sizes */ 783 n_I = 0; 784 if (is_I_layer) { 785 ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr); 786 } 787 economic = PETSC_FALSE; 788 ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr); 789 if (cum != n_I) economic = PETSC_TRUE; 790 ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr); 791 size_active_schur = local_size; 792 793 /* import scaling vector (wrong formulation if we have 3D edges) */ 794 if (scaling && compute_Stilda) { 795 const PetscScalar *array; 796 PetscScalar *array2; 797 const PetscInt *idxs; 798 PetscInt i; 799 800 ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 801 ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr); 802 ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr); 803 ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr); 804 for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]]; 805 ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr); 806 ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr); 807 ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 808 deluxe = PETSC_FALSE; 809 } 810 811 /* size active schurs does not count any dirichlet or vertex dof on the interface */ 812 factor_workaround = PETSC_FALSE; 813 schur_has_vertices = PETSC_FALSE; 814 cum = n_I+size_active_schur; 815 if (sub_schurs->is_dir) { 816 const PetscInt* idxs; 817 PetscInt n_dir; 818 819 ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr); 820 ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 821 ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr); 822 ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 823 cum += n_dir; 824 factor_workaround = PETSC_TRUE; 825 } 826 /* include the primal vertices in the Schur complement */ 827 if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) { 828 PetscInt n_v; 829 830 ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr); 831 if (n_v) { 832 const PetscInt* idxs; 833 834 ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 835 ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr); 836 ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 837 cum += n_v; 838 factor_workaround = PETSC_TRUE; 839 schur_has_vertices = PETSC_TRUE; 840 } 841 } 842 size_schur = cum - n_I; 843 ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr); 844 /* get working mat (TODO: factorize without actually permuting) */ 845 if (cum == n) { 846 ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr); 847 ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr); 848 } else { 849 ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 850 } 851 ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr); 852 ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr); 853 854 /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory 855 n^2 instead of n^1.5 or something. This is a workaround */ 856 if (benign_n) { 857 Vec v; 858 ISLocalToGlobalMapping N_to_reor; 859 IS is_p0,is_p0_p; 860 PetscScalar *cs_AIB,*AIIm1_data; 861 PetscInt sizeA; 862 863 ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr); 864 ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr); 865 ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr); 866 ierr = ISDestroy(&is_p0);CHKERRQ(ierr); 867 ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr); 868 ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr); 869 ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr); 870 ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr); 871 ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 872 ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 873 ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr); 874 /* compute colsum of A_IB restricted to pressures */ 875 for (i=0;i<benign_n;i++) { 876 Vec benign_AIIm1_ones; 877 PetscScalar *array; 878 const PetscInt *idxs; 879 PetscInt j,nz; 880 881 ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr); 882 ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr); 883 ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 884 for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.; 885 ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 886 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr); 887 ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr); 888 ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr); 889 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 890 for (j=0;j<size_schur;j++) { 891 #if defined(PETSC_USE_COMPLEX) 892 cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz)); 893 #else 894 cs_AIB[i*size_schur + j] = array[j+n_I]/nz; 895 #endif 896 } 897 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 898 } 899 ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 900 ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 901 ierr = VecDestroy(&v);CHKERRQ(ierr); 902 ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr); 903 ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 904 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 905 ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr); 906 ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr); 907 ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr); 908 } 909 ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr); 910 ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 911 ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 912 #if defined(PETSC_HAVE_MUMPS) 913 ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr); 914 #else 915 ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr); 916 #endif 917 ierr = PetscOptionsGetString(NULL,((PetscObject)A)->prefix,"-mat_solver_package",solver,256,NULL);CHKERRQ(ierr); 918 919 /* when using the benign subspace trick, the local Schur complements are SPD */ 920 if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE; 921 922 if (n_I) { /* TODO add ordering from options */ 923 IS is_schur; 924 925 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 926 ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 927 } else { 928 ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 929 } 930 ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr); 931 932 /* subsets ordered last */ 933 ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr); 934 ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr); 935 ierr = ISDestroy(&is_schur);CHKERRQ(ierr); 936 937 /* factorization step */ 938 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 939 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 940 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 941 ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr); 942 #endif 943 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 944 S_lower_triangular = PETSC_TRUE; 945 } else { 946 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 947 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 948 ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr); 949 #endif 950 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 951 } 952 ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr); 953 954 /* get explicit Schur Complement computed during numeric factorization */ 955 ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr); 956 ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 957 ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 958 959 /* we can reuse the solvers if we are not using the economic version */ 960 reuse_solvers = (PetscBool)(reuse_solvers && !economic); 961 factor_workaround = (PetscBool)(reuse_solvers && factor_workaround); 962 solver_S = PETSC_TRUE; 963 964 /* update the Schur complement with the change of basis on the pressures */ 965 if (benign_n) { 966 PetscScalar *S_data,*cs_AIB,*AIIm1_data; 967 Vec v; 968 PetscInt sizeA; 969 970 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 971 ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr); 972 ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr); 973 #if defined(PETSC_HAVE_MUMPS) 974 ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr); 975 #endif 976 #if defined(PETSC_HAVE_MKL_PARDISO) 977 ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr); 978 #endif 979 ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 980 ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 981 for (i=0;i<benign_n;i++) { 982 Vec benign_AIIm1_ones; 983 PetscScalar *array,sum = 0.,one = 1.; 984 const PetscInt *idxs; 985 PetscInt j,nz; 986 PetscBLASInt B_k,B_n; 987 988 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr); 989 ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr); 990 ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr); 991 ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr); 992 ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 993 /* p0 dof (eliminated) is excluded from the sum */ 994 for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i]; 995 ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 996 ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr); 997 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 998 /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */ 999 /* cs_AIB already scaled by 1./nz */ 1000 B_k = 1; 1001 ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr); 1002 PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n)); 1003 sum = 1.; 1004 PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n)); 1005 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 1006 /* set p0 entry of AIIm1_ones to zero */ 1007 ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 1008 for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.; 1009 ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 1010 ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr); 1011 } 1012 /* restore defaults */ 1013 #if defined(PETSC_HAVE_MUMPS) 1014 ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr); 1015 #endif 1016 #if defined(PETSC_HAVE_MKL_PARDISO) 1017 ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr); 1018 #endif 1019 ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 1020 ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 1021 ierr = VecDestroy(&v);CHKERRQ(ierr); 1022 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 1023 } 1024 if (!reuse_solvers) { 1025 for (i=0;i<benign_n;i++) { 1026 ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr); 1027 } 1028 ierr = PetscFree(is_p_r);CHKERRQ(ierr); 1029 ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr); 1030 ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr); 1031 } 1032 } else { /* we can't use MatFactor when size_schur == size_of_the_problem */ 1033 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 1034 reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */ 1035 factor_workaround = PETSC_FALSE; 1036 solver_S = PETSC_FALSE; 1037 } 1038 1039 if (reuse_solvers) { 1040 Mat A_II,Afake; 1041 Vec vec1_B; 1042 PCBDDCReuseSolvers msolv_ctx; 1043 PetscInt n_R; 1044 1045 if (sub_schurs->reuse_solver) { 1046 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 1047 } else { 1048 ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr); 1049 } 1050 msolv_ctx = sub_schurs->reuse_solver; 1051 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1052 ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 1053 msolv_ctx->F = F; 1054 ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr); 1055 /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */ 1056 { 1057 PetscScalar *array; 1058 PetscInt n; 1059 1060 ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr); 1061 ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr); 1062 ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr); 1063 ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr); 1064 } 1065 msolv_ctx->has_vertices = schur_has_vertices; 1066 1067 /* interior solver */ 1068 ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr); 1069 ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr); 1070 ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr); 1071 ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr); 1072 ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr); 1073 ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr); 1074 1075 /* correction solver */ 1076 ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr); 1077 ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr); 1078 ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr); 1079 ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr); 1080 ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr); 1081 1082 /* scatter and vecs for Schur complement solver */ 1083 ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr); 1084 ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr); 1085 if (!schur_has_vertices) { 1086 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr); 1087 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 1088 ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr); 1089 msolv_ctx->is_R = is_A_all; 1090 } else { 1091 IS is_B_all; 1092 const PetscInt* idxs; 1093 PetscInt dual,n_v,n; 1094 1095 ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr); 1096 dual = size_schur - n_v; 1097 ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr); 1098 ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr); 1099 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr); 1100 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr); 1101 ierr = ISDestroy(&is_B_all);CHKERRQ(ierr); 1102 ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr); 1103 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 1104 ierr = ISDestroy(&is_B_all);CHKERRQ(ierr); 1105 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr); 1106 ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr); 1107 } 1108 ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr); 1109 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr); 1110 ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1111 ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1112 ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr); 1113 ierr = MatDestroy(&Afake);CHKERRQ(ierr); 1114 ierr = VecDestroy(&vec1_B);CHKERRQ(ierr); 1115 1116 /* communicate benign info to solver context */ 1117 if (benign_n) { 1118 PetscScalar *array; 1119 1120 msolv_ctx->benign_n = benign_n; 1121 msolv_ctx->benign_zerodiag_subs = is_p_r; 1122 ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr); 1123 msolv_ctx->benign_csAIB = cs_AIB_mat; 1124 ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr); 1125 ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr); 1126 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 1127 ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr); 1128 msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat; 1129 } 1130 } else { 1131 if (sub_schurs->reuse_solver) { 1132 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 1133 } 1134 ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr); 1135 } 1136 ierr = MatDestroy(&A);CHKERRQ(ierr); 1137 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 1138 1139 /* Work arrays */ 1140 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 1141 1142 /* matrices for deluxe scaling and adaptive selection */ 1143 if (!sub_schurs->sum_S_Ej_all) { 1144 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1145 ierr = MatSetSizes(sub_schurs->sum_S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr); 1146 ierr = MatSetType(sub_schurs->sum_S_Ej_all,MATAIJ);CHKERRQ(ierr); 1147 ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_all,0,nnz);CHKERRQ(ierr); 1148 } 1149 if (compute_Stilda) { 1150 if (!sub_schurs->sum_S_Ej_tilda_all) { 1151 ierr = MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1152 } 1153 if (!sub_schurs->sum_S_Ej_inv_all && deluxe) { 1154 ierr = MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1155 } 1156 } 1157 1158 /* S_Ej_all */ 1159 cum = cum2 = 0; 1160 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 1161 for (i=0;i<sub_schurs->n_subs;i++) { 1162 PetscInt j; 1163 1164 /* get S_E */ 1165 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1166 if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */ 1167 PetscInt k; 1168 for (k=0;k<subset_size;k++) { 1169 for (j=k;j<subset_size;j++) { 1170 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 1171 work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]); 1172 } 1173 } 1174 } else { /* just copy to workspace */ 1175 PetscInt k; 1176 for (k=0;k<subset_size;k++) { 1177 for (j=0;j<subset_size;j++) { 1178 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 1179 } 1180 } 1181 } 1182 /* insert S_E values */ 1183 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 1184 if (sub_schurs->change) { 1185 Mat change_sub,SEj,T; 1186 1187 /* change basis */ 1188 ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr); 1189 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1190 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1191 Mat T2; 1192 ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr); 1193 ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1194 ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1195 ierr = MatDestroy(&T2);CHKERRQ(ierr); 1196 } else { 1197 ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1198 } 1199 ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1200 ierr = MatDestroy(&T);CHKERRQ(ierr); 1201 ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr); 1202 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1203 } 1204 if (deluxe) { 1205 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 1206 /* if adaptivity is requested, invert S_E blocks */ 1207 if (compute_Stilda) { 1208 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 1209 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1210 if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */ 1211 PetscInt k; 1212 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr)); 1213 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 1214 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr)); 1215 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 1216 for (k=0;k<subset_size;k++) { 1217 for (j=k;j<subset_size;j++) { 1218 work[j*subset_size+k] = work[k*subset_size+j]; 1219 } 1220 } 1221 } else { 1222 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr)); 1223 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 1224 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1225 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 1226 } 1227 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1228 ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 1229 } 1230 } else if (compute_Stilda) { /* not using deluxe */ 1231 Mat SEj; 1232 Vec D; 1233 PetscScalar *array; 1234 1235 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1236 ierr = VecGetArray(Dall,&array);CHKERRQ(ierr); 1237 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr); 1238 ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr); 1239 ierr = VecShift(D,-1.);CHKERRQ(ierr); 1240 ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr); 1241 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1242 ierr = VecDestroy(&D);CHKERRQ(ierr); 1243 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 1244 } 1245 cum += subset_size; 1246 cum2 += subset_size*(size_schur + 1); 1247 } 1248 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 1249 1250 if (solver_S) { 1251 ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1252 } 1253 1254 schur_factor = NULL; 1255 if (compute_Stilda && size_active_schur) { 1256 1257 if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */ 1258 PetscInt j; 1259 for (j=0;j<size_schur;j++) dummy_idx[j] = j; 1260 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 1261 } else { 1262 Mat S_all_inv=NULL; 1263 if (solver_S) { 1264 /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1. 1265 The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */ 1266 if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */ 1267 PetscScalar *data; 1268 PetscInt nd = 0; 1269 1270 if (!sub_schurs->is_posdef) { 1271 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices"); 1272 } 1273 ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr); 1274 ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr); 1275 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1276 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1277 } 1278 1279 /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */ 1280 if (schur_has_vertices) { 1281 Mat M; 1282 PetscScalar *tdata; 1283 PetscInt nv = 0, news; 1284 1285 ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr); 1286 news = size_active_schur + nv; 1287 ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr); 1288 for (i=0;i<size_active_schur;i++) { 1289 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr); 1290 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr); 1291 } 1292 for (i=0;i<nv;i++) { 1293 PetscInt k = i+size_active_schur; 1294 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr); 1295 } 1296 1297 ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr); 1298 ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1299 ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr); 1300 /* save the factors */ 1301 cum = 0; 1302 ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr); 1303 for (i=0;i<size_active_schur;i++) { 1304 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr); 1305 cum += size_active_schur - i; 1306 } 1307 for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)])); 1308 ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr); 1309 /* move back just the active dofs to the Schur complement */ 1310 for (i=0;i<size_active_schur;i++) { 1311 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1312 } 1313 ierr = PetscFree(tdata);CHKERRQ(ierr); 1314 ierr = MatDestroy(&M);CHKERRQ(ierr); 1315 } else { /* we can factorize and invert just the activedofs */ 1316 Mat M; 1317 1318 ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr); 1319 ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr); 1320 ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1321 ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr); 1322 ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr); 1323 ierr = MatDestroy(&M);CHKERRQ(ierr); 1324 for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = 1.0/data[(i+size_active_schur)*(size_schur+1)]; 1325 } 1326 ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr); 1327 } else { /* use MatFactor calls to invert S */ 1328 ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr); 1329 ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr); 1330 } 1331 } else { /* we need to invert explicitly since we are not using MatFactor for S */ 1332 ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr); 1333 S_all_inv = S_all; 1334 ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr); 1335 ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr); 1336 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1337 if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */ 1338 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr)); 1339 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 1340 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr)); 1341 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 1342 } else { 1343 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr)); 1344 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 1345 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1346 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 1347 } 1348 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1349 ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr); 1350 } 1351 /* S_Ej_tilda_all */ 1352 cum = cum2 = 0; 1353 ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr); 1354 for (i=0;i<sub_schurs->n_subs;i++) { 1355 PetscInt j; 1356 1357 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1358 /* get (St^-1)_E */ 1359 /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1 1360 will be properly accessed later during adaptive selection */ 1361 if (S_lower_triangular) { 1362 PetscInt k; 1363 if (sub_schurs->change) { 1364 for (k=0;k<subset_size;k++) { 1365 for (j=k;j<subset_size;j++) { 1366 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 1367 work[j*subset_size+k] = work[k*subset_size+j]; 1368 } 1369 } 1370 } else { 1371 for (k=0;k<subset_size;k++) { 1372 for (j=k;j<subset_size;j++) { 1373 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 1374 } 1375 } 1376 } 1377 } else { 1378 PetscInt k; 1379 for (k=0;k<subset_size;k++) { 1380 for (j=0;j<subset_size;j++) { 1381 work[k*subset_size+j] = S_data[cum2+k*size_schur+j]; 1382 } 1383 } 1384 } 1385 if (sub_schurs->change) { 1386 Mat change_sub,SEj,T; 1387 1388 /* change basis */ 1389 ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr); 1390 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1391 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1392 Mat T2; 1393 ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr); 1394 ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1395 ierr = MatDestroy(&T2);CHKERRQ(ierr); 1396 ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1397 } else { 1398 ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1399 } 1400 ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1401 ierr = MatDestroy(&T);CHKERRQ(ierr); 1402 /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */ 1403 ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr); 1404 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1405 } 1406 for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j; 1407 ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 1408 cum += subset_size; 1409 cum2 += subset_size*(size_schur + 1); 1410 } 1411 ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr); 1412 if (solver_S) { 1413 if (schur_has_vertices) { 1414 ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr); 1415 } else { 1416 ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr); 1417 } 1418 } 1419 ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr); 1420 } 1421 1422 /* move back factors if needed */ 1423 if (schur_has_vertices) { 1424 Mat S_tmp; 1425 PetscInt nd = 0; 1426 1427 if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 1428 ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr); 1429 if (sub_schurs->is_posdef) { 1430 PetscScalar *data; 1431 1432 ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr); 1433 ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1434 1435 if (S_lower_triangular) { 1436 cum = 0; 1437 for (i=0;i<size_active_schur;i++) { 1438 ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr); 1439 cum += size_active_schur-i; 1440 } 1441 } else { 1442 ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1443 } 1444 if (sub_schurs->is_dir) { 1445 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1446 for (i=0;i<nd;i++) { 1447 data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i]; 1448 } 1449 } 1450 /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions, 1451 set the diagonal entry of the Schur factor to a very large value */ 1452 for (i=size_active_schur+nd;i<size_schur;i++) { 1453 data[i*(size_schur+1)] = infty; 1454 } 1455 ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr); 1456 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices"); 1457 ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr); 1458 } 1459 } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */ 1460 PetscScalar *data; 1461 PetscInt nd = 0; 1462 1463 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1464 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1465 } 1466 ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr); 1467 ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr); 1468 for (i=0;i<size_active_schur;i++) { 1469 ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr); 1470 } 1471 for (i=size_active_schur+nd;i<size_schur;i++) { 1472 ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr); 1473 data[i*(size_schur+1)] = infty; 1474 } 1475 ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr); 1476 ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1477 } 1478 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 1479 ierr = PetscFree(schur_factor);CHKERRQ(ierr); 1480 ierr = VecDestroy(&Dall);CHKERRQ(ierr); 1481 } 1482 ierr = PetscFree(nnz);CHKERRQ(ierr); 1483 ierr = MatDestroy(&F);CHKERRQ(ierr); 1484 ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr); 1485 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 1486 ierr = MatDestroy(&A_BB);CHKERRQ(ierr); 1487 ierr = MatDestroy(&A_IB);CHKERRQ(ierr); 1488 ierr = MatDestroy(&A_BI);CHKERRQ(ierr); 1489 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1490 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1491 if (compute_Stilda) { 1492 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1493 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1494 if (deluxe) { 1495 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1496 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1497 } 1498 } 1499 1500 /* Global matrix of all assembled Schur on subsets */ 1501 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 1502 ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 1503 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1504 1505 /* Get local part of (\sum_j S_Ej) */ 1506 ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr); 1507 ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1508 1509 /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */ 1510 if (compute_Stilda) { 1511 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1512 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1513 ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1514 ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1515 if (deluxe) { 1516 ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1517 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1518 ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 1519 ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1520 } else { 1521 PetscScalar *array; 1522 PetscInt cum; 1523 1524 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr); 1525 cum = 0; 1526 for (i=0;i<sub_schurs->n_subs;i++) { 1527 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1528 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 1529 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1530 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 1531 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 1532 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 1533 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 1534 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1535 cum += subset_size*subset_size; 1536 } 1537 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr); 1538 ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1539 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1540 sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all; 1541 } 1542 } 1543 ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr); 1544 1545 /* free workspace */ 1546 ierr = PetscFree(submats);CHKERRQ(ierr); 1547 ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 1548 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 1549 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 1550 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 1551 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 1552 PetscFunctionReturn(0); 1553 } 1554 1555 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc) 1556 { 1557 IS *faces,*edges,*all_cc,vertices; 1558 PetscInt i,n_faces,n_edges,n_all_cc; 1559 PetscBool is_sorted; 1560 PetscErrorCode ierr; 1561 1562 PetscFunctionBegin; 1563 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 1564 if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 1565 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 1566 if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 1567 1568 /* reset any previous data */ 1569 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 1570 1571 /* get index sets for faces and edges (already sorted by global ordering) */ 1572 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1573 n_all_cc = n_faces+n_edges; 1574 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 1575 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 1576 for (i=0;i<n_faces;i++) { 1577 if (copycc) { 1578 ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr); 1579 } else { 1580 ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr); 1581 all_cc[i] = faces[i]; 1582 } 1583 } 1584 for (i=0;i<n_edges;i++) { 1585 if (copycc) { 1586 ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr); 1587 } else { 1588 ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr); 1589 all_cc[n_faces+i] = edges[i]; 1590 } 1591 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 1592 } 1593 ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr); 1594 sub_schurs->is_vertices = vertices; 1595 ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1596 sub_schurs->is_dir = NULL; 1597 ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr); 1598 1599 /* Determine if MatFactor can be used */ 1600 sub_schurs->schur_explicit = PETSC_FALSE; 1601 #if defined(PETSC_HAVE_MUMPS) 1602 sub_schurs->schur_explicit = PETSC_TRUE; 1603 #endif 1604 #if defined(PETSC_HAVE_MKL_PARDISO) 1605 sub_schurs->schur_explicit = PETSC_TRUE; 1606 #endif 1607 1608 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 1609 sub_schurs->is_I = is_I; 1610 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 1611 sub_schurs->is_B = is_B; 1612 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 1613 sub_schurs->l2gmap = graph->l2gmap; 1614 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 1615 sub_schurs->BtoNmap = BtoNmap; 1616 sub_schurs->n_subs = n_all_cc; 1617 sub_schurs->is_subs = all_cc; 1618 if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */ 1619 for (i=0;i<sub_schurs->n_subs;i++) { 1620 ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 1621 } 1622 } 1623 sub_schurs->S_Ej_all = NULL; 1624 sub_schurs->sum_S_Ej_all = NULL; 1625 sub_schurs->sum_S_Ej_inv_all = NULL; 1626 sub_schurs->sum_S_Ej_tilda_all = NULL; 1627 sub_schurs->is_Ej_all = NULL; 1628 PetscFunctionReturn(0); 1629 } 1630 1631 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 1632 { 1633 PCBDDCSubSchurs schurs_ctx; 1634 PetscErrorCode ierr; 1635 1636 PetscFunctionBegin; 1637 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 1638 schurs_ctx->n_subs = 0; 1639 *sub_schurs = schurs_ctx; 1640 PetscFunctionReturn(0); 1641 } 1642 1643 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 1644 { 1645 PetscInt i; 1646 PetscErrorCode ierr; 1647 1648 PetscFunctionBegin; 1649 if (!sub_schurs) PetscFunctionReturn(0); 1650 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 1651 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 1652 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 1653 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 1654 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 1655 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 1656 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1657 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1658 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1659 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 1660 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 1661 ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr); 1662 ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr); 1663 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 1664 for (i=0;i<sub_schurs->n_subs;i++) { 1665 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 1666 } 1667 if (sub_schurs->n_subs) { 1668 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 1669 } 1670 if (sub_schurs->reuse_solver) { 1671 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 1672 } 1673 ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr); 1674 if (sub_schurs->change) { 1675 for (i=0;i<sub_schurs->n_subs;i++) { 1676 ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr); 1677 ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr); 1678 } 1679 } 1680 ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr); 1681 ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr); 1682 sub_schurs->n_subs = 0; 1683 PetscFunctionReturn(0); 1684 } 1685 1686 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs) 1687 { 1688 PetscErrorCode ierr; 1689 1690 PetscFunctionBegin; 1691 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 1692 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 1693 PetscFunctionReturn(0); 1694 } 1695 1696 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 1697 { 1698 PetscInt i,j,n; 1699 PetscErrorCode ierr; 1700 1701 PetscFunctionBegin; 1702 n = 0; 1703 for (i=-n_prev;i<0;i++) { 1704 PetscInt start_dof = queue_tip[i]; 1705 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 1706 PetscInt dof = adjncy[j]; 1707 if (!PetscBTLookup(touched,dof)) { 1708 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 1709 queue_tip[n] = dof; 1710 n++; 1711 } 1712 } 1713 } 1714 *n_added = n; 1715 PetscFunctionReturn(0); 1716 } 1717