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