1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <../src/mat/impls/dense/seq/dense.h> 4 #include <petscblaslapack.h> 5 6 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 7 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*); 8 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec); 9 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec); 10 11 /* if v2 is not present, correction is done in-place */ 12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full) 13 { 14 PetscScalar *array; 15 PetscScalar *array2; 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 if (!ctx->benign_n) PetscFunctionReturn(0); 20 if (sol && full) { 21 PetscInt n_I,size_schur; 22 23 /* get sizes */ 24 ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr); 25 ierr = VecGetSize(v,&n_I);CHKERRQ(ierr); 26 n_I = n_I - size_schur; 27 /* get schur sol from array */ 28 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 29 ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr); 30 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 31 /* apply interior sol correction */ 32 ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr); 33 ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 34 ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr); 35 } 36 if (v2) { 37 PetscInt nl; 38 39 ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 40 ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr); 41 ierr = VecGetArray(v2,&array2);CHKERRQ(ierr); 42 ierr = PetscArraycpy(array2,array,nl);CHKERRQ(ierr); 43 } else { 44 ierr = VecGetArray(v,&array);CHKERRQ(ierr); 45 array2 = array; 46 } 47 if (!sol) { /* change rhs */ 48 PetscInt n; 49 for (n=0;n<ctx->benign_n;n++) { 50 PetscScalar sum = 0.; 51 const PetscInt *cols; 52 PetscInt nz,i; 53 54 ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr); 55 ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 56 for (i=0;i<nz-1;i++) sum += array[cols[i]]; 57 #if defined(PETSC_USE_COMPLEX) 58 sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz)); 59 #else 60 sum = -sum/nz; 61 #endif 62 for (i=0;i<nz-1;i++) array2[cols[i]] += sum; 63 ctx->benign_save_vals[n] = array2[cols[nz-1]]; 64 array2[cols[nz-1]] = sum; 65 ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 66 } 67 } else { 68 PetscInt n; 69 for (n=0;n<ctx->benign_n;n++) { 70 PetscScalar sum = 0.; 71 const PetscInt *cols; 72 PetscInt nz,i; 73 ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr); 74 ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 75 for (i=0;i<nz-1;i++) sum += array[cols[i]]; 76 #if defined(PETSC_USE_COMPLEX) 77 sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz)); 78 #else 79 sum = -sum/nz; 80 #endif 81 for (i=0;i<nz-1;i++) array2[cols[i]] += sum; 82 array2[cols[nz-1]] = ctx->benign_save_vals[n]; 83 ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr); 84 } 85 } 86 if (v2) { 87 ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 88 ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr); 89 } else { 90 ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 91 } 92 if (!sol && full) { 93 Vec usedv; 94 PetscInt n_I,size_schur; 95 96 /* get sizes */ 97 ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr); 98 ierr = VecGetSize(v,&n_I);CHKERRQ(ierr); 99 n_I = n_I - size_schur; 100 /* compute schur rhs correction */ 101 if (v2) { 102 usedv = v2; 103 } else { 104 usedv = v; 105 } 106 /* apply schur rhs correction */ 107 ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr); 108 ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr); 109 ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr); 110 ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr); 111 ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 112 ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 113 } 114 PetscFunctionReturn(0); 115 } 116 117 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full) 118 { 119 PCBDDCReuseSolvers ctx; 120 PetscBool copy = PETSC_FALSE; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr); 125 if (full) { 126 #if defined(PETSC_HAVE_MUMPS) 127 ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr); 128 #endif 129 #if defined(PETSC_HAVE_MKL_PARDISO) 130 ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr); 131 #endif 132 copy = ctx->has_vertices; 133 } else { /* interior solver */ 134 #if defined(PETSC_HAVE_MUMPS) 135 ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr); 136 #endif 137 #if defined(PETSC_HAVE_MKL_PARDISO) 138 ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr); 139 #endif 140 copy = PETSC_TRUE; 141 } 142 /* copy rhs into factored matrix workspace */ 143 if (copy) { 144 PetscInt n; 145 PetscScalar *array,*array_solver; 146 147 ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr); 148 ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 149 ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr); 150 ierr = PetscArraycpy(array_solver,array,n);CHKERRQ(ierr); 151 ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr); 152 ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr); 153 154 ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr); 155 if (transpose) { 156 ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 157 } else { 158 ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr); 159 } 160 ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr); 161 162 /* get back data to caller worskpace */ 163 ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr); 164 ierr = VecGetArray(sol,&array);CHKERRQ(ierr); 165 ierr = PetscArraycpy(array,array_solver,n);CHKERRQ(ierr); 166 ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr); 167 ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr); 168 } else { 169 if (ctx->benign_n) { 170 ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr); 171 if (transpose) { 172 ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr); 173 } else { 174 ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr); 175 } 176 ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr); 177 } else { 178 if (transpose) { 179 ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr); 180 } else { 181 ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr); 182 } 183 } 184 } 185 /* restore defaults */ 186 #if defined(PETSC_HAVE_MUMPS) 187 ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr); 188 #endif 189 #if defined(PETSC_HAVE_MKL_PARDISO) 190 ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr); 191 #endif 192 PetscFunctionReturn(0); 193 } 194 195 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol) 196 { 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol) 205 { 206 PetscErrorCode ierr; 207 208 PetscFunctionBegin; 209 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr); 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 222 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol) 223 { 224 PetscErrorCode ierr; 225 226 PetscFunctionBegin; 227 ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer) 232 { 233 PCBDDCReuseSolvers ctx; 234 PetscBool iascii; 235 PetscErrorCode ierr; 236 237 PetscFunctionBegin; 238 ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr); 239 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 240 if (iascii) { 241 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 242 } 243 ierr = MatView(ctx->F,viewer);CHKERRQ(ierr); 244 if (iascii) { 245 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 246 } 247 PetscFunctionReturn(0); 248 } 249 250 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse) 251 { 252 PetscInt i; 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 ierr = MatDestroy(&reuse->F);CHKERRQ(ierr); 257 ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr); 258 ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr); 259 ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr); 260 ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr); 261 ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr); 262 ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr); 263 ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr); 264 ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr); 265 ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr); 266 for (i=0;i<reuse->benign_n;i++) { 267 ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr); 268 } 269 ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr); 270 ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr); 271 ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr); 272 ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr); 273 ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr); 274 ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 279 { 280 Mat B, C, D, Bd, Cd, AinvBd; 281 KSP ksp; 282 PC pc; 283 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 284 PetscReal fill = 2.0; 285 PetscInt n_I; 286 PetscMPIInt size; 287 PetscErrorCode ierr; 288 289 PetscFunctionBegin; 290 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRMPI(ierr); 291 if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 292 if (reuse == MAT_REUSE_MATRIX) { 293 PetscBool Sdense; 294 295 ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr); 296 if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense"); 297 } 298 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 299 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 300 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 301 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 302 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 303 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 304 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 305 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 306 ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr); 307 if (n_I) { 308 if (!Bdense) { 309 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 310 } else { 311 Bd = B; 312 } 313 314 if (isLU || isILU || isCHOL) { 315 Mat fact; 316 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 317 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 318 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 319 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 320 } else { 321 PetscBool ex = PETSC_TRUE; 322 323 if (ex) { 324 Mat Ainvd; 325 326 ierr = PCComputeOperator(pc, MATDENSE, &Ainvd);CHKERRQ(ierr); 327 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 328 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 329 } else { 330 Vec sol,rhs; 331 PetscScalar *arrayrhs,*arraysol; 332 PetscInt i,nrhs,n; 333 334 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 335 ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr); 336 ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr); 337 ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr); 338 ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr); 339 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 340 for (i=0;i<nrhs;i++) { 341 ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr); 342 ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr); 343 ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr); 344 ierr = VecResetArray(rhs);CHKERRQ(ierr); 345 ierr = VecResetArray(sol);CHKERRQ(ierr); 346 } 347 ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr); 348 ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr); 349 } 350 } 351 if (!Bdense & !issym) { 352 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 353 } 354 355 if (!issym) { 356 if (!Cdense) { 357 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 358 } else { 359 Cd = C; 360 } 361 ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 362 if (!Cdense) { 363 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 364 } 365 } else { 366 ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 367 if (!Bdense) { 368 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 369 } 370 } 371 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 372 } 373 374 if (D) { 375 Mat Dd; 376 PetscBool Ddense; 377 378 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 379 if (!Ddense) { 380 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 381 } else { 382 Dd = D; 383 } 384 if (n_I) { 385 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 386 } else { 387 if (reuse == MAT_INITIAL_MATRIX) { 388 ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr); 389 } else { 390 ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 391 } 392 } 393 if (!Ddense) { 394 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 395 } 396 } else { 397 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 398 } 399 PetscFunctionReturn(0); 400 } 401 402 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) 403 { 404 Mat F,A_II,A_IB,A_BI,A_BB,AE_II; 405 Mat S_all; 406 Vec gstash,lstash; 407 VecScatter sstash; 408 IS is_I,is_I_layer; 409 IS all_subsets,all_subsets_mult,all_subsets_n; 410 PetscScalar *stasharray,*Bwork; 411 PetscInt *nnz,*all_local_idx_N; 412 PetscInt *auxnum1,*auxnum2; 413 PetscInt i,subset_size,max_subset_size; 414 PetscInt n_B,extra,local_size,global_size; 415 PetscInt local_stash_size; 416 PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 417 MPI_Comm comm_n; 418 PetscBool deluxe = PETSC_TRUE; 419 PetscBool use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE; 420 PetscViewer matl_dbg_viewer = NULL; 421 PetscErrorCode ierr; 422 PetscBool flg; 423 424 PetscFunctionBegin; 425 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 426 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 427 if (Ain) { 428 ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr); 429 sub_schurs->A = Ain; 430 } 431 432 ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr); 433 sub_schurs->S = Sin; 434 if (sub_schurs->schur_explicit) { 435 sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A); 436 } 437 438 /* preliminary checks */ 439 if (!sub_schurs->schur_explicit && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO"); 440 441 if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE; 442 443 /* debug (MATLAB) */ 444 if (sub_schurs->debug) { 445 PetscMPIInt size,rank; 446 PetscInt nr,*print_schurs_ranks,print_schurs = PETSC_FALSE; 447 448 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);CHKERRMPI(ierr); 449 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRMPI(ierr); 450 nr = size; 451 ierr = PetscMalloc1(nr,&print_schurs_ranks);CHKERRQ(ierr); 452 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr); 453 ierr = PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);CHKERRQ(ierr); 454 if (!flg) print_schurs = PETSC_TRUE; 455 else { 456 print_schurs = PETSC_FALSE; 457 for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; } 458 } 459 ierr = PetscOptionsEnd();CHKERRQ(ierr); 460 ierr = PetscFree(print_schurs_ranks);CHKERRQ(ierr); 461 if (print_schurs) { 462 char filename[256]; 463 464 ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);CHKERRQ(ierr); 465 ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&matl_dbg_viewer);CHKERRQ(ierr); 466 ierr = PetscViewerPushFormat(matl_dbg_viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 467 } 468 } 469 470 /* restrict work on active processes */ 471 if (sub_schurs->restrict_comm) { 472 PetscSubcomm subcomm; 473 PetscMPIInt color,rank; 474 475 color = 0; 476 if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */ 477 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRMPI(ierr); 478 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr); 479 ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr); 480 ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr); 481 ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr); 482 ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr); 483 if (!sub_schurs->n_subs) { 484 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 } else { 488 ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); 489 } 490 491 /* get Schur complement matrices */ 492 if (!sub_schurs->schur_explicit) { 493 Mat tA_IB,tA_BI,tA_BB; 494 PetscBool isseqsbaij; 495 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr); 496 ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 497 if (isseqsbaij) { 498 ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr); 499 ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr); 500 ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr); 501 } else { 502 ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr); 503 A_BB = tA_BB; 504 ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr); 505 A_IB = tA_IB; 506 ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr); 507 A_BI = tA_BI; 508 } 509 } else { 510 A_II = NULL; 511 A_IB = NULL; 512 A_BI = NULL; 513 A_BB = NULL; 514 } 515 S_all = NULL; 516 517 /* determine interior problems */ 518 ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 519 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 520 PetscBT touched; 521 const PetscInt* idx_B; 522 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 523 524 if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 525 /* get sizes */ 526 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 527 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 528 529 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 530 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 531 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 532 533 /* all boundary dofs must be skipped when adding layers */ 534 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 535 for (j=0;j<n_B;j++) { 536 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 537 } 538 ierr = PetscArraycpy(local_numbering,idx_B,n_B);CHKERRQ(ierr); 539 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 540 541 /* add prescribed number of layers of dofs */ 542 n_local_dofs = n_B; 543 n_prev_added = n_B; 544 for (layer=0;layer<nlayers;layer++) { 545 PetscInt n_added; 546 if (n_local_dofs == n_I+n_B) break; 547 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); 548 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 549 n_prev_added = n_added; 550 n_local_dofs += n_added; 551 if (!n_added) break; 552 } 553 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 554 555 /* IS for I layer dofs in original numbering */ 556 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); 557 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 558 ierr = ISSort(is_I_layer);CHKERRQ(ierr); 559 /* IS for I layer dofs in I numbering */ 560 if (!sub_schurs->schur_explicit) { 561 ISLocalToGlobalMapping ItoNmap; 562 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 563 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr); 564 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 565 566 /* II block */ 567 ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 568 } 569 } else { 570 PetscInt n_I; 571 572 /* IS for I dofs in original numbering */ 573 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 574 is_I_layer = sub_schurs->is_I; 575 576 /* IS for I dofs in I numbering (strided 1) */ 577 if (!sub_schurs->schur_explicit) { 578 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 579 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 580 581 /* II block is the same */ 582 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 583 AE_II = A_II; 584 } 585 } 586 587 /* Get info on subset sizes and sum of all subsets sizes */ 588 max_subset_size = 0; 589 local_size = 0; 590 for (i=0;i<sub_schurs->n_subs;i++) { 591 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 592 max_subset_size = PetscMax(subset_size,max_subset_size); 593 local_size += subset_size; 594 } 595 596 /* Work arrays for local indices */ 597 extra = 0; 598 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 599 if (sub_schurs->schur_explicit && is_I_layer) { 600 ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr); 601 } 602 ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr); 603 if (extra) { 604 const PetscInt *idxs; 605 ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr); 606 ierr = PetscArraycpy(all_local_idx_N,idxs,extra);CHKERRQ(ierr); 607 ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr); 608 } 609 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr); 610 ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr); 611 612 /* Get local indices in local numbering */ 613 local_size = 0; 614 local_stash_size = 0; 615 for (i=0;i<sub_schurs->n_subs;i++) { 616 const PetscInt *idxs; 617 618 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 619 ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 620 /* start (smallest in global ordering) and multiplicity */ 621 auxnum1[i] = idxs[0]; 622 auxnum2[i] = subset_size*subset_size; 623 /* subset indices in local numbering */ 624 ierr = PetscArraycpy(all_local_idx_N+local_size+extra,idxs,subset_size);CHKERRQ(ierr); 625 ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 626 local_size += subset_size; 627 local_stash_size += subset_size*subset_size; 628 } 629 630 /* allocate extra workspace needed only for GETRI or SYTRF */ 631 use_potr = use_sytr = PETSC_FALSE; 632 if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) { 633 use_potr = PETSC_TRUE; 634 } else if (sub_schurs->is_symmetric) { 635 use_sytr = PETSC_TRUE; 636 } 637 if (local_size && !use_potr) { 638 PetscScalar lwork,dummyscalar = 0.; 639 PetscBLASInt dummyint = 0; 640 641 B_lwork = -1; 642 ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 643 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 644 if (use_sytr) { 645 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr)); 646 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr); 647 } else { 648 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr)); 649 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 650 } 651 ierr = PetscFPTrapPop();CHKERRQ(ierr); 652 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 653 ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr); 654 } else { 655 Bwork = NULL; 656 pivots = NULL; 657 } 658 659 /* prepare data for summing up properly schurs on subsets */ 660 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr); 661 ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr); 662 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 663 ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr); 664 ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr); 665 ierr = ISDestroy(&all_subsets);CHKERRQ(ierr); 666 ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr); 667 ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr); 668 if (i != local_stash_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_stash_size); 669 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,local_stash_size,NULL,&lstash);CHKERRQ(ierr); 670 ierr = VecCreateMPI(comm_n,PETSC_DECIDE,global_size,&gstash);CHKERRQ(ierr); 671 ierr = VecScatterCreate(lstash,NULL,gstash,all_subsets_n,&sstash);CHKERRQ(ierr); 672 ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr); 673 674 /* subset indices in local boundary numbering */ 675 if (!sub_schurs->is_Ej_all) { 676 PetscInt *all_local_idx_B; 677 678 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 679 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 680 if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size); 681 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 682 } 683 684 if (change) { 685 ISLocalToGlobalMapping BtoS; 686 IS change_primal_B; 687 IS change_primal_all; 688 689 if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 690 if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 691 ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr); 692 for (i=0;i<sub_schurs->n_subs;i++) { 693 ISLocalToGlobalMapping NtoS; 694 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr); 695 ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr); 696 ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr); 697 } 698 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr); 699 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr); 700 ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr); 701 ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr); 702 ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr); 703 ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr); 704 for (i=0;i<sub_schurs->n_subs;i++) { 705 Mat change_sub; 706 707 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 708 ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr); 709 ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr); 710 if (!sub_schurs->change_with_qr) { 711 ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr); 712 } else { 713 Mat change_subt; 714 ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr); 715 ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr); 716 ierr = MatDestroy(&change_subt);CHKERRQ(ierr); 717 } 718 ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr); 719 ierr = MatDestroy(&change_sub);CHKERRQ(ierr); 720 ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr); 721 ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr); 722 } 723 ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr); 724 } 725 726 /* Local matrix of all local Schur on subsets (transposed) */ 727 if (!sub_schurs->S_Ej_all) { 728 Mat T; 729 PetscScalar *v; 730 PetscInt *ii,*jj; 731 PetscInt cum,i,j,k; 732 733 /* MatSeqAIJSetPreallocation + MatSetValues is slow for these kind of matrices (may have large blocks) 734 Allocate properly a representative matrix and duplicate */ 735 ierr = PetscMalloc3(local_size+1,&ii,local_stash_size,&jj,local_stash_size,&v);CHKERRQ(ierr); 736 ii[0] = 0; 737 cum = 0; 738 for (i=0;i<sub_schurs->n_subs;i++) { 739 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 740 for (j=0;j<subset_size;j++) { 741 const PetscInt row = cum+j; 742 PetscInt col = cum; 743 744 ii[row+1] = ii[row] + subset_size; 745 for (k=ii[row];k<ii[row+1];k++) { 746 jj[k] = col; 747 col++; 748 } 749 } 750 cum += subset_size; 751 } 752 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,local_size,local_size,ii,jj,v,&T);CHKERRQ(ierr); 753 ierr = MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 754 ierr = MatDestroy(&T);CHKERRQ(ierr); 755 ierr = PetscFree3(ii,jj,v);CHKERRQ(ierr); 756 } 757 /* matrices for deluxe scaling and adaptive selection */ 758 if (compute_Stilda) { 759 if (!sub_schurs->sum_S_Ej_tilda_all) { 760 ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 761 } 762 if (!sub_schurs->sum_S_Ej_inv_all && deluxe) { 763 ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 764 } 765 } 766 767 /* Compute Schur complements explicitly */ 768 F = NULL; 769 if (!sub_schurs->schur_explicit) { 770 /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested; 771 it is not efficient, unless the economic version of the scaling is used */ 772 Mat S_Ej_expl; 773 PetscScalar *work; 774 PetscInt j,*dummy_idx; 775 PetscBool Sdense; 776 777 ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 778 local_size = 0; 779 for (i=0;i<sub_schurs->n_subs;i++) { 780 IS is_subset_B; 781 Mat AE_EE,AE_IE,AE_EI,S_Ej; 782 783 /* subsets in original and boundary numbering */ 784 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr); 785 /* EE block */ 786 ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr); 787 /* IE block */ 788 ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr); 789 /* EI block */ 790 if (sub_schurs->is_symmetric) { 791 ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 792 } else if (sub_schurs->is_hermitian) { 793 ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 794 } else { 795 ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr); 796 } 797 ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr); 798 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr); 799 ierr = MatDestroy(&AE_EE);CHKERRQ(ierr); 800 ierr = MatDestroy(&AE_IE);CHKERRQ(ierr); 801 ierr = MatDestroy(&AE_EI);CHKERRQ(ierr); 802 if (AE_II == A_II) { /* we can reuse the same ksp */ 803 KSP ksp; 804 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 805 ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr); 806 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 807 KSP origksp,schurksp; 808 PC origpc,schurpc; 809 KSPType ksp_type; 810 PetscInt n_internal; 811 PetscBool ispcnone; 812 813 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 814 ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr); 815 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 816 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 817 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 818 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 819 ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr); 820 if (!ispcnone) { 821 PCType pc_type; 822 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 823 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 824 } else { 825 ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr); 826 } 827 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 828 if (!n_internal) { /* UMFPACK gives error with 0 sized problems */ 829 MatSolverType solver = NULL; 830 ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr); 831 if (solver) { 832 ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr); 833 } 834 } 835 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 836 } 837 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 838 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 839 ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 840 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 841 if (Sdense) { 842 for (j=0;j<subset_size;j++) { 843 dummy_idx[j]=local_size+j; 844 } 845 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 846 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 847 ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 848 ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr); 849 local_size += subset_size; 850 } 851 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 852 /* free */ 853 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 854 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 855 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 856 } else { 857 Mat A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL; 858 Vec Dall = NULL; 859 IS is_A_all,*is_p_r = NULL; 860 MatType Stype; 861 PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL; 862 PetscScalar *SEj_arr = NULL,*SEjinv_arr = NULL; 863 const PetscScalar *rS_data; 864 PetscInt n,n_I,size_schur,size_active_schur,cum,cum2; 865 PetscBool economic,solver_S,S_lower_triangular = PETSC_FALSE; 866 PetscBool schur_has_vertices,factor_workaround; 867 PetscBool use_cholesky; 868 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 869 PetscBool oldpin; 870 #endif 871 872 /* get sizes */ 873 n_I = 0; 874 if (is_I_layer) { 875 ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr); 876 } 877 economic = PETSC_FALSE; 878 ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr); 879 if (cum != n_I) economic = PETSC_TRUE; 880 ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr); 881 size_active_schur = local_size; 882 883 /* import scaling vector (wrong formulation if we have 3D edges) */ 884 if (scaling && compute_Stilda) { 885 const PetscScalar *array; 886 PetscScalar *array2; 887 const PetscInt *idxs; 888 PetscInt i; 889 890 ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 891 ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr); 892 ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr); 893 ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr); 894 for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]]; 895 ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr); 896 ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr); 897 ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 898 deluxe = PETSC_FALSE; 899 } 900 901 /* size active schurs does not count any dirichlet or vertex dof on the interface */ 902 factor_workaround = PETSC_FALSE; 903 schur_has_vertices = PETSC_FALSE; 904 cum = n_I+size_active_schur; 905 if (sub_schurs->is_dir) { 906 const PetscInt* idxs; 907 PetscInt n_dir; 908 909 ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr); 910 ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 911 ierr = PetscArraycpy(all_local_idx_N+cum,idxs,n_dir);CHKERRQ(ierr); 912 ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 913 cum += n_dir; 914 factor_workaround = PETSC_TRUE; 915 } 916 /* include the primal vertices in the Schur complement */ 917 if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) { 918 PetscInt n_v; 919 920 ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr); 921 if (n_v) { 922 const PetscInt* idxs; 923 924 ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 925 ierr = PetscArraycpy(all_local_idx_N+cum,idxs,n_v);CHKERRQ(ierr); 926 ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 927 cum += n_v; 928 factor_workaround = PETSC_TRUE; 929 schur_has_vertices = PETSC_TRUE; 930 } 931 } 932 size_schur = cum - n_I; 933 ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr); 934 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 935 oldpin = sub_schurs->A->boundtocpu; 936 ierr = MatBindToCPU(sub_schurs->A,PETSC_TRUE);CHKERRQ(ierr); 937 #endif 938 if (cum == n) { 939 ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr); 940 ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr); 941 } else { 942 ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 943 } 944 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 945 ierr = MatBindToCPU(sub_schurs->A,oldpin);CHKERRQ(ierr); 946 #endif 947 ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr); 948 ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr); 949 950 /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory 951 this is a workaround */ 952 if (benign_n) { 953 Vec v,benign_AIIm1_ones; 954 ISLocalToGlobalMapping N_to_reor; 955 IS is_p0,is_p0_p; 956 PetscScalar *cs_AIB,*AIIm1_data; 957 PetscInt sizeA; 958 959 ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr); 960 ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr); 961 ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr); 962 ierr = ISDestroy(&is_p0);CHKERRQ(ierr); 963 ierr = MatCreateVecs(A,&v,&benign_AIIm1_ones);CHKERRQ(ierr); 964 ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr); 965 ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr); 966 ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr); 967 ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 968 ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 969 ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr); 970 /* compute colsum of A_IB restricted to pressures */ 971 for (i=0;i<benign_n;i++) { 972 const PetscScalar *array; 973 const PetscInt *idxs; 974 PetscInt j,nz; 975 976 ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr); 977 ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr); 978 ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 979 for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.; 980 ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 981 ierr = VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);CHKERRQ(ierr); 982 ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr); 983 ierr = VecResetArray(benign_AIIm1_ones);CHKERRQ(ierr); 984 ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 985 for (j=0;j<size_schur;j++) { 986 #if defined(PETSC_USE_COMPLEX) 987 cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz)); 988 #else 989 cs_AIB[i*size_schur + j] = array[j+n_I]/nz; 990 #endif 991 } 992 ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 993 } 994 ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 995 ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 996 ierr = VecDestroy(&v);CHKERRQ(ierr); 997 ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr); 998 ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr); 999 ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 1000 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 1001 ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr); 1002 ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr); 1003 ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr); 1004 } 1005 ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr); 1006 ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 1007 ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 1008 1009 /* for complexes, symmetric and hermitian at the same time implies null imaginary part */ 1010 use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric); 1011 1012 /* when using the benign subspace trick, the local Schur complements are SPD */ 1013 /* MKL_PARDISO does not handle well the computation of a Schur complement from a symmetric indefinite factorization 1014 Use LU and adapt pivoting perturbation (still, solution is not as accurate as with using MUMPS) */ 1015 if (benign_trick) { 1016 sub_schurs->is_posdef = PETSC_TRUE; 1017 ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,&flg);CHKERRQ(ierr); 1018 if (flg) use_cholesky = PETSC_FALSE; 1019 } 1020 1021 if (n_I) { 1022 IS is_schur; 1023 char stype[64]; 1024 PetscBool gpu = PETSC_FALSE; 1025 1026 if (use_cholesky) { 1027 ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 1028 } else { 1029 ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 1030 } 1031 ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr); 1032 #if defined(PETSC_HAVE_MKL_PARDISO) 1033 if (benign_trick) { ierr = MatMkl_PardisoSetCntl(F,10,10);CHKERRQ(ierr); } 1034 #endif 1035 /* subsets ordered last */ 1036 ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr); 1037 ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr); 1038 ierr = ISDestroy(&is_schur);CHKERRQ(ierr); 1039 1040 /* factorization step */ 1041 if (use_cholesky) { 1042 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 1043 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 1044 ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr); 1045 #endif 1046 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 1047 S_lower_triangular = PETSC_TRUE; 1048 } else { 1049 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 1050 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */ 1051 ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr); 1052 #endif 1053 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 1054 } 1055 ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr); 1056 1057 if (matl_dbg_viewer) { 1058 Mat S; 1059 IS is; 1060 1061 ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr); 1062 ierr = MatView(A,matl_dbg_viewer);CHKERRQ(ierr); 1063 ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr); 1064 ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr); 1065 ierr = MatView(S,matl_dbg_viewer);CHKERRQ(ierr); 1066 ierr = MatDestroy(&S);CHKERRQ(ierr); 1067 ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr); 1068 ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr); 1069 ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr); 1070 ierr = ISDestroy(&is);CHKERRQ(ierr); 1071 ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr); 1072 ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr); 1073 ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr); 1074 ierr = ISDestroy(&is);CHKERRQ(ierr); 1075 ierr = PetscObjectSetName((PetscObject)is_A_all,"IA");CHKERRQ(ierr); 1076 ierr = ISView(is_A_all,matl_dbg_viewer);CHKERRQ(ierr); 1077 } 1078 1079 /* get explicit Schur Complement computed during numeric factorization */ 1080 ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr); 1081 ierr = PetscStrncpy(stype,MATSEQDENSE,sizeof(stype));CHKERRQ(ierr); 1082 #if defined(PETSC_HAVE_CUDA) 1083 ierr = PetscObjectTypeCompareAny((PetscObject)A,&gpu,MATSEQAIJVIENNACL,MATSEQAIJCUSPARSE,"");CHKERRQ(ierr); 1084 #endif 1085 if (gpu) { 1086 ierr = PetscStrncpy(stype,MATSEQDENSECUDA,sizeof(stype));CHKERRQ(ierr); 1087 } 1088 ierr = PetscOptionsGetString(NULL,sub_schurs->prefix,"-sub_schurs_schur_mat_type",stype,sizeof(stype),NULL);CHKERRQ(ierr); 1089 ierr = MatConvert(S_all,stype,MAT_INPLACE_MATRIX,&S_all);CHKERRQ(ierr); 1090 ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 1091 ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 1092 ierr = MatGetType(S_all,&Stype);CHKERRQ(ierr); 1093 1094 /* we can reuse the solvers if we are not using the economic version */ 1095 reuse_solvers = (PetscBool)(reuse_solvers && !economic); 1096 factor_workaround = (PetscBool)(reuse_solvers && factor_workaround); 1097 if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur) 1098 reuse_solvers = factor_workaround = PETSC_FALSE; 1099 1100 solver_S = PETSC_TRUE; 1101 1102 /* update the Schur complement with the change of basis on the pressures */ 1103 if (benign_n) { 1104 const PetscScalar *cs_AIB; 1105 PetscScalar *S_data,*AIIm1_data; 1106 Mat S2 = NULL,S3 = NULL; /* dbg */ 1107 PetscScalar *S2_data,*S3_data; /* dbg */ 1108 Vec v,benign_AIIm1_ones; 1109 PetscInt sizeA; 1110 1111 ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr); 1112 ierr = MatCreateVecs(A,&v,&benign_AIIm1_ones);CHKERRQ(ierr); 1113 ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr); 1114 #if defined(PETSC_HAVE_MUMPS) 1115 ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr); 1116 #endif 1117 #if defined(PETSC_HAVE_MKL_PARDISO) 1118 ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr); 1119 #endif 1120 ierr = MatDenseGetArrayRead(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 1121 ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 1122 if (matl_dbg_viewer) { 1123 ierr = MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S2);CHKERRQ(ierr); 1124 ierr = MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S3);CHKERRQ(ierr); 1125 ierr = MatDenseGetArray(S2,&S2_data);CHKERRQ(ierr); 1126 ierr = MatDenseGetArray(S3,&S3_data);CHKERRQ(ierr); 1127 } 1128 for (i=0;i<benign_n;i++) { 1129 PetscScalar *array,sum = 0.,one = 1.,*sums; 1130 const PetscInt *idxs; 1131 PetscInt k,j,nz; 1132 PetscBLASInt B_k,B_n; 1133 1134 ierr = PetscCalloc1(benign_n,&sums);CHKERRQ(ierr); 1135 ierr = VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);CHKERRQ(ierr); 1136 ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr); 1137 ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr); 1138 ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr); 1139 ierr = VecResetArray(benign_AIIm1_ones);CHKERRQ(ierr); 1140 /* p0 dofs (eliminated) are excluded from the sums */ 1141 for (k=0;k<benign_n;k++) { 1142 ierr = ISGetLocalSize(is_p_r[k],&nz);CHKERRQ(ierr); 1143 ierr = ISGetIndices(is_p_r[k],&idxs);CHKERRQ(ierr); 1144 for (j=0;j<nz-1;j++) sums[k] -= AIIm1_data[idxs[j]+sizeA*i]; 1145 ierr = ISRestoreIndices(is_p_r[k],&idxs);CHKERRQ(ierr); 1146 } 1147 ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 1148 if (matl_dbg_viewer) { 1149 Vec vv; 1150 char name[16]; 1151 1152 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array+n_I,&vv);CHKERRQ(ierr); 1153 ierr = PetscSNPrintf(name,sizeof(name),"Pvs%D",i);CHKERRQ(ierr); 1154 ierr = PetscObjectSetName((PetscObject)vv,name);CHKERRQ(ierr); 1155 ierr = VecView(vv,matl_dbg_viewer);CHKERRQ(ierr); 1156 } 1157 /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */ 1158 /* cs_AIB already scaled by 1./nz */ 1159 B_k = 1; 1160 for (k=0;k<benign_n;k++) { 1161 sum = sums[k]; 1162 ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr); 1163 1164 if (PetscAbsScalar(sum) == 0.0) continue; 1165 if (k == i) { 1166 PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n)); 1167 if (matl_dbg_viewer) { 1168 PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n)); 1169 } 1170 } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */ 1171 sum /= 2.0; 1172 PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n)); 1173 if (matl_dbg_viewer) { 1174 PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n)); 1175 } 1176 } 1177 } 1178 sum = 1.; 1179 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)); 1180 if (matl_dbg_viewer) { 1181 PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S2_data,&B_n)); 1182 } 1183 ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr); 1184 /* set p0 entry of AIIm1_ones to zero */ 1185 ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr); 1186 ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 1187 for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.; 1188 ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr); 1189 ierr = PetscFree(sums);CHKERRQ(ierr); 1190 } 1191 ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr); 1192 if (matl_dbg_viewer) { 1193 ierr = MatDenseRestoreArray(S2,&S2_data);CHKERRQ(ierr); 1194 ierr = MatDenseRestoreArray(S3,&S3_data);CHKERRQ(ierr); 1195 } 1196 if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */ 1197 PetscInt k,j; 1198 for (k=0;k<size_schur;k++) { 1199 for (j=k;j<size_schur;j++) { 1200 S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]); 1201 } 1202 } 1203 } 1204 1205 /* restore defaults */ 1206 #if defined(PETSC_HAVE_MUMPS) 1207 ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr); 1208 #endif 1209 #if defined(PETSC_HAVE_MKL_PARDISO) 1210 ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr); 1211 #endif 1212 ierr = MatDenseRestoreArrayRead(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr); 1213 ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr); 1214 ierr = VecDestroy(&v);CHKERRQ(ierr); 1215 ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr); 1216 if (matl_dbg_viewer) { 1217 Mat S; 1218 1219 ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1220 ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr); 1221 ierr = PetscObjectSetName((PetscObject)S,"Sb");CHKERRQ(ierr); 1222 ierr = MatView(S,matl_dbg_viewer);CHKERRQ(ierr); 1223 ierr = MatDestroy(&S);CHKERRQ(ierr); 1224 ierr = PetscObjectSetName((PetscObject)S2,"S2P");CHKERRQ(ierr); 1225 ierr = MatView(S2,matl_dbg_viewer);CHKERRQ(ierr); 1226 ierr = PetscObjectSetName((PetscObject)S3,"S3P");CHKERRQ(ierr); 1227 ierr = MatView(S3,matl_dbg_viewer);CHKERRQ(ierr); 1228 ierr = PetscObjectSetName((PetscObject)cs_AIB_mat,"cs");CHKERRQ(ierr); 1229 ierr = MatView(cs_AIB_mat,matl_dbg_viewer);CHKERRQ(ierr); 1230 ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr); 1231 } 1232 ierr = MatDestroy(&S2);CHKERRQ(ierr); 1233 ierr = MatDestroy(&S3);CHKERRQ(ierr); 1234 } 1235 if (!reuse_solvers) { 1236 for (i=0;i<benign_n;i++) { 1237 ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr); 1238 } 1239 ierr = PetscFree(is_p_r);CHKERRQ(ierr); 1240 ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr); 1241 ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr); 1242 } 1243 } else { /* we can't use MatFactor when size_schur == size_of_the_problem */ 1244 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 1245 ierr = MatGetType(S_all,&Stype);CHKERRQ(ierr); 1246 reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */ 1247 factor_workaround = PETSC_FALSE; 1248 solver_S = PETSC_FALSE; 1249 } 1250 1251 if (reuse_solvers) { 1252 Mat A_II,Afake; 1253 Vec vec1_B; 1254 PCBDDCReuseSolvers msolv_ctx; 1255 PetscInt n_R; 1256 1257 if (sub_schurs->reuse_solver) { 1258 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 1259 } else { 1260 ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr); 1261 } 1262 msolv_ctx = sub_schurs->reuse_solver; 1263 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1264 ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 1265 msolv_ctx->F = F; 1266 ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr); 1267 /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */ 1268 { 1269 PetscScalar *array; 1270 PetscInt n; 1271 1272 ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr); 1273 ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr); 1274 ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr); 1275 ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr); 1276 } 1277 msolv_ctx->has_vertices = schur_has_vertices; 1278 1279 /* interior solver */ 1280 ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr); 1281 ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr); 1282 ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr); 1283 ierr = PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");CHKERRQ(ierr); 1284 ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr); 1285 ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr); 1286 ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr); 1287 ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr); 1288 1289 /* correction solver */ 1290 ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr); 1291 ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr); 1292 ierr = PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");CHKERRQ(ierr); 1293 ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr); 1294 ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr); 1295 ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr); 1296 ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr); 1297 1298 /* scatter and vecs for Schur complement solver */ 1299 ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr); 1300 ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr); 1301 if (!schur_has_vertices) { 1302 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr); 1303 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 1304 ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr); 1305 msolv_ctx->is_R = is_A_all; 1306 } else { 1307 IS is_B_all; 1308 const PetscInt* idxs; 1309 PetscInt dual,n_v,n; 1310 1311 ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr); 1312 dual = size_schur - n_v; 1313 ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr); 1314 ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr); 1315 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr); 1316 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr); 1317 ierr = ISDestroy(&is_B_all);CHKERRQ(ierr); 1318 ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr); 1319 ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr); 1320 ierr = ISDestroy(&is_B_all);CHKERRQ(ierr); 1321 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr); 1322 ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr); 1323 } 1324 ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr); 1325 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr); 1326 ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1327 ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1328 ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr); 1329 ierr = MatDestroy(&Afake);CHKERRQ(ierr); 1330 ierr = VecDestroy(&vec1_B);CHKERRQ(ierr); 1331 1332 /* communicate benign info to solver context */ 1333 if (benign_n) { 1334 PetscScalar *array; 1335 1336 msolv_ctx->benign_n = benign_n; 1337 msolv_ctx->benign_zerodiag_subs = is_p_r; 1338 ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr); 1339 msolv_ctx->benign_csAIB = cs_AIB_mat; 1340 ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr); 1341 ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr); 1342 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr); 1343 ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr); 1344 msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat; 1345 } 1346 } else { 1347 if (sub_schurs->reuse_solver) { 1348 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 1349 } 1350 ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr); 1351 } 1352 ierr = MatDestroy(&A);CHKERRQ(ierr); 1353 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 1354 1355 /* Work arrays */ 1356 ierr = PetscMalloc1(max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 1357 1358 /* S_Ej_all */ 1359 cum = cum2 = 0; 1360 ierr = MatDenseGetArrayRead(S_all,&rS_data);CHKERRQ(ierr); 1361 ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&SEj_arr);CHKERRQ(ierr); 1362 if (compute_Stilda) { 1363 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);CHKERRQ(ierr); 1364 } 1365 for (i=0;i<sub_schurs->n_subs;i++) { 1366 PetscInt j; 1367 1368 /* get S_E */ 1369 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1370 if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */ 1371 PetscInt k; 1372 for (k=0;k<subset_size;k++) { 1373 for (j=k;j<subset_size;j++) { 1374 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j]; 1375 work[j*subset_size+k] = PetscConj(rS_data[cum2+k*size_schur+j]); 1376 } 1377 } 1378 } else { /* just copy to workspace */ 1379 PetscInt k; 1380 for (k=0;k<subset_size;k++) { 1381 for (j=0;j<subset_size;j++) { 1382 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j]; 1383 } 1384 } 1385 } 1386 /* insert S_E values */ 1387 if (sub_schurs->change) { 1388 Mat change_sub,SEj,T; 1389 1390 /* change basis */ 1391 ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr); 1392 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1393 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1394 Mat T2; 1395 ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr); 1396 ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1397 ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1398 ierr = MatDestroy(&T2);CHKERRQ(ierr); 1399 } else { 1400 ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1401 } 1402 ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1403 ierr = MatDestroy(&T);CHKERRQ(ierr); 1404 ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr); 1405 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1406 } 1407 if (deluxe) { 1408 ierr = PetscArraycpy(SEj_arr,work,subset_size*subset_size);CHKERRQ(ierr); 1409 /* if adaptivity is requested, invert S_E blocks */ 1410 if (compute_Stilda) { 1411 Mat M; 1412 const PetscScalar *vals; 1413 PetscBool isdense,isdensecuda; 1414 1415 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&M);CHKERRQ(ierr); 1416 ierr = MatSetOption(M,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 1417 ierr = MatSetOption(M,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 1418 if (!PetscBTLookup(sub_schurs->is_edge,i)) { 1419 ierr = MatSetType(M,Stype);CHKERRQ(ierr); 1420 } 1421 ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&isdense);CHKERRQ(ierr); 1422 ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQDENSECUDA,&isdensecuda);CHKERRQ(ierr); 1423 if (use_cholesky) { 1424 ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr); 1425 } else { 1426 ierr = MatLUFactor(M,NULL,NULL,NULL);CHKERRQ(ierr); 1427 } 1428 if (isdense) { 1429 ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr); 1430 #if defined(PETSC_HAVE_CUDA) 1431 } else if (isdensecuda) { 1432 ierr = MatSeqDenseCUDAInvertFactors_Private(M);CHKERRQ(ierr); 1433 #endif 1434 } else SETERRQ1(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"Not implemented for type %s",Stype); 1435 ierr = MatDenseGetArrayRead(M,&vals);CHKERRQ(ierr); 1436 ierr = PetscArraycpy(SEjinv_arr,vals,subset_size*subset_size);CHKERRQ(ierr); 1437 ierr = MatDenseRestoreArrayRead(M,&vals);CHKERRQ(ierr); 1438 ierr = MatDestroy(&M);CHKERRQ(ierr); 1439 } 1440 } else if (compute_Stilda) { /* not using deluxe */ 1441 Mat SEj; 1442 Vec D; 1443 PetscScalar *array; 1444 1445 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1446 ierr = VecGetArray(Dall,&array);CHKERRQ(ierr); 1447 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr); 1448 ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr); 1449 ierr = VecShift(D,-1.);CHKERRQ(ierr); 1450 ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr); 1451 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1452 ierr = VecDestroy(&D);CHKERRQ(ierr); 1453 ierr = PetscArraycpy(SEj_arr,work,subset_size*subset_size);CHKERRQ(ierr); 1454 } 1455 cum += subset_size; 1456 cum2 += subset_size*(size_schur + 1); 1457 SEj_arr += subset_size*subset_size; 1458 if (SEjinv_arr) SEjinv_arr += subset_size*subset_size; 1459 } 1460 ierr = MatDenseRestoreArrayRead(S_all,&rS_data);CHKERRQ(ierr); 1461 ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&SEj_arr);CHKERRQ(ierr); 1462 if (compute_Stilda) { 1463 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);CHKERRQ(ierr); 1464 } 1465 if (solver_S) { 1466 ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1467 } 1468 1469 /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory 1470 however, preliminary tests indicate using GPUs is still faster in the solve phase */ 1471 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1472 if (reuse_solvers) { 1473 Mat St; 1474 MatFactorSchurStatus st; 1475 1476 flg = PETSC_FALSE; 1477 ierr = PetscOptionsGetBool(NULL,sub_schurs->prefix,"-sub_schurs_schur_pin_to_cpu",&flg,NULL);CHKERRQ(ierr); 1478 ierr = MatFactorGetSchurComplement(F,&St,&st);CHKERRQ(ierr); 1479 ierr = MatBindToCPU(St,flg);CHKERRQ(ierr); 1480 ierr = MatFactorRestoreSchurComplement(F,&St,st);CHKERRQ(ierr); 1481 } 1482 #endif 1483 1484 schur_factor = NULL; 1485 if (compute_Stilda && size_active_schur) { 1486 1487 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);CHKERRQ(ierr); 1488 if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */ 1489 ierr = PetscArraycpy(SEjinv_arr,work,size_schur*size_schur);CHKERRQ(ierr); 1490 } else { 1491 Mat S_all_inv=NULL; 1492 1493 if (solver_S) { 1494 /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1. 1495 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 */ 1496 if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */ 1497 PetscScalar *data; 1498 PetscInt nd = 0; 1499 1500 if (!use_potr) { 1501 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices"); 1502 } 1503 ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr); 1504 ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr); 1505 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1506 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1507 } 1508 1509 /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */ 1510 if (schur_has_vertices) { 1511 Mat M; 1512 PetscScalar *tdata; 1513 PetscInt nv = 0, news; 1514 1515 ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr); 1516 news = size_active_schur + nv; 1517 ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr); 1518 for (i=0;i<size_active_schur;i++) { 1519 ierr = PetscArraycpy(tdata+i*(news+1),data+i*(size_schur+1),size_active_schur-i);CHKERRQ(ierr); 1520 ierr = PetscArraycpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv);CHKERRQ(ierr); 1521 } 1522 for (i=0;i<nv;i++) { 1523 PetscInt k = i+size_active_schur; 1524 ierr = PetscArraycpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),nv-i);CHKERRQ(ierr); 1525 } 1526 1527 ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr); 1528 ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1529 ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr); 1530 /* save the factors */ 1531 cum = 0; 1532 ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr); 1533 for (i=0;i<size_active_schur;i++) { 1534 ierr = PetscArraycpy(schur_factor+cum,tdata+i*(news+1),size_active_schur-i);CHKERRQ(ierr); 1535 cum += size_active_schur - i; 1536 } 1537 for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)])); 1538 ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr); 1539 /* move back just the active dofs to the Schur complement */ 1540 for (i=0;i<size_active_schur;i++) { 1541 ierr = PetscArraycpy(data+i*size_schur,tdata+i*news,size_active_schur);CHKERRQ(ierr); 1542 } 1543 ierr = PetscFree(tdata);CHKERRQ(ierr); 1544 ierr = MatDestroy(&M);CHKERRQ(ierr); 1545 } else { /* we can factorize and invert just the activedofs */ 1546 Mat M; 1547 PetscScalar *aux; 1548 1549 ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr); 1550 for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)]; 1551 ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr); 1552 ierr = MatDenseSetLDA(M,size_schur);CHKERRQ(ierr); 1553 ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1554 ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr); 1555 ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr); 1556 ierr = MatDestroy(&M);CHKERRQ(ierr); 1557 ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr); 1558 ierr = MatZeroEntries(M);CHKERRQ(ierr); 1559 ierr = MatDestroy(&M);CHKERRQ(ierr); 1560 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr); 1561 ierr = MatDenseSetLDA(M,size_schur);CHKERRQ(ierr); 1562 ierr = MatZeroEntries(M);CHKERRQ(ierr); 1563 ierr = MatDestroy(&M);CHKERRQ(ierr); 1564 for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i]; 1565 ierr = PetscFree(aux);CHKERRQ(ierr); 1566 } 1567 ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr); 1568 } else { /* use MatFactor calls to invert S */ 1569 ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr); 1570 ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr); 1571 } 1572 } else { /* we need to invert explicitly since we are not using MatFactor for S */ 1573 ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr); 1574 S_all_inv = S_all; 1575 ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr); 1576 ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr); 1577 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1578 if (use_potr) { 1579 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr)); 1580 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 1581 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr)); 1582 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 1583 } else if (use_sytr) { 1584 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1585 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr); 1586 PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr)); 1587 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr); 1588 } else { 1589 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr)); 1590 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 1591 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1592 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 1593 } 1594 ierr = PetscLogFlops(1.0*size_schur*size_schur*size_schur);CHKERRQ(ierr); 1595 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1596 ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr); 1597 } 1598 /* S_Ej_tilda_all */ 1599 cum = cum2 = 0; 1600 ierr = MatDenseGetArrayRead(S_all_inv,&rS_data);CHKERRQ(ierr); 1601 for (i=0;i<sub_schurs->n_subs;i++) { 1602 PetscInt j; 1603 1604 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1605 /* get (St^-1)_E */ 1606 /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1 1607 will be properly accessed later during adaptive selection */ 1608 if (S_lower_triangular) { 1609 PetscInt k; 1610 if (sub_schurs->change) { 1611 for (k=0;k<subset_size;k++) { 1612 for (j=k;j<subset_size;j++) { 1613 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j]; 1614 work[j*subset_size+k] = work[k*subset_size+j]; 1615 } 1616 } 1617 } else { 1618 for (k=0;k<subset_size;k++) { 1619 for (j=k;j<subset_size;j++) { 1620 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j]; 1621 } 1622 } 1623 } 1624 } else { 1625 PetscInt k; 1626 for (k=0;k<subset_size;k++) { 1627 for (j=0;j<subset_size;j++) { 1628 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j]; 1629 } 1630 } 1631 } 1632 if (sub_schurs->change) { 1633 Mat change_sub,SEj,T; 1634 1635 /* change basis */ 1636 ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr); 1637 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr); 1638 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1639 Mat T2; 1640 ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr); 1641 ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1642 ierr = MatDestroy(&T2);CHKERRQ(ierr); 1643 ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1644 } else { 1645 ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr); 1646 } 1647 ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1648 ierr = MatDestroy(&T);CHKERRQ(ierr); 1649 /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */ 1650 ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr); 1651 ierr = MatDestroy(&SEj);CHKERRQ(ierr); 1652 } 1653 ierr = PetscArraycpy(SEjinv_arr,work,subset_size*subset_size);CHKERRQ(ierr); 1654 cum += subset_size; 1655 cum2 += subset_size*(size_schur + 1); 1656 SEjinv_arr += subset_size*subset_size; 1657 } 1658 ierr = MatDenseRestoreArrayRead(S_all_inv,&rS_data);CHKERRQ(ierr); 1659 if (solver_S) { 1660 if (schur_has_vertices) { 1661 ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr); 1662 } else { 1663 ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr); 1664 } 1665 } 1666 ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr); 1667 } 1668 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);CHKERRQ(ierr); 1669 1670 /* move back factors if needed */ 1671 if (schur_has_vertices) { 1672 Mat S_tmp; 1673 PetscInt nd = 0; 1674 1675 if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen"); 1676 ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr); 1677 if (use_potr) { 1678 PetscScalar *data; 1679 1680 ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr); 1681 ierr = PetscArrayzero(data,size_schur*size_schur);CHKERRQ(ierr); 1682 1683 if (S_lower_triangular) { 1684 cum = 0; 1685 for (i=0;i<size_active_schur;i++) { 1686 ierr = PetscArraycpy(data+i*(size_schur+1),schur_factor+cum,size_active_schur-i);CHKERRQ(ierr); 1687 cum += size_active_schur-i; 1688 } 1689 } else { 1690 ierr = PetscArraycpy(data,schur_factor,size_schur*size_schur);CHKERRQ(ierr); 1691 } 1692 if (sub_schurs->is_dir) { 1693 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1694 for (i=0;i<nd;i++) { 1695 data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i]; 1696 } 1697 } 1698 /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions, 1699 set the diagonal entry of the Schur factor to a very large value */ 1700 for (i=size_active_schur+nd;i<size_schur;i++) { 1701 data[i*(size_schur+1)] = infty; 1702 } 1703 ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr); 1704 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices"); 1705 ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr); 1706 } 1707 } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */ 1708 PetscScalar *data; 1709 PetscInt nd = 0; 1710 1711 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1712 ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr); 1713 } 1714 ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr); 1715 ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr); 1716 for (i=0;i<size_active_schur;i++) { 1717 ierr = PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);CHKERRQ(ierr); 1718 } 1719 for (i=size_active_schur+nd;i<size_schur;i++) { 1720 ierr = PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);CHKERRQ(ierr); 1721 data[i*(size_schur+1)] = infty; 1722 } 1723 ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr); 1724 ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1725 } 1726 ierr = PetscFree(work);CHKERRQ(ierr); 1727 ierr = PetscFree(schur_factor);CHKERRQ(ierr); 1728 ierr = VecDestroy(&Dall);CHKERRQ(ierr); 1729 } 1730 ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr); 1731 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 1732 ierr = MatDestroy(&A_BB);CHKERRQ(ierr); 1733 ierr = MatDestroy(&A_IB);CHKERRQ(ierr); 1734 ierr = MatDestroy(&A_BI);CHKERRQ(ierr); 1735 ierr = MatDestroy(&F);CHKERRQ(ierr); 1736 1737 ierr = PetscMalloc1(sub_schurs->n_subs,&nnz);CHKERRQ(ierr); 1738 for (i=0;i<sub_schurs->n_subs;i++) { 1739 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);CHKERRQ(ierr); 1740 } 1741 ierr = ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);CHKERRQ(ierr); 1742 ierr = MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr); 1743 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1744 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1745 if (compute_Stilda) { 1746 ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr); 1747 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1748 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1749 if (deluxe) { 1750 ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr); 1751 ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1752 ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1753 } 1754 } 1755 ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr); 1756 1757 /* Get local part of (\sum_j S_Ej) */ 1758 if (!sub_schurs->sum_S_Ej_all) { 1759 ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1760 } 1761 ierr = VecSet(gstash,0.0);CHKERRQ(ierr); 1762 ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&stasharray);CHKERRQ(ierr); 1763 ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr); 1764 ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1765 ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1766 ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&stasharray);CHKERRQ(ierr); 1767 ierr = VecResetArray(lstash);CHKERRQ(ierr); 1768 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&stasharray);CHKERRQ(ierr); 1769 ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr); 1770 ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1771 ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1772 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&stasharray);CHKERRQ(ierr); 1773 ierr = VecResetArray(lstash);CHKERRQ(ierr); 1774 1775 /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */ 1776 if (compute_Stilda) { 1777 ierr = VecSet(gstash,0.0);CHKERRQ(ierr); 1778 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);CHKERRQ(ierr); 1779 ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr); 1780 ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1781 ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1782 ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1783 ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1784 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);CHKERRQ(ierr); 1785 ierr = VecResetArray(lstash);CHKERRQ(ierr); 1786 if (deluxe) { 1787 ierr = VecSet(gstash,0.0);CHKERRQ(ierr); 1788 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);CHKERRQ(ierr); 1789 ierr = VecPlaceArray(lstash,stasharray);CHKERRQ(ierr); 1790 ierr = VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1791 ierr = VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1792 ierr = VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1793 ierr = VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1794 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);CHKERRQ(ierr); 1795 ierr = VecResetArray(lstash);CHKERRQ(ierr); 1796 } else { 1797 PetscScalar *array; 1798 PetscInt cum; 1799 1800 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr); 1801 cum = 0; 1802 for (i=0;i<sub_schurs->n_subs;i++) { 1803 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1804 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 1805 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1806 if (use_potr) { 1807 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 1808 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 1809 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 1810 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 1811 } else if (use_sytr) { 1812 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1813 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr); 1814 PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr)); 1815 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr); 1816 } else { 1817 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 1818 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 1819 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 1820 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 1821 } 1822 ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr); 1823 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1824 cum += subset_size*subset_size; 1825 } 1826 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr); 1827 ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 1828 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 1829 sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all; 1830 } 1831 } 1832 ierr = VecDestroy(&lstash);CHKERRQ(ierr); 1833 ierr = VecDestroy(&gstash);CHKERRQ(ierr); 1834 ierr = VecScatterDestroy(&sstash);CHKERRQ(ierr); 1835 1836 if (matl_dbg_viewer) { 1837 PetscInt cum; 1838 1839 if (sub_schurs->S_Ej_all) { 1840 ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr); 1841 ierr = MatView(sub_schurs->S_Ej_all,matl_dbg_viewer);CHKERRQ(ierr); 1842 } 1843 if (sub_schurs->sum_S_Ej_all) { 1844 ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr); 1845 ierr = MatView(sub_schurs->sum_S_Ej_all,matl_dbg_viewer);CHKERRQ(ierr); 1846 } 1847 if (sub_schurs->sum_S_Ej_inv_all) { 1848 ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr); 1849 ierr = MatView(sub_schurs->sum_S_Ej_inv_all,matl_dbg_viewer);CHKERRQ(ierr); 1850 } 1851 if (sub_schurs->sum_S_Ej_tilda_all) { 1852 ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr); 1853 ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,matl_dbg_viewer);CHKERRQ(ierr); 1854 } 1855 for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 1856 IS is; 1857 char name[16]; 1858 1859 ierr = PetscSNPrintf(name,sizeof(name),"IE%D",i);CHKERRQ(ierr); 1860 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 1861 ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr); 1862 ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr); 1863 ierr = ISView(is,matl_dbg_viewer);CHKERRQ(ierr); 1864 ierr = ISDestroy(&is);CHKERRQ(ierr); 1865 cum += subset_size; 1866 } 1867 } 1868 1869 /* free workspace */ 1870 ierr = PetscViewerDestroy(&matl_dbg_viewer);CHKERRQ(ierr); 1871 ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 1872 ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 1873 PetscFunctionReturn(0); 1874 } 1875 1876 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc) 1877 { 1878 IS *faces,*edges,*all_cc,vertices; 1879 PetscInt i,n_faces,n_edges,n_all_cc; 1880 PetscBool is_sorted,ispardiso,ismumps; 1881 PetscErrorCode ierr; 1882 1883 PetscFunctionBegin; 1884 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 1885 if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 1886 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 1887 if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 1888 1889 /* reset any previous data */ 1890 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 1891 1892 /* get index sets for faces and edges (already sorted by global ordering) */ 1893 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1894 n_all_cc = n_faces+n_edges; 1895 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 1896 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 1897 for (i=0;i<n_faces;i++) { 1898 if (copycc) { 1899 ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr); 1900 } else { 1901 ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr); 1902 all_cc[i] = faces[i]; 1903 } 1904 } 1905 for (i=0;i<n_edges;i++) { 1906 if (copycc) { 1907 ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr); 1908 } else { 1909 ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr); 1910 all_cc[n_faces+i] = edges[i]; 1911 } 1912 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 1913 } 1914 ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr); 1915 sub_schurs->is_vertices = vertices; 1916 ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 1917 sub_schurs->is_dir = NULL; 1918 ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr); 1919 1920 /* Determine if MatFactor can be used */ 1921 ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr); 1922 #if defined(PETSC_HAVE_MUMPS) 1923 ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,sizeof(sub_schurs->mat_solver_type));CHKERRQ(ierr); 1924 #elif defined(PETSC_HAVE_MKL_PARDISO) 1925 ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,sizeof(sub_schurs->mat_solver_type));CHKERRQ(ierr); 1926 #else 1927 ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,sizeof(sub_schurs->mat_solver_type));CHKERRQ(ierr); 1928 #endif 1929 #if defined(PETSC_USE_COMPLEX) 1930 sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */ 1931 #else 1932 sub_schurs->is_hermitian = PETSC_TRUE; 1933 #endif 1934 sub_schurs->is_posdef = PETSC_TRUE; 1935 sub_schurs->is_symmetric = PETSC_TRUE; 1936 sub_schurs->debug = PETSC_FALSE; 1937 sub_schurs->restrict_comm = PETSC_FALSE; 1938 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr); 1939 ierr = PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,sizeof(sub_schurs->mat_solver_type),NULL);CHKERRQ(ierr); 1940 ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr); 1941 ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr); 1942 ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr); 1943 ierr = PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);CHKERRQ(ierr); 1944 ierr = PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);CHKERRQ(ierr); 1945 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1946 ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMUMPS,&ismumps);CHKERRQ(ierr); 1947 ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,&ispardiso);CHKERRQ(ierr); 1948 sub_schurs->schur_explicit = (PetscBool)(ispardiso || ismumps); 1949 1950 /* for reals, symmetric and hermitian are synonims */ 1951 #if !defined(PETSC_USE_COMPLEX) 1952 sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian); 1953 sub_schurs->is_hermitian = sub_schurs->is_symmetric; 1954 #endif 1955 1956 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 1957 sub_schurs->is_I = is_I; 1958 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 1959 sub_schurs->is_B = is_B; 1960 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 1961 sub_schurs->l2gmap = graph->l2gmap; 1962 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 1963 sub_schurs->BtoNmap = BtoNmap; 1964 sub_schurs->n_subs = n_all_cc; 1965 sub_schurs->is_subs = all_cc; 1966 sub_schurs->S_Ej_all = NULL; 1967 sub_schurs->sum_S_Ej_all = NULL; 1968 sub_schurs->sum_S_Ej_inv_all = NULL; 1969 sub_schurs->sum_S_Ej_tilda_all = NULL; 1970 sub_schurs->is_Ej_all = NULL; 1971 PetscFunctionReturn(0); 1972 } 1973 1974 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 1975 { 1976 PCBDDCSubSchurs schurs_ctx; 1977 PetscErrorCode ierr; 1978 1979 PetscFunctionBegin; 1980 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 1981 schurs_ctx->n_subs = 0; 1982 *sub_schurs = schurs_ctx; 1983 PetscFunctionReturn(0); 1984 } 1985 1986 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 1987 { 1988 PetscInt i; 1989 PetscErrorCode ierr; 1990 1991 PetscFunctionBegin; 1992 if (!sub_schurs) PetscFunctionReturn(0); 1993 ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr); 1994 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 1995 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 1996 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 1997 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 1998 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 1999 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 2000 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 2001 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 2002 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 2003 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 2004 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 2005 ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr); 2006 ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr); 2007 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 2008 for (i=0;i<sub_schurs->n_subs;i++) { 2009 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 2010 } 2011 if (sub_schurs->n_subs) { 2012 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 2013 } 2014 if (sub_schurs->reuse_solver) { 2015 ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr); 2016 } 2017 ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr); 2018 if (sub_schurs->change) { 2019 for (i=0;i<sub_schurs->n_subs;i++) { 2020 ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr); 2021 ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr); 2022 } 2023 } 2024 ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr); 2025 ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr); 2026 sub_schurs->n_subs = 0; 2027 PetscFunctionReturn(0); 2028 } 2029 2030 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs) 2031 { 2032 PetscErrorCode ierr; 2033 2034 PetscFunctionBegin; 2035 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 2036 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 2037 PetscFunctionReturn(0); 2038 } 2039 2040 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 2041 { 2042 PetscInt i,j,n; 2043 PetscErrorCode ierr; 2044 2045 PetscFunctionBegin; 2046 n = 0; 2047 for (i=-n_prev;i<0;i++) { 2048 PetscInt start_dof = queue_tip[i]; 2049 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 2050 PetscInt dof = adjncy[j]; 2051 if (!PetscBTLookup(touched,dof)) { 2052 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 2053 queue_tip[n] = dof; 2054 n++; 2055 } 2056 } 2057 } 2058 *n_added = n; 2059 PetscFunctionReturn(0); 2060 } 2061