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