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