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