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