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