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 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCBDDCComputeExplicitSchur" 10 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 11 { 12 Mat B, C, D, Bd, Cd, AinvBd; 13 KSP ksp; 14 PC pc; 15 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 16 PetscReal fill = 2.0; 17 PetscMPIInt size; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 22 if (size != 1) { 23 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 24 } 25 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 26 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 27 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 28 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 29 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 30 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 31 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 32 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 33 if (!Bdense) { 34 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 35 } else { 36 Bd = B; 37 } 38 39 if (isLU || isILU || isCHOL) { 40 Mat fact; 41 42 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 43 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 44 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 45 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 46 } else { 47 Mat Ainvd; 48 49 ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 50 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 51 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 52 } 53 if (!Bdense & !issym) { 54 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 55 } 56 57 if (!issym) { 58 if (!Cdense) { 59 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 60 } else { 61 Cd = C; 62 } 63 ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 64 if (!Cdense) { 65 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 66 } 67 } else { 68 ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 69 if (!Bdense) { 70 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 71 } 72 } 73 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 74 75 if (D) { 76 Mat Dd; 77 PetscBool Ddense; 78 79 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 80 if (!Ddense) { 81 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 82 } else { 83 Dd = D; 84 } 85 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 86 if (!Ddense) { 87 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 88 } 89 } else { 90 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 91 } 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCBDDCSubSchursSetUp" 97 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool deluxe, PetscBool invert_Stildas, PetscBool use_edges, PetscBool use_faces) 98 { 99 Mat A_II,A_IB,A_BI,A_BB; 100 Mat AE_II,*AE_IE,*AE_EI,*AE_EE; 101 Mat S_all,global_schur_subsets,work_mat; 102 Mat S_Ej_tilda_all,S_Ej_inv_all; 103 ISLocalToGlobalMapping l2gmap_subsets; 104 IS is_I,*is_subset_B,temp_is; 105 PetscInt *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N; 106 PetscInt i,subset_size,max_subset_size_seq; 107 PetscInt extra,local_size,global_size; 108 PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 109 PetscScalar *Bwork; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 /* get Schur complement matrices */ 114 if (!sub_schurs->use_mumps) { 115 if (compute_Stilda) { 116 SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS"); 117 } 118 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr); 119 ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B, 120 sub_schurs->n_subs,&AE_IE, 121 sub_schurs->n_subs,&AE_EI, 122 sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr); 123 } else { 124 is_subset_B = NULL; 125 AE_IE = NULL; 126 AE_EI = NULL; 127 AE_EE = NULL; 128 } 129 130 /* determine interior problems */ 131 ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 132 ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 133 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 134 PetscBT touched; 135 const PetscInt* idx_B; 136 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 137 138 if (xadj == NULL || adjncy == NULL) { 139 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 140 } 141 /* get sizes */ 142 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 143 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 144 145 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 146 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 147 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 148 149 /* all boundary dofs must be skipped when adding layers */ 150 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 151 for (j=0;j<n_B;j++) { 152 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 153 } 154 ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 155 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 156 157 /* add prescribed number of layers of dofs */ 158 n_local_dofs = n_B; 159 n_prev_added = n_B; 160 for (layer=0;layer<nlayers;layer++) { 161 PetscInt n_added; 162 if (n_local_dofs == n_I+n_B) break; 163 if (n_local_dofs > n_I+n_B) { 164 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); 165 } 166 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 167 n_prev_added = n_added; 168 n_local_dofs += n_added; 169 if (!n_added) break; 170 } 171 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 172 173 /* IS for I layer dofs in original numbering */ 174 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_I_layer);CHKERRQ(ierr); 175 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 176 ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr); 177 /* IS for I layer dofs in I numbering */ 178 if (!sub_schurs->use_mumps) { 179 ISLocalToGlobalMapping ItoNmap; 180 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 181 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr); 182 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 183 184 /* II block */ 185 ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 186 } 187 } else { 188 PetscInt n_I; 189 190 /* IS for I dofs in original numbering */ 191 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 192 sub_schurs->is_I_layer = sub_schurs->is_I; 193 194 /* IS for I dofs in I numbering (strided 1) */ 195 if (!sub_schurs->use_mumps) { 196 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 197 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 198 199 /* II block is the same */ 200 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 201 AE_II = A_II; 202 } 203 } 204 /* Get info on subset sizes and sum of all subsets sizes */ 205 max_subset_size_seq = 0; 206 local_size = 0; 207 for (i=0;i<sub_schurs->n_subs_seq;i++) { 208 PetscInt j = sub_schurs->index_sequential[i]; 209 ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr); 210 max_subset_size_seq = PetscMax(subset_size,max_subset_size_seq); 211 local_size += subset_size; 212 } 213 214 /* Work arrays for local indices */ 215 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 216 extra = 0; 217 if (sub_schurs->use_mumps) { 218 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr); 219 } 220 ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 221 if (extra) { 222 const PetscInt *idxs; 223 ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 224 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 225 ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 226 } 227 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 228 229 /* Get local indices in local numbering */ 230 local_size = 0; 231 for (i=0;i<sub_schurs->n_subs_seq;i++) { 232 PetscInt j; 233 const PetscInt *idxs; 234 235 PetscInt local_problem_index = sub_schurs->index_sequential[i]; 236 ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr); 237 ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr); 238 /* subset indices in local numbering */ 239 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 240 ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr); 241 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 242 local_size += subset_size; 243 } 244 245 S_all = NULL; 246 sub_schurs->is_hermitian = PETSC_FALSE; 247 sub_schurs->is_posdef = PETSC_FALSE; 248 if (sub_schurs->A) { 249 ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr); 250 if (sub_schurs->is_hermitian) { 251 PetscScalar val; 252 Vec vec1,vec2; 253 254 ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr); 255 ierr = VecSetRandom(vec1,NULL); 256 ierr = VecCopy(vec1,vec2);CHKERRQ(ierr); 257 ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr); 258 ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr); 259 if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE; 260 ierr = VecDestroy(&vec1);CHKERRQ(ierr); 261 ierr = VecDestroy(&vec2);CHKERRQ(ierr); 262 } 263 } 264 265 if (sub_schurs->n_subs) { 266 if (!sub_schurs->use_mumps) { 267 /* subsets in original and boundary numbering */ 268 for (i=0;i<sub_schurs->n_subs;i++) { 269 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr); 270 } 271 272 /* EE block */ 273 for (i=0;i<sub_schurs->n_subs;i++) { 274 ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr); 275 } 276 277 /* IE block */ 278 for (i=0;i<sub_schurs->n_subs;i++) { 279 ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr); 280 } 281 282 /* EI block */ 283 for (i=0;i<sub_schurs->n_subs;i++) { 284 if (sub_schurs->is_hermitian) { 285 ierr = MatCreateTranspose(AE_IE[i],&AE_EI[i]);CHKERRQ(ierr); 286 } else { 287 ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr); 288 } 289 } 290 291 /* setup Schur complements on subset */ 292 for (i=0;i<sub_schurs->n_subs;i++) { 293 ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 294 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 295 ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr); 296 ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr); 297 ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr); 298 if (AE_II == A_II) { /* we can reuse the same ksp */ 299 KSP ksp; 300 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 301 ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr); 302 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 303 KSP origksp,schurksp; 304 PC origpc,schurpc; 305 KSPType ksp_type; 306 PCType pc_type; 307 PetscInt n_internal; 308 309 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 310 ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr); 311 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 312 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 313 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 314 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 315 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 316 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 317 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 318 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 319 MatSolverPackage solver=NULL; 320 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 321 if (solver) { 322 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 323 } 324 } 325 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 326 } 327 } 328 /* free */ 329 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 330 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 331 for (i=0;i<sub_schurs->n_subs;i++) { 332 ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr); 333 } 334 ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr); 335 #if defined(PETSC_HAVE_MUMPS) 336 } else { 337 Mat A,F; 338 IS is_A_all; 339 PetscInt *idxs_schur,n_I; 340 341 /* get working mat */ 342 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr); 343 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 344 ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 345 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 346 347 if (n_I) { 348 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 349 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 350 } else { 351 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 352 } 353 354 /* subsets ordered last */ 355 ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr); 356 for (i=0;i<local_size;i++) { 357 idxs_schur[i] = n_I+i+1; 358 } 359 ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr); 360 ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 361 362 /* factorization step */ 363 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 364 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 365 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 366 } else { 367 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 368 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 369 } 370 371 /* get explicit Schur Complement computed during numeric factorization */ 372 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 373 ierr = MatDestroy(&F);CHKERRQ(ierr); 374 } else { 375 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 376 } 377 ierr = MatDestroy(&A);CHKERRQ(ierr); 378 #endif 379 } 380 } else { 381 ierr = PetscFree(nnz);CHKERRQ(ierr); 382 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 383 } 384 385 if (!sub_schurs->n_subs_seq_g) { 386 PetscFunctionReturn(0); 387 } 388 389 /* subset indices in local boundary numbering */ 390 if (!sub_schurs->is_Ej_all) { 391 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 392 if (subset_size != local_size) { 393 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 394 } 395 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 396 } 397 398 /* Local matrix of all local Schur on subsets transposed */ 399 if (!sub_schurs->S_Ej_all) { 400 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 401 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 402 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 403 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 404 } else { 405 ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr); 406 } 407 408 S_Ej_tilda_all = 0; 409 S_Ej_inv_all = 0; 410 Bwork = NULL; 411 pivots = NULL; 412 if (sub_schurs->n_subs && deluxe && compute_Stilda && !sub_schurs->is_hermitian) { /* workspace needed only for GETRI */ 413 PetscScalar lwork; 414 415 B_lwork = -1; 416 ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr); 417 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 418 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr)); 419 ierr = PetscFPTrapPop();CHKERRQ(ierr); 420 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 421 ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr); 422 ierr = PetscMalloc2(B_lwork,&Bwork,max_subset_size_seq,&pivots);CHKERRQ(ierr); 423 } 424 425 ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 426 if (!sub_schurs->use_mumps) { 427 if (sub_schurs->n_subs_seq) { 428 PetscScalar *work; 429 PetscInt *dummy_idx; 430 431 /* Work arrays */ 432 ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr); 433 434 /* Loop on local problems end compute Schur complements explicitly */ 435 local_size = 0; 436 for (i=0;i<sub_schurs->n_subs_seq;i++) { 437 Mat S_Ej_expl; 438 PetscInt j,lpi = sub_schurs->index_sequential[i]; 439 PetscBool Sdense; 440 441 ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr); 442 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 443 ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 444 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 445 if (Sdense) { 446 for (j=0;j<subset_size;j++) { 447 dummy_idx[j]=local_size+j; 448 } 449 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 450 } else { 451 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 452 } 453 local_size += subset_size; 454 ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr); 455 sub_schurs->S_Ej[lpi] = S_Ej_expl; 456 } 457 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 458 } 459 } else { 460 PetscInt *dummy_idx,n_all; 461 PetscScalar *work; 462 463 if (compute_Stilda) { 464 ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr); 465 ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 466 ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 467 ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 468 if (deluxe) { 469 ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr); 470 ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 471 ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr); 472 ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr); 473 } 474 } 475 476 if (S_all) { /* multilevel */ 477 ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr); 478 } else { 479 n_all = 0; 480 } 481 local_size = 0; 482 483 /* Work arrays */ 484 if (compute_Stilda) { 485 ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr); 486 } else { 487 ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr); 488 } 489 490 for (i=0;i<sub_schurs->n_subs;i++) { 491 IS is_E; 492 PetscScalar *vals; 493 PetscInt j; 494 495 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 496 for (j=0;j<subset_size;j++) { 497 dummy_idx[j]=local_size+j; 498 } 499 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr); 500 if (sub_schurs->S_Ej[i]) { 501 ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 502 } else { 503 ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 504 } 505 ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 506 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 507 ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 508 509 if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) { 510 Mat Stilda; 511 512 ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr); 513 if (sub_schurs->n_subs == 1) { 514 ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr); 515 Stilda = sub_schurs->S_Ej[i]; 516 } else { 517 KSP ksp; 518 PC pc; 519 Mat S_EF,S_FE,S_FF,Stilda_impl; 520 IS is_F; 521 PetscScalar eps=1.e-8; 522 PetscInt nc,cum; 523 PetscBool chop=PETSC_FALSE; 524 525 ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr); 526 nc = n_all - subset_size; 527 cum = 0; 528 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,nc,work+cum,&S_FF);CHKERRQ(ierr); 529 cum += nc*nc; 530 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,subset_size,work+cum,&S_FE);CHKERRQ(ierr); 531 cum += nc*subset_size; 532 ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr); 533 ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_REUSE_MATRIX,&S_FE);CHKERRQ(ierr); 534 if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */ 535 ierr = MatCreateTranspose(S_FE,&S_EF);CHKERRQ(ierr); 536 } else { 537 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,nc,work+cum,&S_EF);CHKERRQ(ierr); 538 cum += nc*subset_size; 539 ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_REUSE_MATRIX,&S_EF);CHKERRQ(ierr); 540 } 541 ierr = ISDestroy(&is_F);CHKERRQ(ierr); 542 if (chop) { 543 ierr = MatChop(S_FF,eps);CHKERRQ(ierr); 544 ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr); 545 } 546 ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr); 547 ierr = MatDestroy(&S_FF);CHKERRQ(ierr); 548 ierr = MatDestroy(&S_FE);CHKERRQ(ierr); 549 ierr = MatDestroy(&S_EF);CHKERRQ(ierr); 550 ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr); 551 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 552 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 553 if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 554 ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr); 555 } else { 556 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 557 } 558 if (!chop) { 559 ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr); 560 } else { 561 ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr); 562 } 563 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 564 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work+cum,&Stilda);CHKERRQ(ierr); 565 ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr); 566 ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr); 567 } 568 /* 569 PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB); 570 ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 571 PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF); 572 */ 573 ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr); 574 if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */ 575 PetscScalar *vals2; 576 577 ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr); 578 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 579 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 580 if (invert_Stildas) { /* Stildas can be singular */ 581 if (!sub_schurs->is_hermitian) { 582 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr)); 583 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 584 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 585 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 586 } else { 587 PetscInt j,k; 588 589 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr)); 590 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 591 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr)); 592 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 593 for (j=0;j<B_N;j++) { 594 for (k=j+1;k<B_N;k++) { 595 vals[k*B_N+j] = vals[j*B_N+k]; 596 } 597 } 598 } 599 } 600 if (!sub_schurs->is_hermitian) { 601 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr)); 602 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 603 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 604 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 605 } else { 606 PetscInt j,k; 607 608 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals2,&B_N,&B_ierr)); 609 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 610 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals2,&B_N,&B_ierr)); 611 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 612 for (j=0;j<B_N;j++) { 613 for (k=j+1;k<B_N;k++) { 614 vals2[k*B_N+j] = vals2[j*B_N+k]; 615 } 616 } 617 } 618 ierr = PetscFPTrapPop();CHKERRQ(ierr); 619 ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr); 620 ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr); 621 } 622 ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 623 ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr); 624 ierr = MatDestroy(&Stilda);CHKERRQ(ierr); 625 } 626 ierr = ISDestroy(&is_E);CHKERRQ(ierr); 627 local_size += subset_size; 628 } 629 ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 630 } 631 ierr = PetscFree(nnz);CHKERRQ(ierr); 632 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 633 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 634 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 635 if (S_Ej_tilda_all) { 636 ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 637 ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 638 } 639 if (S_Ej_inv_all) { 640 ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 641 ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 642 } 643 644 /* Global matrix of all assembled Schur on subsets */ 645 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->l2gmap,local_size,all_local_idx_N+extra,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr); 646 ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 647 ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 648 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 649 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 650 ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 651 ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 652 ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 653 ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 654 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 655 /* Get local part of (\sum_j S_Ej) */ 656 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 657 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr); 658 if (sub_schurs->sum_S_Ej_all) { 659 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 660 } else { 661 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 662 } 663 664 if (S_Ej_tilda_all) { 665 ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr); 666 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 667 if (sub_schurs->sum_S_Ej_tilda_all) { 668 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 669 } else { 670 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 671 } 672 } 673 if (S_Ej_inv_all) { 674 PetscInt cum; 675 PetscScalar *array,*array2; 676 ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr); 677 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 678 if (sub_schurs->sum_S_Ej_inv_all) { 679 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 680 } else { 681 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 682 } 683 /* invert blocks of sum_S_Ej_inv_all */ 684 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 685 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 686 cum = 0; 687 for (i=0;i<sub_schurs->n_subs;i++) { 688 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 689 if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) { 690 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 691 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 692 if (!sub_schurs->is_hermitian) { 693 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 694 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 695 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 696 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 697 } else { 698 PetscInt j,k; 699 700 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 701 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 702 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 703 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 704 for (j=0;j<B_N;j++) { 705 for (k=j+1;k<B_N;k++) { 706 array[k*B_N+j+cum] = array[j*B_N+k+cum]; 707 } 708 } 709 } 710 if (invert_Stildas) { /* Stildas can be singular */ 711 if (!sub_schurs->is_hermitian) { 712 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr)); 713 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 714 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 715 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 716 } else { 717 PetscInt j,k; 718 719 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr)); 720 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 721 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr)); 722 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 723 for (j=0;j<B_N;j++) { 724 for (k=j+1;k<B_N;k++) { 725 array2[k*B_N+j+cum] = array2[j*B_N+k+cum]; 726 } 727 } 728 } 729 } 730 ierr = PetscFPTrapPop();CHKERRQ(ierr); 731 } 732 cum += subset_size*subset_size; 733 } 734 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 735 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 736 } 737 738 ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 739 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 740 ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr); 741 ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr); 742 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 743 ierr = ISDestroy(&temp_is);CHKERRQ(ierr); 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "PCBDDCSubSchursInit" 749 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold) 750 { 751 IS *faces,*edges,*all_cc,vertices; 752 PetscInt *index_sequential,*index_parallel; 753 PetscInt *auxlocal_sequential,*auxlocal_parallel; 754 PetscInt *auxglobal_sequential,*auxglobal_parallel; 755 PetscInt *auxmapping; 756 PetscInt i,max_subset_size; 757 PetscInt n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems; 758 PetscInt n_faces,n_edges,n_all_cc; 759 PetscBool is_sorted; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 764 if (!is_sorted) { 765 SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 766 } 767 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 768 if (!is_sorted) { 769 SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 770 } 771 772 /* reset any previous data */ 773 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 774 775 /* get index sets for faces and edges */ 776 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 777 n_all_cc = n_faces+n_edges; 778 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 779 ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 780 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 781 for (i=0;i<n_faces;i++) { 782 all_cc[i] = faces[i]; 783 } 784 for (i=0;i<n_edges;i++) { 785 all_cc[n_faces+i] = edges[i]; 786 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 787 } 788 ierr = PetscFree(faces);CHKERRQ(ierr); 789 ierr = PetscFree(edges);CHKERRQ(ierr); 790 791 /* map interface's subsets */ 792 max_subset_size = 0; 793 for (i=0;i<n_all_cc;i++) { 794 PetscInt subset_size; 795 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 796 max_subset_size = PetscMax(max_subset_size,subset_size); 797 } 798 ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr); 799 ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential, 800 graph->ncc,&auxlocal_parallel);CHKERRQ(ierr); 801 ierr = PetscMalloc2(graph->ncc,&index_sequential, 802 graph->ncc,&index_parallel);CHKERRQ(ierr); 803 804 /* if threshold is negative uses all sequential problems (possibly using MUMPS) */ 805 sub_schurs->use_mumps = PETSC_FALSE; 806 if (seqthreshold < 0) { 807 seqthreshold = max_subset_size; 808 #if defined(PETSC_HAVE_MUMPS) 809 sub_schurs->use_mumps = !!A; 810 #endif 811 } 812 813 814 /* determine which problem has to be solved in parallel or sequentially */ 815 n_local_sequential_problems = 0; 816 n_local_parallel_problems = 0; 817 for (i=0;i<n_all_cc;i++) { 818 PetscInt subset_size,j,min_loc = 0; 819 const PetscInt *idxs; 820 821 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 822 ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr); 823 ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr); 824 for (j=1;j<subset_size;j++) { 825 if (auxmapping[j]<auxmapping[min_loc]) { 826 min_loc = j; 827 } 828 } 829 if (subset_size > seqthreshold) { 830 index_parallel[n_local_parallel_problems] = i; 831 auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc]; 832 n_local_parallel_problems++; 833 } else { 834 index_sequential[n_local_sequential_problems] = i; 835 auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc]; 836 n_local_sequential_problems++; 837 } 838 ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr); 839 } 840 841 /* Number parallel problems */ 842 auxglobal_parallel = 0; 843 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr); 844 845 /* Number sequential problems */ 846 auxglobal_sequential = 0; 847 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr); 848 849 /* update info in sub_schurs */ 850 if (A) { 851 PetscBool isseqaij; 852 853 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 854 if (isseqaij) { 855 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 856 sub_schurs->A = A; 857 } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */ 858 ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr); 859 } 860 } 861 ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr); 862 sub_schurs->S = S; 863 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 864 sub_schurs->is_I = is_I; 865 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 866 sub_schurs->is_B = is_B; 867 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 868 sub_schurs->l2gmap = graph->l2gmap; 869 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 870 sub_schurs->BtoNmap = BtoNmap; 871 sub_schurs->n_subs_seq = n_local_sequential_problems; 872 sub_schurs->n_subs_par = n_local_parallel_problems; 873 sub_schurs->n_subs_seq_g = n_sequential_problems; 874 sub_schurs->n_subs_par_g = n_parallel_problems; 875 sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par; 876 sub_schurs->is_subs = all_cc; 877 if (!sub_schurs->use_mumps) { /* for adaptive selection */ 878 for (i=0;i<sub_schurs->n_subs;i++) { 879 ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 880 } 881 } 882 sub_schurs->is_Ej_com = vertices; 883 sub_schurs->index_sequential = index_sequential; 884 sub_schurs->index_parallel = index_parallel; 885 sub_schurs->auxglobal_sequential = auxglobal_sequential; 886 sub_schurs->auxglobal_parallel = auxglobal_parallel; 887 888 889 /* allocate space for schur complements */ 890 ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr); 891 sub_schurs->S_Ej_all = NULL; 892 sub_schurs->sum_S_Ej_all = NULL; 893 sub_schurs->sum_S_Ej_inv_all = NULL; 894 sub_schurs->sum_S_Ej_tilda_all = NULL; 895 sub_schurs->is_Ej_all = NULL; 896 897 /* free workspace */ 898 ierr = PetscFree(auxmapping);CHKERRQ(ierr); 899 ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr); 900 901 PetscFunctionReturn(0); 902 } 903 904 #undef __FUNCT__ 905 #define __FUNCT__ "PCBDDCSubSchursCreate" 906 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 907 { 908 PCBDDCSubSchurs schurs_ctx; 909 PetscErrorCode ierr; 910 911 PetscFunctionBegin; 912 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 913 schurs_ctx->n_subs = 0; 914 *sub_schurs = schurs_ctx; 915 PetscFunctionReturn(0); 916 } 917 918 #undef __FUNCT__ 919 #define __FUNCT__ "PCBDDCSubSchursDestroy" 920 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 921 { 922 PetscErrorCode ierr; 923 924 PetscFunctionBegin; 925 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 926 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "PCBDDCSubSchursReset" 932 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 933 { 934 PetscInt i; 935 PetscErrorCode ierr; 936 937 PetscFunctionBegin; 938 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 939 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 940 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 941 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 942 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 943 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 944 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 945 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 946 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 947 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 948 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 949 ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr); 950 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 951 ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 952 for (i=0;i<sub_schurs->n_subs;i++) { 953 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 954 ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 955 } 956 ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 957 if (sub_schurs->n_subs) { 958 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 959 ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr); 960 } 961 ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr); 962 ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr); 963 ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr); 964 sub_schurs->n_subs = 0; 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 970 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 971 { 972 PetscInt i,j,n; 973 PetscErrorCode ierr; 974 975 PetscFunctionBegin; 976 n = 0; 977 for (i=-n_prev;i<0;i++) { 978 PetscInt start_dof = queue_tip[i]; 979 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 980 PetscInt dof = adjncy[j]; 981 if (!PetscBTLookup(touched,dof)) { 982 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 983 queue_tip[n] = dof; 984 n++; 985 } 986 } 987 } 988 *n_added = n; 989 PetscFunctionReturn(0); 990 } 991