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