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