134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h> 234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 334a97f8cSStefano Zampini 43202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 53202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S); 63202ece2SStefano Zampini 73202ece2SStefano Zampini #undef __FUNCT__ 83202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur" 93202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S) 103202ece2SStefano Zampini { 113202ece2SStefano Zampini Mat B, C, D, Bd, Cd, AinvBd; 123202ece2SStefano Zampini KSP ksp; 133202ece2SStefano Zampini PC pc; 143202ece2SStefano Zampini PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 153202ece2SStefano Zampini PetscReal fill = 2.0; 163202ece2SStefano Zampini PetscMPIInt size; 173202ece2SStefano Zampini PetscErrorCode ierr; 183202ece2SStefano Zampini 193202ece2SStefano Zampini PetscFunctionBegin; 203202ece2SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 213202ece2SStefano Zampini if (size != 1) { 223202ece2SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 233202ece2SStefano Zampini } 243202ece2SStefano Zampini ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 253202ece2SStefano Zampini ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 263202ece2SStefano Zampini ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 273202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 283202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 293202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 303202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 313202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 323202ece2SStefano Zampini if (!Bdense) { 333202ece2SStefano Zampini ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 343202ece2SStefano Zampini } else { 353202ece2SStefano Zampini Bd = B; 363202ece2SStefano Zampini } 373202ece2SStefano Zampini 383202ece2SStefano Zampini if (isLU || isILU || isCHOL) { 393202ece2SStefano Zampini Mat fact; 403202ece2SStefano Zampini 413202ece2SStefano Zampini ierr = KSPSetUp(ksp);CHKERRQ(ierr); 423202ece2SStefano Zampini ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 433202ece2SStefano Zampini ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 443202ece2SStefano Zampini ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 453202ece2SStefano Zampini } else { 463202ece2SStefano Zampini Mat Ainvd; 473202ece2SStefano Zampini 483202ece2SStefano Zampini ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 493202ece2SStefano Zampini ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 503202ece2SStefano Zampini ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 513202ece2SStefano Zampini } 523202ece2SStefano Zampini if (!Bdense) { 533202ece2SStefano Zampini ierr = MatDestroy(&Bd);CHKERRQ(ierr); 543202ece2SStefano Zampini } 553202ece2SStefano Zampini if (!Cdense) { 563202ece2SStefano Zampini ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 573202ece2SStefano Zampini } else { 583202ece2SStefano Zampini Cd = C; 593202ece2SStefano Zampini } 603202ece2SStefano Zampini 613202ece2SStefano Zampini ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr); 623202ece2SStefano Zampini ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 633202ece2SStefano Zampini if (!Cdense) { 643202ece2SStefano Zampini ierr = MatDestroy(&Cd);CHKERRQ(ierr); 653202ece2SStefano Zampini } 663202ece2SStefano Zampini 673202ece2SStefano Zampini if (D) { 683202ece2SStefano Zampini Mat Dd; 693202ece2SStefano Zampini PetscBool Ddense; 703202ece2SStefano Zampini 713202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 723202ece2SStefano Zampini if (!Ddense) { 733202ece2SStefano Zampini ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 743202ece2SStefano Zampini } else { 753202ece2SStefano Zampini Dd = D; 763202ece2SStefano Zampini } 773202ece2SStefano Zampini ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 783202ece2SStefano Zampini if (!Ddense) { 793202ece2SStefano Zampini ierr = MatDestroy(&Dd);CHKERRQ(ierr); 803202ece2SStefano Zampini } 813202ece2SStefano Zampini } else { 823202ece2SStefano Zampini ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 833202ece2SStefano Zampini } 843202ece2SStefano Zampini PetscFunctionReturn(0); 853202ece2SStefano Zampini } 8634a97f8cSStefano Zampini 8734a97f8cSStefano Zampini #undef __FUNCT__ 881580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp" 891580ed26SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers) 90b1b3d7a2SStefano Zampini { 91b1b3d7a2SStefano Zampini Mat A_II,A_IB,A_BI,A_BB; 92b1b3d7a2SStefano Zampini Mat AE_II,*AE_IE,*AE_EI,*AE_EE; 93883469d8SStefano Zampini Mat S_all,global_schur_subsets,work_mat; 94*a1337663SStefano Zampini Mat S_Ej_tilda_all; 955db18549SStefano Zampini ISLocalToGlobalMapping l2gmap_subsets; 965db18549SStefano Zampini IS is_I,*is_subset_B,temp_is; 97883469d8SStefano Zampini PetscInt *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N; 985db18549SStefano Zampini PetscInt i,subset_size,max_subset_size; 99883469d8SStefano Zampini PetscInt extra,local_size,global_size; 100*a1337663SStefano Zampini PetscBool compute_Stilda=PETSC_FALSE; 101b1b3d7a2SStefano Zampini PetscErrorCode ierr; 102b1b3d7a2SStefano Zampini 103b1b3d7a2SStefano Zampini PetscFunctionBegin; 104b1b3d7a2SStefano Zampini /* allocate space for schur complements */ 1052b973097SStefano Zampini ierr = PetscMalloc5(sub_schurs->n_subs,&sub_schurs->is_AEj_B, 106b1b3d7a2SStefano Zampini sub_schurs->n_subs,&sub_schurs->S_Ej, 1072b973097SStefano Zampini sub_schurs->n_subs,&sub_schurs->S_Ej_tilda, 108b1b3d7a2SStefano Zampini sub_schurs->n_subs,&sub_schurs->work1, 109b1b3d7a2SStefano Zampini sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr); 110b1b3d7a2SStefano Zampini 111b1b3d7a2SStefano Zampini /* get Schur complement matrices */ 112883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 113b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr); 114b1b3d7a2SStefano Zampini ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B, 115b1b3d7a2SStefano Zampini sub_schurs->n_subs,&AE_IE, 116b1b3d7a2SStefano Zampini sub_schurs->n_subs,&AE_EI, 117b1b3d7a2SStefano Zampini sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr); 118b1b3d7a2SStefano Zampini } 119b1b3d7a2SStefano Zampini 120b1b3d7a2SStefano Zampini /* determine interior problems */ 121b1b3d7a2SStefano Zampini if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */ 122b1b3d7a2SStefano Zampini PetscBT touched; 123b1b3d7a2SStefano Zampini const PetscInt* idx_B; 124b1b3d7a2SStefano Zampini PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 125b1b3d7a2SStefano Zampini 126b1b3d7a2SStefano Zampini /* get sizes */ 127b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 128b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 129b1b3d7a2SStefano Zampini 130b1b3d7a2SStefano Zampini ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 131b1b3d7a2SStefano Zampini ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 132b1b3d7a2SStefano Zampini ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 133b1b3d7a2SStefano Zampini 134b1b3d7a2SStefano Zampini /* all boundary dofs must be skipped when adding layers */ 135b1b3d7a2SStefano Zampini ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 136b1b3d7a2SStefano Zampini for (j=0;j<n_B;j++) { 137b1b3d7a2SStefano Zampini ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 138b1b3d7a2SStefano Zampini } 139b1b3d7a2SStefano Zampini ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 140b1b3d7a2SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 141b1b3d7a2SStefano Zampini 142b1b3d7a2SStefano Zampini /* add prescribed number of layers of dofs */ 143b1b3d7a2SStefano Zampini n_local_dofs = n_B; 144b1b3d7a2SStefano Zampini n_prev_added = n_B; 145b1b3d7a2SStefano Zampini for (layer=0;layer<nlayers;layer++) { 146b1b3d7a2SStefano Zampini PetscInt n_added; 147b1b3d7a2SStefano Zampini if (n_local_dofs == n_I+n_B) break; 148b1b3d7a2SStefano Zampini if (n_local_dofs > n_I+n_B) { 149b1b3d7a2SStefano Zampini 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); 150b1b3d7a2SStefano Zampini } 151b1b3d7a2SStefano Zampini ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 152b1b3d7a2SStefano Zampini n_prev_added = n_added; 153b1b3d7a2SStefano Zampini n_local_dofs += n_added; 154b1b3d7a2SStefano Zampini if (!n_added) break; 155b1b3d7a2SStefano Zampini } 156b1b3d7a2SStefano Zampini ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 157b1b3d7a2SStefano Zampini 158883469d8SStefano Zampini /* IS for I layer dofs in original numbering */ 15968270318SStefano Zampini 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); 160b1b3d7a2SStefano Zampini ierr = PetscFree(local_numbering);CHKERRQ(ierr); 16168270318SStefano Zampini ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr); 162883469d8SStefano Zampini /* IS for I layer dofs in I numbering */ 163883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 164b1b3d7a2SStefano Zampini ISLocalToGlobalMapping ItoNmap; 165b1b3d7a2SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 16668270318SStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr); 167b1b3d7a2SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 168b1b3d7a2SStefano Zampini 169b1b3d7a2SStefano Zampini /* II block */ 170b1b3d7a2SStefano Zampini ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 171b1b3d7a2SStefano Zampini } 172b1b3d7a2SStefano Zampini } else { 173b1b3d7a2SStefano Zampini PetscInt n_I; 174b1b3d7a2SStefano Zampini 175b1b3d7a2SStefano Zampini /* IS for I dofs in original numbering */ 176b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 17768270318SStefano Zampini sub_schurs->is_I_layer = sub_schurs->is_I; 178b1b3d7a2SStefano Zampini 179b1b3d7a2SStefano Zampini /* IS for I dofs in I numbering (strided 1) */ 180883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 181b1b3d7a2SStefano Zampini ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 182b1b3d7a2SStefano Zampini ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 183b1b3d7a2SStefano Zampini 184b1b3d7a2SStefano Zampini /* II block is the same */ 185b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 186b1b3d7a2SStefano Zampini AE_II = A_II; 187b1b3d7a2SStefano Zampini } 188b1b3d7a2SStefano Zampini } 189b1b3d7a2SStefano Zampini 190b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 191b1b3d7a2SStefano Zampini ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 192b1b3d7a2SStefano Zampini ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 193883469d8SStefano Zampini } 194883469d8SStefano Zampini 195883469d8SStefano Zampini /* Get info on subset sizes and sum of all subsets sizes */ 196883469d8SStefano Zampini max_subset_size = 0; 197883469d8SStefano Zampini local_size = 0; 198883469d8SStefano Zampini for (i=0;i<sub_schurs->n_subs_seq;i++) { 199883469d8SStefano Zampini PetscInt j = sub_schurs->index_sequential[i]; 200883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr); 201883469d8SStefano Zampini max_subset_size = PetscMax(subset_size,max_subset_size); 202883469d8SStefano Zampini local_size += subset_size; 203883469d8SStefano Zampini } 204883469d8SStefano Zampini 205883469d8SStefano Zampini /* Work arrays for local indices */ 206883469d8SStefano Zampini ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 207883469d8SStefano Zampini extra = 0; 208883469d8SStefano Zampini if (sub_schurs->use_mumps) { 209883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr); 210883469d8SStefano Zampini } 211883469d8SStefano Zampini ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 212883469d8SStefano Zampini if (extra) { 213883469d8SStefano Zampini const PetscInt *idxs; 214883469d8SStefano Zampini ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 215883469d8SStefano Zampini ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 216883469d8SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 217883469d8SStefano Zampini } 218883469d8SStefano Zampini ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 219883469d8SStefano Zampini 220883469d8SStefano Zampini /* Get local indices in local numbering */ 221883469d8SStefano Zampini local_size = 0; 222883469d8SStefano Zampini for (i=0;i<sub_schurs->n_subs_seq;i++) { 223883469d8SStefano Zampini PetscInt j; 224883469d8SStefano Zampini const PetscInt *idxs; 225883469d8SStefano Zampini 226883469d8SStefano Zampini PetscInt local_problem_index = sub_schurs->index_sequential[i]; 227883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr); 228883469d8SStefano Zampini ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr); 229883469d8SStefano Zampini /* subset indices in local numbering */ 230883469d8SStefano Zampini ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 231883469d8SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr); 232883469d8SStefano Zampini for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 233883469d8SStefano Zampini local_size += subset_size; 234883469d8SStefano Zampini } 235883469d8SStefano Zampini 236883469d8SStefano Zampini S_all = NULL; 237883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 238883469d8SStefano Zampini /* subsets in original and boundary numbering */ 239883469d8SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 2405db18549SStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr); 241b1b3d7a2SStefano Zampini } 242b1b3d7a2SStefano Zampini 243b1b3d7a2SStefano Zampini /* EE block */ 244b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 245b1b3d7a2SStefano Zampini ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr); 246b1b3d7a2SStefano Zampini } 247b1b3d7a2SStefano Zampini /* IE block */ 248b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 249b1b3d7a2SStefano Zampini ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr); 250b1b3d7a2SStefano Zampini } 251b1b3d7a2SStefano Zampini /* EI block */ 252b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 253b1b3d7a2SStefano Zampini ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr); 254b1b3d7a2SStefano Zampini } 255b1b3d7a2SStefano Zampini 256b1b3d7a2SStefano Zampini /* setup Schur complements on subset */ 257b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 258b1b3d7a2SStefano Zampini ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 259b1b3d7a2SStefano Zampini if (AE_II == A_II) { /* we can reuse the same ksp */ 260b1b3d7a2SStefano Zampini KSP ksp; 261b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 262b1b3d7a2SStefano Zampini ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr); 263b1b3d7a2SStefano Zampini } else { /* build new ksp object which inherits ksp and pc types from the original one */ 264b1b3d7a2SStefano Zampini KSP origksp,schurksp; 265b1b3d7a2SStefano Zampini PC origpc,schurpc; 266b1b3d7a2SStefano Zampini KSPType ksp_type; 267b1b3d7a2SStefano Zampini PCType pc_type; 268b1b3d7a2SStefano Zampini PetscInt n_internal; 269b1b3d7a2SStefano Zampini 270b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 271b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr); 272b1b3d7a2SStefano Zampini ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 273b1b3d7a2SStefano Zampini ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 274b1b3d7a2SStefano Zampini ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 275b1b3d7a2SStefano Zampini ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 276b1b3d7a2SStefano Zampini ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 277b1b3d7a2SStefano Zampini ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 278b1b3d7a2SStefano Zampini ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 279b1b3d7a2SStefano Zampini if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 280b1b3d7a2SStefano Zampini MatSolverPackage solver=NULL; 281b1b3d7a2SStefano Zampini ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 282b1b3d7a2SStefano Zampini if (solver) { 283b1b3d7a2SStefano Zampini ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 284b1b3d7a2SStefano Zampini } 285b1b3d7a2SStefano Zampini } 286b1b3d7a2SStefano Zampini ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 287b1b3d7a2SStefano Zampini } 288b1b3d7a2SStefano Zampini } 289b1b3d7a2SStefano Zampini /* free */ 290b1b3d7a2SStefano Zampini ierr = ISDestroy(&is_I);CHKERRQ(ierr); 291b1b3d7a2SStefano Zampini ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 292b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 293b1b3d7a2SStefano Zampini ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr); 294b1b3d7a2SStefano Zampini ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr); 295b1b3d7a2SStefano Zampini ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr); 296b1b3d7a2SStefano Zampini ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr); 297b1b3d7a2SStefano Zampini } 298b1b3d7a2SStefano Zampini ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr); 299883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 300883469d8SStefano Zampini } else { 301883469d8SStefano Zampini Mat A,F; 302883469d8SStefano Zampini IS is_A_all; 303883469d8SStefano Zampini PetscBool is_symmetric; 304883469d8SStefano Zampini PetscInt *idxs_schur,n_I; 305883469d8SStefano Zampini 306883469d8SStefano Zampini /* get working mat */ 307883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr); 308883469d8SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 309883469d8SStefano Zampini ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr); 310883469d8SStefano Zampini ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 311883469d8SStefano Zampini 312883469d8SStefano Zampini ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr); 313883469d8SStefano Zampini if (is_symmetric) { 314883469d8SStefano Zampini ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 315883469d8SStefano Zampini } else { 316883469d8SStefano Zampini ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 317883469d8SStefano Zampini } 318883469d8SStefano Zampini 319883469d8SStefano Zampini /* subsets ordered last */ 320883469d8SStefano Zampini ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr); 321883469d8SStefano Zampini for (i=0;i<local_size;i++) { 322883469d8SStefano Zampini idxs_schur[i] = n_I+i+1; 323883469d8SStefano Zampini } 324883469d8SStefano Zampini ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr); 325883469d8SStefano Zampini ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 326883469d8SStefano Zampini 327883469d8SStefano Zampini /* factorization step */ 328883469d8SStefano Zampini if (is_symmetric) { 329883469d8SStefano Zampini ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 330883469d8SStefano Zampini ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 331883469d8SStefano Zampini } else { 332883469d8SStefano Zampini ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 333883469d8SStefano Zampini ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 334883469d8SStefano Zampini } 335883469d8SStefano Zampini 336883469d8SStefano Zampini /* get explicit Schur Complement computed during numeric factorization */ 337883469d8SStefano Zampini ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 338883469d8SStefano Zampini 339883469d8SStefano Zampini /* free workspace */ 340883469d8SStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 341883469d8SStefano Zampini ierr = MatDestroy(&F);CHKERRQ(ierr); 342883469d8SStefano Zampini 343883469d8SStefano Zampini #endif 344b1b3d7a2SStefano Zampini } 3455db18549SStefano Zampini 3469221af80SStefano Zampini /* TODO: just for compatibility with the previous version, needs to be fixed */ 3479221af80SStefano Zampini for (i=0;i<sub_schurs->n_subs_par;i++) { 3489221af80SStefano Zampini PetscInt j = sub_schurs->index_parallel[i]; 3499221af80SStefano Zampini ierr = MatCreateVecs(sub_schurs->S_Ej[j],&sub_schurs->work1[j],&sub_schurs->work2[j]);CHKERRQ(ierr); 3509221af80SStefano Zampini } 3519221af80SStefano Zampini for (i=0;i<sub_schurs->n_subs_seq;i++) { 3529221af80SStefano Zampini sub_schurs->work1[sub_schurs->index_sequential[i]] = 0; 3539221af80SStefano Zampini sub_schurs->work2[sub_schurs->index_sequential[i]] = 0; 3549221af80SStefano Zampini } 3559221af80SStefano Zampini 356*a1337663SStefano Zampini /* Stildas are not computed without mumps */ 357*a1337663SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 358*a1337663SStefano Zampini sub_schurs->S_Ej_tilda[i] = 0; 359*a1337663SStefano Zampini } 3605db18549SStefano Zampini if (!sub_schurs->n_subs_seq_g) { 3619221af80SStefano Zampini sub_schurs->S_Ej_all = 0; 3629221af80SStefano Zampini sub_schurs->sum_S_Ej_all = 0; 3639221af80SStefano Zampini sub_schurs->is_Ej_all = 0; 3645db18549SStefano Zampini PetscFunctionReturn(0); 3655db18549SStefano Zampini } 3665db18549SStefano Zampini 3675db18549SStefano Zampini /* subset indices in local boundary numbering */ 368883469d8SStefano Zampini ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 3695db18549SStefano Zampini if (subset_size != local_size) { 3705db18549SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 3715db18549SStefano Zampini } 3725db18549SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 3735db18549SStefano Zampini 3743202ece2SStefano Zampini /* Local matrix of all local Schur on subsets transposed */ 3755db18549SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 3765db18549SStefano Zampini ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 3775db18549SStefano Zampini ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 3785db18549SStefano Zampini ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 379*a1337663SStefano Zampini 380*a1337663SStefano Zampini sub_schurs->sum_S_Ej_tilda_all = 0; 381*a1337663SStefano Zampini S_Ej_tilda_all = 0; 3825db18549SStefano Zampini 383883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 3845db18549SStefano Zampini PetscScalar *fill_vals; 3855db18549SStefano Zampini PetscInt *dummy_idx; 3865db18549SStefano Zampini 3875db18549SStefano Zampini /* Work arrays */ 3883202ece2SStefano Zampini ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 3895db18549SStefano Zampini 3903202ece2SStefano Zampini /* Loop on local problems end compute Schur complements explicitly */ 3915db18549SStefano Zampini local_size = 0; 3925db18549SStefano Zampini for (i=0;i<sub_schurs->n_subs_seq;i++) { 3933202ece2SStefano Zampini Mat S_Ej_expl; 3943202ece2SStefano Zampini PetscInt j,lpi = sub_schurs->index_sequential[i]; 3953202ece2SStefano Zampini PetscBool Sdense; 3965db18549SStefano Zampini 3973202ece2SStefano Zampini ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr); 3983202ece2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr); 3993202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 4003202ece2SStefano Zampini if (Sdense) { 4015db18549SStefano Zampini for (j=0;j<subset_size;j++) { 4025db18549SStefano Zampini dummy_idx[j]=local_size+j; 4035db18549SStefano Zampini } 4043202ece2SStefano Zampini ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 4055db18549SStefano Zampini ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr); 4063202ece2SStefano Zampini ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 4073202ece2SStefano Zampini } else { 4083202ece2SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 4095db18549SStefano Zampini } 4103202ece2SStefano Zampini local_size += subset_size; 41165d8bf0aSStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr); 41265d8bf0aSStefano Zampini sub_schurs->S_Ej[lpi] = S_Ej_expl; 4133202ece2SStefano Zampini } 4142b973097SStefano Zampini /* Stildas are not computed without mumps */ 4152b973097SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 4162b973097SStefano Zampini sub_schurs->S_Ej_tilda[i] = 0; 4172b973097SStefano Zampini } 4183202ece2SStefano Zampini ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 419883469d8SStefano Zampini } else { 420883469d8SStefano Zampini PetscInt *dummy_idx,n_all; 421*a1337663SStefano Zampini PetscBool is_symmetric; 422883469d8SStefano Zampini 423*a1337663SStefano Zampini compute_Stilda=PETSC_TRUE; 424*a1337663SStefano Zampini if (compute_Stilda) { 425*a1337663SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr); 426*a1337663SStefano Zampini ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 427*a1337663SStefano Zampini ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 428*a1337663SStefano Zampini ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 429*a1337663SStefano Zampini } 430883469d8SStefano Zampini /* Work arrays */ 431883469d8SStefano Zampini ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 432883469d8SStefano Zampini 43365d8bf0aSStefano Zampini /* Work arrays */ 43465d8bf0aSStefano Zampini ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr); 43565d8bf0aSStefano Zampini local_size = 0; 43665d8bf0aSStefano Zampini ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr); 43765d8bf0aSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 43865d8bf0aSStefano Zampini IS is_E; 43965d8bf0aSStefano Zampini PetscScalar *vals; 44065d8bf0aSStefano Zampini PetscInt j; 44165d8bf0aSStefano Zampini 44265d8bf0aSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[i],&subset_size);CHKERRQ(ierr); 443*a1337663SStefano Zampini for (j=0;j<subset_size;j++) { 444*a1337663SStefano Zampini dummy_idx[j]=local_size+j; 445*a1337663SStefano Zampini } 44665d8bf0aSStefano Zampini ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr); 447*a1337663SStefano Zampini ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 448*a1337663SStefano Zampini ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 449*a1337663SStefano Zampini ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 450*a1337663SStefano Zampini ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 451*a1337663SStefano Zampini 45265d8bf0aSStefano Zampini if (compute_Stilda) { 45365d8bf0aSStefano Zampini Mat Stilda; 45465d8bf0aSStefano Zampini if (sub_schurs->n_subs == 1) { 45565d8bf0aSStefano Zampini ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr); 45665d8bf0aSStefano Zampini Stilda = sub_schurs->S_Ej[i]; 45765d8bf0aSStefano Zampini } else { 45865d8bf0aSStefano Zampini KSP ksp; 45965d8bf0aSStefano Zampini PC pc; 46065d8bf0aSStefano Zampini Mat S_EF,S_FE,S_FF,Stilda_impl; 46165d8bf0aSStefano Zampini IS is_F; 4622b973097SStefano Zampini PetscScalar eps=1.e-8; 4632b973097SStefano Zampini PetscBool chop=PETSC_FALSE; 46465d8bf0aSStefano Zampini 46565d8bf0aSStefano Zampini ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr); 46665d8bf0aSStefano Zampini ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr); 46765d8bf0aSStefano Zampini ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr); 46865d8bf0aSStefano Zampini ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr); 469*a1337663SStefano Zampini ierr = ISDestroy(&is_F);CHKERRQ(ierr); 4702b973097SStefano Zampini if (chop) { 4712b973097SStefano Zampini ierr = MatChop(S_FF,eps);CHKERRQ(ierr); 4722b973097SStefano Zampini ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr); 4732b973097SStefano Zampini } 474*a1337663SStefano Zampini ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr); 475*a1337663SStefano Zampini ierr = MatDestroy(&S_FF);CHKERRQ(ierr); 476*a1337663SStefano Zampini ierr = MatDestroy(&S_FE);CHKERRQ(ierr); 477*a1337663SStefano Zampini ierr = MatDestroy(&S_EF);CHKERRQ(ierr); 47865d8bf0aSStefano Zampini ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr); 47965d8bf0aSStefano Zampini ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 48065d8bf0aSStefano Zampini ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 48165d8bf0aSStefano Zampini if (is_symmetric) { 48265d8bf0aSStefano Zampini ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr); 48365d8bf0aSStefano Zampini } else { 48465d8bf0aSStefano Zampini ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 48565d8bf0aSStefano Zampini } 4862b973097SStefano Zampini if (!chop) { 48765d8bf0aSStefano Zampini ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr); 4882b973097SStefano Zampini } else { 4892b973097SStefano Zampini ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr); 4902b973097SStefano Zampini } 49165d8bf0aSStefano Zampini ierr = KSPSetUp(ksp);CHKERRQ(ierr); 49265d8bf0aSStefano Zampini ierr = PCBDDCComputeExplicitSchur(Stilda_impl,&Stilda);CHKERRQ(ierr); 49365d8bf0aSStefano Zampini ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr); 49465d8bf0aSStefano Zampini } 4952b973097SStefano Zampini /* 4962b973097SStefano Zampini PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB); 4972b973097SStefano Zampini ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 4982b973097SStefano Zampini PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF); 4992b973097SStefano Zampini */ 5002b973097SStefano Zampini sub_schurs->S_Ej_tilda[i]= Stilda; 501*a1337663SStefano Zampini ierr = MatDenseGetArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr); 502*a1337663SStefano Zampini ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 503*a1337663SStefano Zampini ierr = MatDenseRestoreArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr); 5042b973097SStefano Zampini } else { 5052b973097SStefano Zampini sub_schurs->S_Ej_tilda[i] = 0; 50665d8bf0aSStefano Zampini } 50765d8bf0aSStefano Zampini ierr = ISDestroy(&is_E);CHKERRQ(ierr); 508883469d8SStefano Zampini local_size += subset_size; 509883469d8SStefano Zampini } 510883469d8SStefano Zampini ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 5115db18549SStefano Zampini } 512*a1337663SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 513*a1337663SStefano Zampini ierr = MatDestroy(&S_all);CHKERRQ(ierr); 5145db18549SStefano Zampini ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5155db18549SStefano Zampini ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 516*a1337663SStefano Zampini if (compute_Stilda) { 517*a1337663SStefano Zampini ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 518*a1337663SStefano Zampini ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 519*a1337663SStefano Zampini } 520*a1337663SStefano Zampini 5215db18549SStefano Zampini /* Global matrix of all assembled Schur on subsets */ 522883469d8SStefano Zampini 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); 5235db18549SStefano Zampini ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 5245db18549SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 5255db18549SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 5265db18549SStefano Zampini ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 5275db18549SStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 5285db18549SStefano Zampini /* Get local part of (\sum_j S_Ej) */ 529883469d8SStefano Zampini ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 530883469d8SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr); 531883469d8SStefano Zampini ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 5323202ece2SStefano Zampini 533*a1337663SStefano Zampini if (compute_Stilda) { 534*a1337663SStefano Zampini ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr); 535*a1337663SStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 536*a1337663SStefano Zampini ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 5373202ece2SStefano Zampini } 5383202ece2SStefano Zampini 539*a1337663SStefano Zampini ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 540*a1337663SStefano Zampini ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr); 5413202ece2SStefano Zampini ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 5425db18549SStefano Zampini ierr = ISDestroy(&temp_is);CHKERRQ(ierr); 543b1b3d7a2SStefano Zampini PetscFunctionReturn(0); 544b1b3d7a2SStefano Zampini } 545b1b3d7a2SStefano Zampini 546b1b3d7a2SStefano Zampini #undef __FUNCT__ 547b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit" 5485db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold) 549b1b3d7a2SStefano Zampini { 550b1b3d7a2SStefano Zampini IS *faces,*edges,*all_cc; 551b1b3d7a2SStefano Zampini PetscInt *index_sequential,*index_parallel; 552b1b3d7a2SStefano Zampini PetscInt *auxlocal_sequential,*auxlocal_parallel; 553b1b3d7a2SStefano Zampini PetscInt *auxglobal_sequential,*auxglobal_parallel; 554b1b3d7a2SStefano Zampini PetscInt *auxmapping;//,*idxs; 555b1b3d7a2SStefano Zampini PetscInt i,max_subset_size; 556b1b3d7a2SStefano Zampini PetscInt n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems; 557b1b3d7a2SStefano Zampini PetscInt n_faces,n_edges,n_all_cc; 558b1b3d7a2SStefano Zampini PetscBool is_sorted; 559b1b3d7a2SStefano Zampini PetscErrorCode ierr; 560b1b3d7a2SStefano Zampini 561b1b3d7a2SStefano Zampini PetscFunctionBegin; 562b1b3d7a2SStefano Zampini ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 563b1b3d7a2SStefano Zampini if (!is_sorted) { 564b1b3d7a2SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 565b1b3d7a2SStefano Zampini } 566b1b3d7a2SStefano Zampini ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 567b1b3d7a2SStefano Zampini if (!is_sorted) { 568b1b3d7a2SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 569b1b3d7a2SStefano Zampini } 570b1b3d7a2SStefano Zampini 571b1b3d7a2SStefano Zampini /* reset any previous data */ 572b1b3d7a2SStefano Zampini ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 573b1b3d7a2SStefano Zampini 574b1b3d7a2SStefano Zampini /* get index sets for faces and edges */ 575b1b3d7a2SStefano Zampini ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr); 576b1b3d7a2SStefano Zampini n_all_cc = n_faces+n_edges; 577b1b3d7a2SStefano Zampini ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 578b1b3d7a2SStefano Zampini for (i=0;i<n_faces;i++) { 579b1b3d7a2SStefano Zampini all_cc[i] = faces[i]; 580b1b3d7a2SStefano Zampini } 581b1b3d7a2SStefano Zampini for (i=0;i<n_edges;i++) { 582b1b3d7a2SStefano Zampini all_cc[n_faces+i] = edges[i]; 583b1b3d7a2SStefano Zampini } 584b1b3d7a2SStefano Zampini ierr = PetscFree(faces);CHKERRQ(ierr); 585b1b3d7a2SStefano Zampini ierr = PetscFree(edges);CHKERRQ(ierr); 586b1b3d7a2SStefano Zampini 587b1b3d7a2SStefano Zampini /* map interface's subsets */ 588b1b3d7a2SStefano Zampini max_subset_size = 0; 589b1b3d7a2SStefano Zampini for (i=0;i<n_all_cc;i++) { 590b1b3d7a2SStefano Zampini PetscInt subset_size; 591b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 592b1b3d7a2SStefano Zampini max_subset_size = PetscMax(max_subset_size,subset_size); 593b1b3d7a2SStefano Zampini } 594b1b3d7a2SStefano Zampini ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr); 595b1b3d7a2SStefano Zampini ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential, 596b1b3d7a2SStefano Zampini graph->ncc,&auxlocal_parallel);CHKERRQ(ierr); 597b1b3d7a2SStefano Zampini ierr = PetscMalloc2(graph->ncc,&index_sequential, 598b1b3d7a2SStefano Zampini graph->ncc,&index_parallel);CHKERRQ(ierr); 599b1b3d7a2SStefano Zampini 600883469d8SStefano Zampini /* if threshold is negative uses all sequential problems (possibly using MUMPS) */ 601883469d8SStefano Zampini sub_schurs->use_mumps = PETSC_FALSE; 602883469d8SStefano Zampini if (seqthreshold < 0) { 603883469d8SStefano Zampini seqthreshold = max_subset_size; 604883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 605883469d8SStefano Zampini sub_schurs->use_mumps = !!A; 606883469d8SStefano Zampini #endif 607883469d8SStefano Zampini } 608883469d8SStefano Zampini 609b1b3d7a2SStefano Zampini 610b1b3d7a2SStefano Zampini /* determine which problem has to be solved in parallel or sequentially */ 611b1b3d7a2SStefano Zampini n_local_sequential_problems = 0; 612b1b3d7a2SStefano Zampini n_local_parallel_problems = 0; 613b1b3d7a2SStefano Zampini for (i=0;i<n_all_cc;i++) { 614b1b3d7a2SStefano Zampini PetscInt subset_size,j,min_loc = 0; 615b1b3d7a2SStefano Zampini const PetscInt *idxs; 616b1b3d7a2SStefano Zampini 617b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 618b1b3d7a2SStefano Zampini ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr); 619b1b3d7a2SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr); 620b1b3d7a2SStefano Zampini for (j=1;j<subset_size;j++) { 621b1b3d7a2SStefano Zampini if (auxmapping[j]<auxmapping[min_loc]) { 622b1b3d7a2SStefano Zampini min_loc = j; 623b1b3d7a2SStefano Zampini } 624b1b3d7a2SStefano Zampini } 625b1b3d7a2SStefano Zampini if (subset_size > seqthreshold) { 626b1b3d7a2SStefano Zampini index_parallel[n_local_parallel_problems] = i; 627b1b3d7a2SStefano Zampini auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc]; 628b1b3d7a2SStefano Zampini n_local_parallel_problems++; 629b1b3d7a2SStefano Zampini } else { 630b1b3d7a2SStefano Zampini index_sequential[n_local_sequential_problems] = i; 631b1b3d7a2SStefano Zampini auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc]; 632b1b3d7a2SStefano Zampini n_local_sequential_problems++; 633b1b3d7a2SStefano Zampini } 634b1b3d7a2SStefano Zampini ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr); 635b1b3d7a2SStefano Zampini } 636b1b3d7a2SStefano Zampini 637b1b3d7a2SStefano Zampini /* Number parallel problems */ 638b1b3d7a2SStefano Zampini auxglobal_parallel = 0; 639b1b3d7a2SStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr); 640b1b3d7a2SStefano Zampini 641b1b3d7a2SStefano Zampini /* Number sequential problems */ 642b1b3d7a2SStefano Zampini auxglobal_sequential = 0; 643b1b3d7a2SStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr); 644b1b3d7a2SStefano Zampini 645b1b3d7a2SStefano Zampini /* update info in sub_schurs */ 646883469d8SStefano Zampini if (sub_schurs->use_mumps && A) { 6471e9c79c2SStefano Zampini ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6481e9c79c2SStefano Zampini sub_schurs->A = A; 6491e9c79c2SStefano Zampini } 650b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr); 651b1b3d7a2SStefano Zampini sub_schurs->S = S; 652b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 653b1b3d7a2SStefano Zampini sub_schurs->is_I = is_I; 654b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 655b1b3d7a2SStefano Zampini sub_schurs->is_B = is_B; 6565db18549SStefano Zampini ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 6575db18549SStefano Zampini sub_schurs->l2gmap = graph->l2gmap; 6585db18549SStefano Zampini ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 6595db18549SStefano Zampini sub_schurs->BtoNmap = BtoNmap; 660b1b3d7a2SStefano Zampini sub_schurs->n_subs_seq = n_local_sequential_problems; 661b1b3d7a2SStefano Zampini sub_schurs->n_subs_par = n_local_parallel_problems; 662b1b3d7a2SStefano Zampini sub_schurs->n_subs_seq_g = n_sequential_problems; 663b1b3d7a2SStefano Zampini sub_schurs->n_subs_par_g = n_parallel_problems; 664b1b3d7a2SStefano Zampini sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par; 665b1b3d7a2SStefano Zampini sub_schurs->is_subs = all_cc; 666b1b3d7a2SStefano Zampini sub_schurs->index_sequential = index_sequential; 667b1b3d7a2SStefano Zampini sub_schurs->index_parallel = index_parallel; 668b1b3d7a2SStefano Zampini sub_schurs->auxglobal_sequential = auxglobal_sequential; 669b1b3d7a2SStefano Zampini sub_schurs->auxglobal_parallel = auxglobal_parallel; 670b1b3d7a2SStefano Zampini 671b1b3d7a2SStefano Zampini /* free workspace */ 672b1b3d7a2SStefano Zampini ierr = PetscFree(auxmapping);CHKERRQ(ierr); 673b1b3d7a2SStefano Zampini ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr); 674b1b3d7a2SStefano Zampini 675b1b3d7a2SStefano Zampini PetscFunctionReturn(0); 676b1b3d7a2SStefano Zampini } 677b1b3d7a2SStefano Zampini 678b1b3d7a2SStefano Zampini #undef __FUNCT__ 67934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate" 68034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 68134a97f8cSStefano Zampini { 68234a97f8cSStefano Zampini PCBDDCSubSchurs schurs_ctx; 68334a97f8cSStefano Zampini PetscErrorCode ierr; 68434a97f8cSStefano Zampini 68534a97f8cSStefano Zampini PetscFunctionBegin; 68634a97f8cSStefano Zampini ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 6875ff63025SStefano Zampini schurs_ctx->n_subs = 0; 68834a97f8cSStefano Zampini *sub_schurs = schurs_ctx; 68934a97f8cSStefano Zampini PetscFunctionReturn(0); 69034a97f8cSStefano Zampini } 69134a97f8cSStefano Zampini 69234a97f8cSStefano Zampini #undef __FUNCT__ 69334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy" 69434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 69534a97f8cSStefano Zampini { 69634a97f8cSStefano Zampini PetscErrorCode ierr; 69734a97f8cSStefano Zampini 69834a97f8cSStefano Zampini PetscFunctionBegin; 69934a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 70034a97f8cSStefano Zampini ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 70134a97f8cSStefano Zampini PetscFunctionReturn(0); 70234a97f8cSStefano Zampini } 70334a97f8cSStefano Zampini 70434a97f8cSStefano Zampini #undef __FUNCT__ 70534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset" 70634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 70734a97f8cSStefano Zampini { 70834a97f8cSStefano Zampini PetscInt i; 70934a97f8cSStefano Zampini PetscErrorCode ierr; 71034a97f8cSStefano Zampini 71134a97f8cSStefano Zampini PetscFunctionBegin; 7121e9c79c2SStefano Zampini ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 713b1b3d7a2SStefano Zampini ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 714b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 715b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 7165db18549SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 7175db18549SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 71841c3ba1bSStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 71941c3ba1bSStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 720*a1337663SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 7215db18549SStefano Zampini ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 72234a97f8cSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 723b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 72434a97f8cSStefano Zampini ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 72534a97f8cSStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 7262b973097SStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr); 7279221af80SStefano Zampini ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr); 7289221af80SStefano Zampini ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr); 72934a97f8cSStefano Zampini } 73068270318SStefano Zampini ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 7315ff63025SStefano Zampini if (sub_schurs->n_subs) { 732b1b3d7a2SStefano Zampini ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 7332b973097SStefano Zampini ierr = PetscFree5(sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->S_Ej_tilda,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr); 734b1b3d7a2SStefano Zampini ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr); 735b1b3d7a2SStefano Zampini ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr); 736b1b3d7a2SStefano Zampini ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr); 7375ff63025SStefano Zampini } 73834a97f8cSStefano Zampini sub_schurs->n_subs = 0; 73934a97f8cSStefano Zampini PetscFunctionReturn(0); 74034a97f8cSStefano Zampini } 74134a97f8cSStefano Zampini 74234a97f8cSStefano Zampini #undef __FUNCT__ 74334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 7442a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 74534a97f8cSStefano Zampini { 74634a97f8cSStefano Zampini PetscInt i,j,n; 74734a97f8cSStefano Zampini PetscErrorCode ierr; 74834a97f8cSStefano Zampini 74934a97f8cSStefano Zampini PetscFunctionBegin; 75034a97f8cSStefano Zampini n = 0; 75134a97f8cSStefano Zampini for (i=-n_prev;i<0;i++) { 75234a97f8cSStefano Zampini PetscInt start_dof = queue_tip[i]; 75334a97f8cSStefano Zampini for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 75434a97f8cSStefano Zampini PetscInt dof = adjncy[j]; 75534a97f8cSStefano Zampini if (!PetscBTLookup(touched,dof)) { 75634a97f8cSStefano Zampini ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 75734a97f8cSStefano Zampini queue_tip[n] = dof; 75834a97f8cSStefano Zampini n++; 75934a97f8cSStefano Zampini } 76034a97f8cSStefano Zampini } 76134a97f8cSStefano Zampini } 76234a97f8cSStefano Zampini *n_added = n; 76334a97f8cSStefano Zampini PetscFunctionReturn(0); 76434a97f8cSStefano Zampini } 765