xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 5ff630250620d1ec1a1edb70acde3521db631784)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
334a97f8cSStefano Zampini 
434a97f8cSStefano Zampini static PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
534a97f8cSStefano Zampini 
634a97f8cSStefano Zampini #undef __FUNCT__
734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
934a97f8cSStefano Zampini {
1034a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
1134a97f8cSStefano Zampini   PetscErrorCode  ierr;
1234a97f8cSStefano Zampini 
1334a97f8cSStefano Zampini   PetscFunctionBegin;
1434a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
15*5ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
1634a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
1734a97f8cSStefano Zampini   PetscFunctionReturn(0);
1834a97f8cSStefano Zampini }
1934a97f8cSStefano Zampini 
2034a97f8cSStefano Zampini #undef __FUNCT__
2134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
2234a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
2334a97f8cSStefano Zampini {
2434a97f8cSStefano Zampini   PetscErrorCode ierr;
2534a97f8cSStefano Zampini 
2634a97f8cSStefano Zampini   PetscFunctionBegin;
2734a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
2834a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
2934a97f8cSStefano Zampini   PetscFunctionReturn(0);
3034a97f8cSStefano Zampini }
3134a97f8cSStefano Zampini 
3234a97f8cSStefano Zampini #undef __FUNCT__
3334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
3434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
3534a97f8cSStefano Zampini {
3634a97f8cSStefano Zampini   PetscInt       i;
3734a97f8cSStefano Zampini   PetscErrorCode ierr;
3834a97f8cSStefano Zampini 
3934a97f8cSStefano Zampini   PetscFunctionBegin;
4034a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4134a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_I[i]);CHKERRQ(ierr);
4234a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
4334a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
4434a97f8cSStefano Zampini     ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
4534a97f8cSStefano Zampini     ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
4634a97f8cSStefano Zampini   }
47*5ff63025SStefano Zampini   if (sub_schurs->n_subs) {
4834a97f8cSStefano Zampini     ierr = PetscFree5(sub_schurs->is_AEj_I,sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
49*5ff63025SStefano Zampini   }
5034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
5134a97f8cSStefano Zampini   PetscFunctionReturn(0);
5234a97f8cSStefano Zampini }
5334a97f8cSStefano Zampini 
5434a97f8cSStefano Zampini #undef __FUNCT__
5534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
562a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
5734a97f8cSStefano Zampini {
5834a97f8cSStefano Zampini   PetscInt       i,j,n;
5934a97f8cSStefano Zampini   PetscErrorCode ierr;
6034a97f8cSStefano Zampini 
6134a97f8cSStefano Zampini   PetscFunctionBegin;
6234a97f8cSStefano Zampini   n = 0;
6334a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
6434a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
6534a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
6634a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
6734a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
6834a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
6934a97f8cSStefano Zampini         queue_tip[n] = dof;
7034a97f8cSStefano Zampini         n++;
7134a97f8cSStefano Zampini       }
7234a97f8cSStefano Zampini     }
7334a97f8cSStefano Zampini   }
7434a97f8cSStefano Zampini   *n_added = n;
7534a97f8cSStefano Zampini   PetscFunctionReturn(0);
7634a97f8cSStefano Zampini }
7734a97f8cSStefano Zampini 
7834a97f8cSStefano Zampini #undef __FUNCT__
7934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
807821e9e7SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat S, IS is_A_I, IS is_A_B, PetscInt ncc, IS is_cc[], PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers)
8134a97f8cSStefano Zampini {
8234a97f8cSStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
832a155e38SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
842a155e38SStefano Zampini   IS                     is_I,*is_subset_B;
852a155e38SStefano Zampini   ISLocalToGlobalMapping BtoNmap;
862a155e38SStefano Zampini   PetscInt               i;
8734a97f8cSStefano Zampini   PetscBool              is_sorted;
8834a97f8cSStefano Zampini   PetscErrorCode         ierr;
8934a97f8cSStefano Zampini 
9034a97f8cSStefano Zampini   PetscFunctionBegin;
9134a97f8cSStefano Zampini   ierr = ISSorted(is_A_I,&is_sorted);CHKERRQ(ierr);
9234a97f8cSStefano Zampini   if (!is_sorted) {
9334a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
9434a97f8cSStefano Zampini   }
9534a97f8cSStefano Zampini   ierr = ISSorted(is_A_B,&is_sorted);CHKERRQ(ierr);
9634a97f8cSStefano Zampini   if (!is_sorted) {
9734a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
9834a97f8cSStefano Zampini   }
9934a97f8cSStefano Zampini 
10034a97f8cSStefano Zampini   /* get Schur complement matrices */
10134a97f8cSStefano Zampini   ierr = MatSchurComplementGetSubMatrices(S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
10234a97f8cSStefano Zampini 
10334a97f8cSStefano Zampini   /* allocate space for schur complements */
10434a97f8cSStefano Zampini   ierr = PetscMalloc5(ncc,&sub_schurs->is_AEj_I,ncc,&sub_schurs->is_AEj_B,ncc,&sub_schurs->S_Ej,ncc,&sub_schurs->work1,ncc,&sub_schurs->work2);CHKERRQ(ierr);
1052a155e38SStefano Zampini   ierr = PetscMalloc4(ncc,&is_subset_B,ncc,&AE_IE,ncc,&AE_EI,ncc,&AE_EE);CHKERRQ(ierr);
10634a97f8cSStefano Zampini   sub_schurs->n_subs = ncc;
10734a97f8cSStefano Zampini 
1082a155e38SStefano Zampini   /* maps */
1092a155e38SStefano Zampini   if (sub_schurs->n_subs && nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
1102a155e38SStefano Zampini     ISLocalToGlobalMapping ItoNmap;
1112a155e38SStefano Zampini     PetscBT                touched;
11234a97f8cSStefano Zampini     const PetscInt*        idx_B;
1132a155e38SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
1142a155e38SStefano Zampini 
1152a155e38SStefano Zampini     /* get sizes */
1162a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_I,&n_I);CHKERRQ(ierr);
1172a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_B,&n_B);CHKERRQ(ierr);
1182a155e38SStefano Zampini 
1192a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is_A_I,&ItoNmap);CHKERRQ(ierr);
1202a155e38SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
1212a155e38SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
1222a155e38SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
12334a97f8cSStefano Zampini 
12434a97f8cSStefano Zampini     /* all boundary dofs must be skipped when adding layers */
12534a97f8cSStefano Zampini     ierr = ISGetIndices(is_A_B,&idx_B);CHKERRQ(ierr);
12634a97f8cSStefano Zampini     for (j=0;j<n_B;j++) {
12734a97f8cSStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
12834a97f8cSStefano Zampini     }
1292a155e38SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
13034a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_A_B,&idx_B);CHKERRQ(ierr);
13134a97f8cSStefano Zampini 
13234a97f8cSStefano Zampini     /* add next layers of dofs */
1332a155e38SStefano Zampini     n_local_dofs = n_B;
1342a155e38SStefano Zampini     n_prev_added = n_B;
13534a97f8cSStefano Zampini     for (layer=0;layer<nlayers;layer++) {
13634a97f8cSStefano Zampini       PetscInt n_added;
1372a155e38SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
1382a155e38SStefano Zampini       if (n_local_dofs > n_I+n_B) {
1392a155e38SStefano 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);
14034a97f8cSStefano Zampini       }
14134a97f8cSStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
14234a97f8cSStefano Zampini       n_prev_added = n_added;
14334a97f8cSStefano Zampini       n_local_dofs += n_added;
1449b0e569cSStefano Zampini       if (!n_added) break;
14534a97f8cSStefano Zampini     }
1462a155e38SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
14734a97f8cSStefano Zampini 
14834a97f8cSStefano Zampini     /* IS for I dofs in original numbering and in I numbering */
1492a155e38SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)ItoNmap),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
1502a155e38SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
1512a155e38SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
1522a155e38SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_I[0],&is_I);CHKERRQ(ierr);
1532a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
15434a97f8cSStefano Zampini 
15534a97f8cSStefano Zampini     /* II block */
15634a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
1572a155e38SStefano Zampini   } else {
1582a155e38SStefano Zampini     PetscInt n_I;
1592a155e38SStefano Zampini 
16034a97f8cSStefano Zampini     /* IS for I dofs in original numbering */
16134a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)is_A_I);CHKERRQ(ierr);
1622a155e38SStefano Zampini     sub_schurs->is_AEj_I[0] = is_A_I;
16334a97f8cSStefano Zampini 
1642a155e38SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
1652a155e38SStefano Zampini     ierr = ISGetSize(is_A_I,&n_I);CHKERRQ(ierr);
16634a97f8cSStefano Zampini     ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_I),n_I,0,1,&is_I);CHKERRQ(ierr);
16734a97f8cSStefano Zampini 
16834a97f8cSStefano Zampini     /* II block is the same */
16934a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
17034a97f8cSStefano Zampini     AE_II = A_II;
17134a97f8cSStefano Zampini   }
17234a97f8cSStefano Zampini 
1732a155e38SStefano Zampini   /* TODO: just for compatibility with the previous version, needs to be fixed */
1742a155e38SStefano Zampini   for (i=1;i<sub_schurs->n_subs;i++) {
1752a155e38SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
1762a155e38SStefano Zampini     sub_schurs->is_AEj_I[i] = sub_schurs->is_AEj_I[0];
1772a155e38SStefano Zampini   }
17834a97f8cSStefano Zampini 
1792a155e38SStefano Zampini   /* subsets in original and boundary numbering */
1802a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is_A_B,&BtoNmap);CHKERRQ(ierr);
1812a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1822a155e38SStefano Zampini     ierr = ISDuplicate(is_cc[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
1832a155e38SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
1842a155e38SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
1852a155e38SStefano Zampini   }
1862a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&BtoNmap);CHKERRQ(ierr);
1872a155e38SStefano Zampini 
1882a155e38SStefano Zampini   /* EE block */
1892a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1902a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
1912a155e38SStefano Zampini   }
1922a155e38SStefano Zampini   /* IE block */
1932a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1942a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
1952a155e38SStefano Zampini   }
19634a97f8cSStefano Zampini   /* EI block */
1972a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1982a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
1992a155e38SStefano Zampini   }
20034a97f8cSStefano Zampini 
20134a97f8cSStefano Zampini   /* setup Schur complements on subset */
2022a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
2032a155e38SStefano Zampini     ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
2048a26ef87SStefano Zampini     ierr = MatCreateVecs(sub_schurs->S_Ej[i],&sub_schurs->work1[i],&sub_schurs->work2[i]);CHKERRQ(ierr);
20534a97f8cSStefano Zampini     if (AE_II == A_II) { /* we can reuse the same ksp */
20634a97f8cSStefano Zampini       KSP ksp;
20734a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&ksp);CHKERRQ(ierr);
20834a97f8cSStefano Zampini       ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
20934a97f8cSStefano Zampini     } else { /* build new ksp object which inherits ksp and pc types from the original one */
21034a97f8cSStefano Zampini       KSP      origksp,schurksp;
21134a97f8cSStefano Zampini       PC       origpc,schurpc;
21234a97f8cSStefano Zampini       KSPType  ksp_type;
21334a97f8cSStefano Zampini       PCType   pc_type;
21434a97f8cSStefano Zampini       PetscInt n_internal;
21534a97f8cSStefano Zampini 
21634a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&origksp);CHKERRQ(ierr);
21734a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
21834a97f8cSStefano Zampini       ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
21934a97f8cSStefano Zampini       ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
22034a97f8cSStefano Zampini       ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
22134a97f8cSStefano Zampini       ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
22234a97f8cSStefano Zampini       ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
22334a97f8cSStefano Zampini       ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
22434a97f8cSStefano Zampini       ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
22534a97f8cSStefano Zampini       if (n_internal) { /* UMFPACK gives error with 0 sized problems */
22634a97f8cSStefano Zampini         MatSolverPackage solver=NULL;
22734a97f8cSStefano Zampini         ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
22834a97f8cSStefano Zampini         if (solver) {
22934a97f8cSStefano Zampini           ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
23034a97f8cSStefano Zampini         }
23134a97f8cSStefano Zampini       }
23234a97f8cSStefano Zampini       ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
23334a97f8cSStefano Zampini     }
23434a97f8cSStefano Zampini   }
23534a97f8cSStefano Zampini   /* free */
2362a155e38SStefano Zampini   ierr = ISDestroy(&is_I);CHKERRQ(ierr);
2372a155e38SStefano Zampini   ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
2382a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
2392a155e38SStefano Zampini     ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
2402a155e38SStefano Zampini     ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
2412a155e38SStefano Zampini     ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
2422a155e38SStefano Zampini     ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
2432a155e38SStefano Zampini   }
2442a155e38SStefano Zampini   ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
24534a97f8cSStefano Zampini   PetscFunctionReturn(0);
24634a97f8cSStefano Zampini }
247