134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h> 234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 308122e43SStefano Zampini #include <petscblaslapack.h> 434a97f8cSStefano Zampini 53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 65ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*); 73202ece2SStefano Zampini 83202ece2SStefano Zampini #undef __FUNCT__ 93202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur" 105ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 113202ece2SStefano Zampini { 123202ece2SStefano Zampini Mat B, C, D, Bd, Cd, AinvBd; 133202ece2SStefano Zampini KSP ksp; 143202ece2SStefano Zampini PC pc; 153202ece2SStefano Zampini PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 163202ece2SStefano Zampini PetscReal fill = 2.0; 17f11841e3SStefano Zampini PetscInt n_I; 183202ece2SStefano Zampini PetscMPIInt size; 193202ece2SStefano Zampini PetscErrorCode ierr; 203202ece2SStefano Zampini 213202ece2SStefano Zampini PetscFunctionBegin; 223202ece2SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 233202ece2SStefano Zampini if (size != 1) { 243202ece2SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 253202ece2SStefano Zampini } 26f11841e3SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 27f11841e3SStefano Zampini PetscBool Sdense; 28f11841e3SStefano Zampini 29f11841e3SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr); 30f11841e3SStefano Zampini if (!Sdense) { 31f11841e3SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense"); 32f11841e3SStefano Zampini } 33f11841e3SStefano Zampini } 343202ece2SStefano Zampini ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 353202ece2SStefano Zampini ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 363202ece2SStefano Zampini ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 373202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 383202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 393202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 403202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 413202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 42f11841e3SStefano Zampini ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr); 43f11841e3SStefano Zampini if (n_I) { 443202ece2SStefano Zampini if (!Bdense) { 453202ece2SStefano Zampini ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 463202ece2SStefano Zampini } else { 473202ece2SStefano Zampini Bd = B; 483202ece2SStefano Zampini } 493202ece2SStefano Zampini 503202ece2SStefano Zampini if (isLU || isILU || isCHOL) { 513202ece2SStefano Zampini Mat fact; 523202ece2SStefano Zampini ierr = KSPSetUp(ksp);CHKERRQ(ierr); 533202ece2SStefano Zampini ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 543202ece2SStefano Zampini ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 553202ece2SStefano Zampini ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 563202ece2SStefano Zampini } else { 5707b1e237SStefano Zampini PetscBool ex = PETSC_TRUE; 5807b1e237SStefano Zampini 5907b1e237SStefano Zampini if (ex) { 603202ece2SStefano Zampini Mat Ainvd; 613202ece2SStefano Zampini 623202ece2SStefano Zampini ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 633202ece2SStefano Zampini ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 643202ece2SStefano Zampini ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 6507b1e237SStefano Zampini } else { 6607b1e237SStefano Zampini Vec sol,rhs; 6707b1e237SStefano Zampini PetscScalar *arrayrhs,*arraysol; 6807b1e237SStefano Zampini PetscInt i,nrhs,n; 6907b1e237SStefano Zampini 7007b1e237SStefano Zampini ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 7107b1e237SStefano Zampini ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr); 7207b1e237SStefano Zampini ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr); 7307b1e237SStefano Zampini ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr); 7407b1e237SStefano Zampini ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr); 7507b1e237SStefano Zampini ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 7607b1e237SStefano Zampini for (i=0;i<nrhs;i++) { 7707b1e237SStefano Zampini ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr); 7807b1e237SStefano Zampini ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr); 7907b1e237SStefano Zampini ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr); 8007b1e237SStefano Zampini ierr = VecResetArray(rhs);CHKERRQ(ierr); 8107b1e237SStefano Zampini ierr = VecResetArray(sol);CHKERRQ(ierr); 8207b1e237SStefano Zampini } 8307b1e237SStefano Zampini ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr); 8407b1e237SStefano Zampini ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr); 8507b1e237SStefano Zampini } 863202ece2SStefano Zampini } 875ec10c6aSStefano Zampini if (!Bdense & !issym) { 883202ece2SStefano Zampini ierr = MatDestroy(&Bd);CHKERRQ(ierr); 893202ece2SStefano Zampini } 905ec10c6aSStefano Zampini 915ec10c6aSStefano Zampini if (!issym) { 923202ece2SStefano Zampini if (!Cdense) { 933202ece2SStefano Zampini ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 943202ece2SStefano Zampini } else { 953202ece2SStefano Zampini Cd = C; 963202ece2SStefano Zampini } 975ec10c6aSStefano Zampini ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 983202ece2SStefano Zampini if (!Cdense) { 993202ece2SStefano Zampini ierr = MatDestroy(&Cd);CHKERRQ(ierr); 1003202ece2SStefano Zampini } 1015ec10c6aSStefano Zampini } else { 1025ec10c6aSStefano Zampini ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr); 1035ec10c6aSStefano Zampini if (!Bdense) { 1045ec10c6aSStefano Zampini ierr = MatDestroy(&Bd);CHKERRQ(ierr); 1055ec10c6aSStefano Zampini } 1065ec10c6aSStefano Zampini } 1075ec10c6aSStefano Zampini ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 108f11841e3SStefano Zampini } 1093202ece2SStefano Zampini 1103202ece2SStefano Zampini if (D) { 1113202ece2SStefano Zampini Mat Dd; 1123202ece2SStefano Zampini PetscBool Ddense; 1133202ece2SStefano Zampini 1143202ece2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 1153202ece2SStefano Zampini if (!Ddense) { 1163202ece2SStefano Zampini ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 1173202ece2SStefano Zampini } else { 1183202ece2SStefano Zampini Dd = D; 1193202ece2SStefano Zampini } 120f11841e3SStefano Zampini if (n_I) { 1213202ece2SStefano Zampini ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 122f11841e3SStefano Zampini } else { 123f11841e3SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 124f11841e3SStefano Zampini ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr); 125f11841e3SStefano Zampini } else { 126f11841e3SStefano Zampini ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 127f11841e3SStefano Zampini } 128f11841e3SStefano Zampini } 1293202ece2SStefano Zampini if (!Ddense) { 1303202ece2SStefano Zampini ierr = MatDestroy(&Dd);CHKERRQ(ierr); 1313202ece2SStefano Zampini } 1323202ece2SStefano Zampini } else { 1333202ece2SStefano Zampini ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 1343202ece2SStefano Zampini } 1353202ece2SStefano Zampini PetscFunctionReturn(0); 1363202ece2SStefano Zampini } 13734a97f8cSStefano Zampini 13834a97f8cSStefano Zampini #undef __FUNCT__ 1391580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp" 1405a95e1ceSStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool invert_Stildas, PetscBool use_edges, PetscBool use_faces) 141b1b3d7a2SStefano Zampini { 1425a95e1ceSStefano Zampini Mat A_II,A_IB,A_BI,A_BB,AE_II; 143d2627357SStefano Zampini Mat S_all,S_all_inv; 144d2627357SStefano Zampini Mat global_schur_subsets,work_mat; 14508122e43SStefano Zampini Mat S_Ej_tilda_all,S_Ej_inv_all; 1465db18549SStefano Zampini ISLocalToGlobalMapping l2gmap_subsets; 1475a95e1ceSStefano Zampini IS is_I,temp_is; 148d648f858SStefano Zampini PetscInt *nnz,*all_local_idx_G,*all_local_idx_N; 1495a95e1ceSStefano Zampini PetscInt i,subset_size,max_subset_size; 150883469d8SStefano Zampini PetscInt extra,local_size,global_size; 15108122e43SStefano Zampini PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 15206a4b1faSStefano Zampini PetscScalar *Bwork; 1535a95e1ceSStefano Zampini PetscSubcomm subcomm; 1545a95e1ceSStefano Zampini PetscMPIInt color,rank; 1555a95e1ceSStefano Zampini MPI_Comm comm_n; 156b1b3d7a2SStefano Zampini PetscErrorCode ierr; 157b1b3d7a2SStefano Zampini 158b1b3d7a2SStefano Zampini PetscFunctionBegin; 1595a95e1ceSStefano Zampini /* preliminary checks */ 1605a95e1ceSStefano Zampini if (!sub_schurs->use_mumps && compute_Stilda) { 1615a95e1ceSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS"); 1625a95e1ceSStefano Zampini } 1635a95e1ceSStefano Zampini /* determine if we are dealing with hermitian positive definite problems */ 1645a95e1ceSStefano Zampini sub_schurs->is_hermitian = PETSC_FALSE; 1655a95e1ceSStefano Zampini sub_schurs->is_posdef = PETSC_FALSE; 1665a95e1ceSStefano Zampini if (sub_schurs->A) { 1675a95e1ceSStefano Zampini PetscInt lsize; 1685a95e1ceSStefano Zampini 1695a95e1ceSStefano Zampini ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr); 1705a95e1ceSStefano Zampini if (lsize) { 1715a95e1ceSStefano Zampini ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr); 1725a95e1ceSStefano Zampini if (sub_schurs->is_hermitian) { 1735a95e1ceSStefano Zampini PetscScalar val; 1745a95e1ceSStefano Zampini Vec vec1,vec2; 1755a95e1ceSStefano Zampini 1765a95e1ceSStefano Zampini ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr); 1775a95e1ceSStefano Zampini ierr = VecSetRandom(vec1,NULL); 1785a95e1ceSStefano Zampini ierr = VecCopy(vec1,vec2);CHKERRQ(ierr); 1795a95e1ceSStefano Zampini ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr); 1805a95e1ceSStefano Zampini ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr); 1815a95e1ceSStefano Zampini if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE; 1825a95e1ceSStefano Zampini ierr = VecDestroy(&vec1);CHKERRQ(ierr); 1835a95e1ceSStefano Zampini ierr = VecDestroy(&vec2);CHKERRQ(ierr); 1845a95e1ceSStefano Zampini } 1855a95e1ceSStefano Zampini } else { 1865a95e1ceSStefano Zampini sub_schurs->is_hermitian = PETSC_TRUE; 1875a95e1ceSStefano Zampini sub_schurs->is_posdef = PETSC_TRUE; 1885a95e1ceSStefano Zampini } 1895a95e1ceSStefano Zampini if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) { 1905a95e1ceSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported"); 1915a95e1ceSStefano Zampini } 1925a95e1ceSStefano Zampini } 1935a95e1ceSStefano Zampini /* restrict work on active processes */ 1945a95e1ceSStefano Zampini color = 0; 1955a95e1ceSStefano Zampini if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */ 1965a95e1ceSStefano Zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr); 1975a95e1ceSStefano Zampini ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr); 1985a95e1ceSStefano Zampini ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr); 1995a95e1ceSStefano Zampini ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr); 2005a95e1ceSStefano Zampini ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr); 2015a95e1ceSStefano Zampini ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr); 2025a95e1ceSStefano Zampini if (!sub_schurs->n_subs) { 2035a95e1ceSStefano Zampini ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 2045a95e1ceSStefano Zampini PetscFunctionReturn(0); 2055a95e1ceSStefano Zampini } 2065a95e1ceSStefano Zampini 207b1b3d7a2SStefano Zampini /* get Schur complement matrices */ 208883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 209f11841e3SStefano Zampini PetscBool isseqaij; 210b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr); 211f11841e3SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 212f11841e3SStefano Zampini if (!isseqaij) { 213f11841e3SStefano Zampini ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr); 214f11841e3SStefano Zampini ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr); 215f11841e3SStefano Zampini ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr); 216f11841e3SStefano Zampini } 217a58a30b4SStefano Zampini } else { 2185a95e1ceSStefano Zampini A_II = NULL; 2195a95e1ceSStefano Zampini A_IB = NULL; 2205a95e1ceSStefano Zampini A_BI = NULL; 2215a95e1ceSStefano Zampini A_BB = NULL; 222b1b3d7a2SStefano Zampini } 2235a95e1ceSStefano Zampini S_all = NULL; 2245a95e1ceSStefano Zampini S_all_inv = NULL; 2255a95e1ceSStefano Zampini S_Ej_tilda_all = NULL; 2265a95e1ceSStefano Zampini S_Ej_inv_all = NULL; 227b1b3d7a2SStefano Zampini 228b1b3d7a2SStefano Zampini /* determine interior problems */ 229b96c3477SStefano Zampini ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 2303dc780c3SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr); 2313dc780c3SStefano Zampini if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 232b1b3d7a2SStefano Zampini PetscBT touched; 233b1b3d7a2SStefano Zampini const PetscInt* idx_B; 234b1b3d7a2SStefano Zampini PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 235b1b3d7a2SStefano Zampini 2363dc780c3SStefano Zampini if (xadj == NULL || adjncy == NULL) { 2373dc780c3SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency"); 2383dc780c3SStefano Zampini } 239b1b3d7a2SStefano Zampini /* get sizes */ 240b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 241b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 242b1b3d7a2SStefano Zampini 243b1b3d7a2SStefano Zampini ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 244b1b3d7a2SStefano Zampini ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 245b1b3d7a2SStefano Zampini ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 246b1b3d7a2SStefano Zampini 247b1b3d7a2SStefano Zampini /* all boundary dofs must be skipped when adding layers */ 248b1b3d7a2SStefano Zampini ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 249b1b3d7a2SStefano Zampini for (j=0;j<n_B;j++) { 250b1b3d7a2SStefano Zampini ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 251b1b3d7a2SStefano Zampini } 252b1b3d7a2SStefano Zampini ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 253b1b3d7a2SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 254b1b3d7a2SStefano Zampini 255b1b3d7a2SStefano Zampini /* add prescribed number of layers of dofs */ 256b1b3d7a2SStefano Zampini n_local_dofs = n_B; 257b1b3d7a2SStefano Zampini n_prev_added = n_B; 258b1b3d7a2SStefano Zampini for (layer=0;layer<nlayers;layer++) { 259b1b3d7a2SStefano Zampini PetscInt n_added; 260b1b3d7a2SStefano Zampini if (n_local_dofs == n_I+n_B) break; 261b1b3d7a2SStefano Zampini if (n_local_dofs > n_I+n_B) { 262b1b3d7a2SStefano 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); 263b1b3d7a2SStefano Zampini } 264b1b3d7a2SStefano Zampini ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 265b1b3d7a2SStefano Zampini n_prev_added = n_added; 266b1b3d7a2SStefano Zampini n_local_dofs += n_added; 267b1b3d7a2SStefano Zampini if (!n_added) break; 268b1b3d7a2SStefano Zampini } 269b1b3d7a2SStefano Zampini ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 270b1b3d7a2SStefano Zampini 271883469d8SStefano Zampini /* IS for I layer dofs in original numbering */ 27268270318SStefano 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); 273b1b3d7a2SStefano Zampini ierr = PetscFree(local_numbering);CHKERRQ(ierr); 27468270318SStefano Zampini ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr); 275883469d8SStefano Zampini /* IS for I layer dofs in I numbering */ 276883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 277b1b3d7a2SStefano Zampini ISLocalToGlobalMapping ItoNmap; 278b1b3d7a2SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 27968270318SStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr); 280b1b3d7a2SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 281b1b3d7a2SStefano Zampini 282b1b3d7a2SStefano Zampini /* II block */ 283b1b3d7a2SStefano Zampini ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 284b1b3d7a2SStefano Zampini } 285b1b3d7a2SStefano Zampini } else { 286b1b3d7a2SStefano Zampini PetscInt n_I; 287b1b3d7a2SStefano Zampini 288b1b3d7a2SStefano Zampini /* IS for I dofs in original numbering */ 289b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 29068270318SStefano Zampini sub_schurs->is_I_layer = sub_schurs->is_I; 291b1b3d7a2SStefano Zampini 292b1b3d7a2SStefano Zampini /* IS for I dofs in I numbering (strided 1) */ 293883469d8SStefano Zampini if (!sub_schurs->use_mumps) { 294b1b3d7a2SStefano Zampini ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 295b1b3d7a2SStefano Zampini ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 296b1b3d7a2SStefano Zampini 297b1b3d7a2SStefano Zampini /* II block is the same */ 298b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 299b1b3d7a2SStefano Zampini AE_II = A_II; 300b1b3d7a2SStefano Zampini } 301b1b3d7a2SStefano Zampini } 3025a95e1ceSStefano Zampini 303883469d8SStefano Zampini /* Get info on subset sizes and sum of all subsets sizes */ 3045a95e1ceSStefano Zampini max_subset_size = 0; 305883469d8SStefano Zampini local_size = 0; 3065a95e1ceSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 3075a95e1ceSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 3085a95e1ceSStefano Zampini max_subset_size = PetscMax(subset_size,max_subset_size); 309883469d8SStefano Zampini local_size += subset_size; 310883469d8SStefano Zampini } 311883469d8SStefano Zampini 312883469d8SStefano Zampini /* Work arrays for local indices */ 313883469d8SStefano Zampini extra = 0; 314883469d8SStefano Zampini if (sub_schurs->use_mumps) { 315883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr); 316883469d8SStefano Zampini } 317883469d8SStefano Zampini ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 318883469d8SStefano Zampini if (extra) { 319883469d8SStefano Zampini const PetscInt *idxs; 320883469d8SStefano Zampini ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 321883469d8SStefano Zampini ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 322883469d8SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 323883469d8SStefano Zampini } 324883469d8SStefano Zampini ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 325883469d8SStefano Zampini 326883469d8SStefano Zampini /* Get local indices in local numbering */ 327883469d8SStefano Zampini local_size = 0; 3285a95e1ceSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 329883469d8SStefano Zampini PetscInt j; 330883469d8SStefano Zampini const PetscInt *idxs; 331883469d8SStefano Zampini 3325a95e1ceSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 3335a95e1ceSStefano Zampini ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 334883469d8SStefano Zampini /* subset indices in local numbering */ 335883469d8SStefano Zampini ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 3365a95e1ceSStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 337883469d8SStefano Zampini for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 338883469d8SStefano Zampini local_size += subset_size; 339883469d8SStefano Zampini } 340883469d8SStefano Zampini 3415a95e1ceSStefano Zampini /* allocate extra workspace needed only for GETRI */ 342d2627357SStefano Zampini Bwork = NULL; 343d2627357SStefano Zampini pivots = NULL; 3445a95e1ceSStefano Zampini if (sub_schurs->n_subs && !sub_schurs->is_hermitian) { 345d2627357SStefano Zampini PetscScalar lwork; 346d2627357SStefano Zampini 347d2627357SStefano Zampini B_lwork = -1; 348d2627357SStefano Zampini ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 349d2627357SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 350d2627357SStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr)); 351d2627357SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 352d2627357SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 353d2627357SStefano Zampini ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 354d2627357SStefano Zampini ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr); 355d2627357SStefano Zampini } 356d2627357SStefano Zampini 357d2627357SStefano Zampini /* prepare parallel matrices for summing up properly schurs on subsets */ 3585a95e1ceSStefano Zampini ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,local_size,all_local_idx_N+extra,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr); 3595a95e1ceSStefano Zampini ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 3605a95e1ceSStefano Zampini ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 361d2627357SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 362d2627357SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 363d2627357SStefano Zampini ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 364d2627357SStefano Zampini ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 3652972d61bSStefano Zampini 3665a95e1ceSStefano Zampini /* subset indices in local boundary numbering */ 3675a95e1ceSStefano Zampini if (!sub_schurs->is_Ej_all) { 3685a95e1ceSStefano Zampini PetscInt *all_local_idx_B; 3695a95e1ceSStefano Zampini 3705a95e1ceSStefano Zampini ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 3715a95e1ceSStefano Zampini ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 3725a95e1ceSStefano Zampini if (subset_size != local_size) { 3735a95e1ceSStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 3745a95e1ceSStefano Zampini } 3755a95e1ceSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 376b1b3d7a2SStefano Zampini } 377b1b3d7a2SStefano Zampini 3785a95e1ceSStefano Zampini /* Local matrix of all local Schur on subsets (transposed) */ 3795a95e1ceSStefano Zampini if (!sub_schurs->S_Ej_all) { 3805a95e1ceSStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 3815a95e1ceSStefano Zampini ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 3825a95e1ceSStefano Zampini ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 3835a95e1ceSStefano Zampini ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 384aa83b6aeSStefano Zampini } else { 3855a95e1ceSStefano Zampini ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr); 386aa83b6aeSStefano Zampini } 387b1b3d7a2SStefano Zampini 3885a95e1ceSStefano Zampini /* Compute Schur complements explicitly */ 3895a95e1ceSStefano Zampini ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 3905a95e1ceSStefano Zampini if (!sub_schurs->use_mumps) { 3915a95e1ceSStefano Zampini Mat S_Ej_expl; 3925a95e1ceSStefano Zampini PetscScalar *work; 3935a95e1ceSStefano Zampini PetscInt j,*dummy_idx; 3945a95e1ceSStefano Zampini PetscBool Sdense; 3955a95e1ceSStefano Zampini 3965a95e1ceSStefano Zampini ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 3975a95e1ceSStefano Zampini local_size = 0; 398b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 3995a95e1ceSStefano Zampini IS is_subset_B; 4005a95e1ceSStefano Zampini Mat AE_EE,AE_IE,AE_EI,S_Ej; 4015a95e1ceSStefano Zampini 4025a95e1ceSStefano Zampini /* subsets in original and boundary numbering */ 4035a95e1ceSStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr); 4045a95e1ceSStefano Zampini /* EE block */ 4055a95e1ceSStefano Zampini ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr); 4065a95e1ceSStefano Zampini /* IE block */ 4075a95e1ceSStefano Zampini ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr); 4085a95e1ceSStefano Zampini /* EI block */ 4095a95e1ceSStefano Zampini if (sub_schurs->is_hermitian) { 4105a95e1ceSStefano Zampini ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr); 4115a95e1ceSStefano Zampini } else { 4125a95e1ceSStefano Zampini ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr); 4135a95e1ceSStefano Zampini } 4145a95e1ceSStefano Zampini ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr); 4155a95e1ceSStefano Zampini ierr = MatDestroy(&AE_EE);CHKERRQ(ierr); 4165a95e1ceSStefano Zampini ierr = MatDestroy(&AE_IE);CHKERRQ(ierr); 4175a95e1ceSStefano Zampini ierr = MatDestroy(&AE_EI);CHKERRQ(ierr); 418b1b3d7a2SStefano Zampini if (AE_II == A_II) { /* we can reuse the same ksp */ 419b1b3d7a2SStefano Zampini KSP ksp; 420b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 4215a95e1ceSStefano Zampini ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr); 422b1b3d7a2SStefano Zampini } else { /* build new ksp object which inherits ksp and pc types from the original one */ 423b1b3d7a2SStefano Zampini KSP origksp,schurksp; 424b1b3d7a2SStefano Zampini PC origpc,schurpc; 425b1b3d7a2SStefano Zampini KSPType ksp_type; 426b1b3d7a2SStefano Zampini PetscInt n_internal; 4275a95e1ceSStefano Zampini PetscBool ispcnone; 428b1b3d7a2SStefano Zampini 429b1b3d7a2SStefano Zampini ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 4305a95e1ceSStefano Zampini ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr); 431b1b3d7a2SStefano Zampini ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 432b1b3d7a2SStefano Zampini ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 433b1b3d7a2SStefano Zampini ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 434b1b3d7a2SStefano Zampini ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 4355a95e1ceSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr); 4365a95e1ceSStefano Zampini if (!ispcnone) { 4375a95e1ceSStefano Zampini PCType pc_type; 438b1b3d7a2SStefano Zampini ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 439b1b3d7a2SStefano Zampini ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 4405a95e1ceSStefano Zampini } else { 4415a95e1ceSStefano Zampini ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr); 4425a95e1ceSStefano Zampini } 443b1b3d7a2SStefano Zampini ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 444b1b3d7a2SStefano Zampini if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 445b1b3d7a2SStefano Zampini MatSolverPackage solver=NULL; 446b1b3d7a2SStefano Zampini ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 447b1b3d7a2SStefano Zampini if (solver) { 448b1b3d7a2SStefano Zampini ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 449b1b3d7a2SStefano Zampini } 450b1b3d7a2SStefano Zampini } 451b1b3d7a2SStefano Zampini ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 452b1b3d7a2SStefano Zampini } 4535a95e1ceSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 4545a95e1ceSStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr); 4555a95e1ceSStefano Zampini ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr); 4565a95e1ceSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 4575a95e1ceSStefano Zampini if (Sdense) { 4585a95e1ceSStefano Zampini for (j=0;j<subset_size;j++) { 4595a95e1ceSStefano Zampini dummy_idx[j]=local_size+j; 460b1b3d7a2SStefano Zampini } 4615a95e1ceSStefano Zampini ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 4625a95e1ceSStefano Zampini } else { 4635a95e1ceSStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 4645a95e1ceSStefano Zampini } 4655a95e1ceSStefano Zampini ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 4665a95e1ceSStefano Zampini local_size += subset_size; 4675a95e1ceSStefano Zampini } 4685a95e1ceSStefano Zampini ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 469b1b3d7a2SStefano Zampini /* free */ 470b1b3d7a2SStefano Zampini ierr = ISDestroy(&is_I);CHKERRQ(ierr); 471b1b3d7a2SStefano Zampini ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 4725a95e1ceSStefano Zampini ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 473883469d8SStefano Zampini } else { 474883469d8SStefano Zampini Mat A,F; 475883469d8SStefano Zampini IS is_A_all; 4765a95e1ceSStefano Zampini PetscScalar *work; 4775a95e1ceSStefano Zampini PetscInt *idxs_schur,n_I,*dummy_idx; 478883469d8SStefano Zampini 479883469d8SStefano Zampini /* get working mat */ 480883469d8SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr); 481883469d8SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 482d648f858SStefano Zampini ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr); 483883469d8SStefano Zampini ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 484883469d8SStefano Zampini 48508122e43SStefano Zampini if (n_I) { 4869ab7bb16SStefano Zampini if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 487883469d8SStefano Zampini ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 488883469d8SStefano Zampini } else { 489883469d8SStefano Zampini ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 490883469d8SStefano Zampini } 491883469d8SStefano Zampini 492883469d8SStefano Zampini /* subsets ordered last */ 493883469d8SStefano Zampini ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr); 494883469d8SStefano Zampini for (i=0;i<local_size;i++) { 495883469d8SStefano Zampini idxs_schur[i] = n_I+i+1; 496883469d8SStefano Zampini } 4975a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS) 498883469d8SStefano Zampini ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr); 4995a95e1ceSStefano Zampini #endif 500883469d8SStefano Zampini ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 501883469d8SStefano Zampini 502883469d8SStefano Zampini /* factorization step */ 5039ab7bb16SStefano Zampini if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 504883469d8SStefano Zampini ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 505883469d8SStefano Zampini ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 506883469d8SStefano Zampini } else { 507883469d8SStefano Zampini ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 508883469d8SStefano Zampini ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 509883469d8SStefano Zampini } 510883469d8SStefano Zampini 511883469d8SStefano Zampini /* get explicit Schur Complement computed during numeric factorization */ 5125a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS) 513883469d8SStefano Zampini ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 5145a95e1ceSStefano Zampini #endif 515883469d8SStefano Zampini ierr = MatDestroy(&F);CHKERRQ(ierr); 51608122e43SStefano Zampini } else { 51708122e43SStefano Zampini ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 51808122e43SStefano Zampini } 51908122e43SStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 520d2627357SStefano Zampini ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 5215db18549SStefano Zampini 5225a95e1ceSStefano Zampini if (compute_Stilda) { /* TODO PICKUP BETTER NAMES */ 523a1337663SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr); 524a1337663SStefano Zampini ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 525a1337663SStefano Zampini ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 526a1337663SStefano Zampini ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 52708122e43SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr); 52808122e43SStefano Zampini ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 52908122e43SStefano Zampini ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr); 53008122e43SStefano Zampini ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr); 531a1337663SStefano Zampini } 53206a4b1faSStefano Zampini 5335a95e1ceSStefano Zampini /* Get St^-1 */ 5345a95e1ceSStefano Zampini if (compute_Stilda) { 535d2627357SStefano Zampini PetscScalar *vals; 5365a95e1ceSStefano Zampini ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr); 537d2627357SStefano Zampini ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr); 538d2627357SStefano Zampini ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr); 539d2627357SStefano Zampini if (!sub_schurs->is_hermitian) { 540d2627357SStefano Zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr)); 541d2627357SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 542d2627357SStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 543d2627357SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 544d2627357SStefano Zampini } else { 545d2627357SStefano Zampini PetscInt j,k; 546d2627357SStefano Zampini 547d2627357SStefano Zampini PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr)); 548d2627357SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 549d2627357SStefano Zampini PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr)); 550d2627357SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 551d2627357SStefano Zampini for (j=0;j<B_N;j++) { 552d2627357SStefano Zampini for (k=j+1;k<B_N;k++) { 553d2627357SStefano Zampini vals[k*B_N+j] = vals[j*B_N+k]; 554d2627357SStefano Zampini } 555d2627357SStefano Zampini } 556d2627357SStefano Zampini } 557d2627357SStefano Zampini ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr); 558d2627357SStefano Zampini } 559d2627357SStefano Zampini 560*9087bf02SStefano Zampini /* Work arrays */ 561*9087bf02SStefano Zampini if (sub_schurs->n_subs == 1) { 562*9087bf02SStefano Zampini ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 563*9087bf02SStefano Zampini } else { 564*9087bf02SStefano Zampini ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr); 565*9087bf02SStefano Zampini } 566*9087bf02SStefano Zampini 5675a95e1ceSStefano Zampini local_size = 0; 56865d8bf0aSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 5695a95e1ceSStefano Zampini Mat S_Ej; 57065d8bf0aSStefano Zampini IS is_E; 57165d8bf0aSStefano Zampini PetscInt j; 57265d8bf0aSStefano Zampini 5735a95e1ceSStefano Zampini /* get S_E */ 574b96c3477SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 575*9087bf02SStefano Zampini if (sub_schurs->n_subs == 1) { 576*9087bf02SStefano Zampini ierr = MatDenseGetArray(S_all,&work);CHKERRQ(ierr); 577*9087bf02SStefano Zampini S_Ej = NULL; 578*9087bf02SStefano Zampini is_E = NULL; 579*9087bf02SStefano Zampini } else { 5805a95e1ceSStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,local_size,1,&is_E);CHKERRQ(ierr); 5815a95e1ceSStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej);CHKERRQ(ierr); 5825a95e1ceSStefano Zampini ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr); 583*9087bf02SStefano Zampini } 5845a95e1ceSStefano Zampini /* insert S_E values */ 585a1337663SStefano Zampini for (j=0;j<subset_size;j++) { 586a1337663SStefano Zampini dummy_idx[j]=local_size+j; 587a1337663SStefano Zampini } 5885a95e1ceSStefano Zampini ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 589a1337663SStefano Zampini 5905a95e1ceSStefano Zampini /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */ 591d2627357SStefano Zampini if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) { 5925a95e1ceSStefano Zampini /* get S_E^-1 */ 59308122e43SStefano Zampini ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 59408122e43SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 5952972d61bSStefano Zampini if (!sub_schurs->is_hermitian) { 5965a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr)); 59708122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 5985a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 59908122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 6002972d61bSStefano Zampini } else { 6012972d61bSStefano Zampini PetscInt j,k; 6022972d61bSStefano Zampini 6035a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr)); 6042972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 6055a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr)); 6062972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 6072972d61bSStefano Zampini for (j=0;j<B_N;j++) { 6082972d61bSStefano Zampini for (k=j+1;k<B_N;k++) { 6095a95e1ceSStefano Zampini work[k*B_N+j] = work[j*B_N+k]; 6102972d61bSStefano Zampini } 6112972d61bSStefano Zampini } 6122972d61bSStefano Zampini } 61308122e43SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 6145a95e1ceSStefano Zampini ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 6155a95e1ceSStefano Zampini 6165a95e1ceSStefano Zampini /* get St_E^-1 */ 617*9087bf02SStefano Zampini if (sub_schurs->n_subs == 1) { 618*9087bf02SStefano Zampini ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr); 619*9087bf02SStefano Zampini ierr = MatDenseGetArray(S_all_inv,&work);CHKERRQ(ierr); 620*9087bf02SStefano Zampini } else { 6215a95e1ceSStefano Zampini ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr); 622*9087bf02SStefano Zampini } 6235a95e1ceSStefano Zampini ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr); 624*9087bf02SStefano Zampini if (sub_schurs->n_subs == 1) { 625*9087bf02SStefano Zampini ierr = MatDenseRestoreArray(S_all_inv,&work);CHKERRQ(ierr); 626*9087bf02SStefano Zampini } 627*9087bf02SStefano Zampini ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr); 628*9087bf02SStefano Zampini } else if (sub_schurs->n_subs == 1) { 629*9087bf02SStefano Zampini ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr); 63008122e43SStefano Zampini } 6315a95e1ceSStefano Zampini ierr = MatDestroy(&S_Ej);CHKERRQ(ierr); 63265d8bf0aSStefano Zampini ierr = ISDestroy(&is_E);CHKERRQ(ierr); 633883469d8SStefano Zampini local_size += subset_size; 634883469d8SStefano Zampini } 635*9087bf02SStefano Zampini if (sub_schurs->n_subs == 1) { 636*9087bf02SStefano Zampini ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 637*9087bf02SStefano Zampini } else { 6385ec10c6aSStefano Zampini ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr); 6395db18549SStefano Zampini } 640*9087bf02SStefano Zampini } 641a1337663SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 642a1337663SStefano Zampini ierr = MatDestroy(&S_all);CHKERRQ(ierr); 643d2627357SStefano Zampini ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr); 6445db18549SStefano Zampini ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6455db18549SStefano Zampini ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6465a95e1ceSStefano Zampini if (compute_Stilda) { 647a1337663SStefano Zampini ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 648a1337663SStefano Zampini ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64908122e43SStefano Zampini ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65008122e43SStefano Zampini ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65108122e43SStefano Zampini } 652a1337663SStefano Zampini 6535db18549SStefano Zampini /* Global matrix of all assembled Schur on subsets */ 6545db18549SStefano Zampini ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 6553927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 6563927de2eSStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 6575a95e1ceSStefano Zampini 6585db18549SStefano Zampini /* Get local part of (\sum_j S_Ej) */ 6595a95e1ceSStefano Zampini ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr); 660d648f858SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 661d648f858SStefano Zampini ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 66208122e43SStefano Zampini 6635a95e1ceSStefano Zampini /* Compute explicitly (\sum_j S_Ej)^-1 */ 6645a95e1ceSStefano Zampini { 6655a95e1ceSStefano Zampini Mat tmpmat; 6665a95e1ceSStefano Zampini PetscScalar *array; 6675a95e1ceSStefano Zampini PetscInt cum; 6685a95e1ceSStefano Zampini 6695a95e1ceSStefano Zampini ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 6705a95e1ceSStefano Zampini cum = 0; 6715a95e1ceSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 6725a95e1ceSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 6735a95e1ceSStefano Zampini ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 6745a95e1ceSStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 6755a95e1ceSStefano Zampini if (!sub_schurs->is_hermitian) { 6765a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 6775a95e1ceSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 6785a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 6795a95e1ceSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 6805a95e1ceSStefano Zampini } else { 6815a95e1ceSStefano Zampini PetscInt j,k; 6825a95e1ceSStefano Zampini 6835a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 6845a95e1ceSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 6855a95e1ceSStefano Zampini PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 6865a95e1ceSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 6875a95e1ceSStefano Zampini for (j=0;j<B_N;j++) { 6885a95e1ceSStefano Zampini for (k=j+1;k<B_N;k++) { 6895a95e1ceSStefano Zampini array[k*B_N+j+cum] = array[j*B_N+k+cum]; 6905a95e1ceSStefano Zampini } 6915a95e1ceSStefano Zampini } 6925a95e1ceSStefano Zampini } 6935a95e1ceSStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 6945a95e1ceSStefano Zampini cum += subset_size*subset_size; 6955a95e1ceSStefano Zampini } 6965a95e1ceSStefano Zampini ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr); 6975a95e1ceSStefano Zampini ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr); 6985a95e1ceSStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 6995a95e1ceSStefano Zampini sub_schurs->S_Ej_all = tmpmat; 7005a95e1ceSStefano Zampini } 7015a95e1ceSStefano Zampini 7025a95e1ceSStefano Zampini /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) and invert them */ 7035a95e1ceSStefano Zampini if (compute_Stilda) { 7045a95e1ceSStefano Zampini PetscInt cum; 7055a95e1ceSStefano Zampini PetscScalar *array,*array2; 7065a95e1ceSStefano Zampini 707a1337663SStefano Zampini ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr); 708a1337663SStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 709d648f858SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 710d648f858SStefano Zampini ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 71108122e43SStefano Zampini ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr); 71208122e43SStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 713d648f858SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 714d648f858SStefano Zampini ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 7155a95e1ceSStefano Zampini /* invert blocks */ 71608122e43SStefano Zampini ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 71708122e43SStefano Zampini ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 71808122e43SStefano Zampini cum = 0; 71908122e43SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 72008122e43SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 72108122e43SStefano Zampini if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) { 72208122e43SStefano Zampini ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 72308122e43SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7242972d61bSStefano Zampini if (!sub_schurs->is_hermitian) { 72508122e43SStefano Zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 72608122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 72706a4b1faSStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 72808122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 7292972d61bSStefano Zampini } else { 7302972d61bSStefano Zampini PetscInt j,k; 7312972d61bSStefano Zampini 7322972d61bSStefano Zampini PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr)); 7332972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 7342972d61bSStefano Zampini PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr)); 7352972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 7362972d61bSStefano Zampini for (j=0;j<B_N;j++) { 7372972d61bSStefano Zampini for (k=j+1;k<B_N;k++) { 7382972d61bSStefano Zampini array[k*B_N+j+cum] = array[j*B_N+k+cum]; 7392972d61bSStefano Zampini } 7402972d61bSStefano Zampini } 7412972d61bSStefano Zampini } 7429552c7c7SStefano Zampini if (invert_Stildas) { /* Stildas can be singular */ 7432972d61bSStefano Zampini if (!sub_schurs->is_hermitian) { 74408122e43SStefano Zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr)); 74508122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 74606a4b1faSStefano Zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 74708122e43SStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 7482972d61bSStefano Zampini } else { 7492972d61bSStefano Zampini PetscInt j,k; 7502972d61bSStefano Zampini 7512972d61bSStefano Zampini PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr)); 7522972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 7532972d61bSStefano Zampini PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr)); 7542972d61bSStefano Zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 7552972d61bSStefano Zampini for (j=0;j<B_N;j++) { 7562972d61bSStefano Zampini for (k=j+1;k<B_N;k++) { 7572972d61bSStefano Zampini array2[k*B_N+j+cum] = array2[j*B_N+k+cum]; 7582972d61bSStefano Zampini } 7592972d61bSStefano Zampini } 7602972d61bSStefano Zampini } 7619552c7c7SStefano Zampini } 76208122e43SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 76308122e43SStefano Zampini } 76408122e43SStefano Zampini cum += subset_size*subset_size; 76508122e43SStefano Zampini } 76608122e43SStefano Zampini ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 76708122e43SStefano Zampini ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 76808122e43SStefano Zampini } 7693202ece2SStefano Zampini 7705a95e1ceSStefano Zampini /* free workspace */ 77106a4b1faSStefano Zampini ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr); 772a1337663SStefano Zampini ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 773a1337663SStefano Zampini ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr); 77408122e43SStefano Zampini ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr); 7753202ece2SStefano Zampini ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 7765db18549SStefano Zampini ierr = ISDestroy(&temp_is);CHKERRQ(ierr); 7775a95e1ceSStefano Zampini ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr); 778b1b3d7a2SStefano Zampini PetscFunctionReturn(0); 779b1b3d7a2SStefano Zampini } 780b1b3d7a2SStefano Zampini 781b1b3d7a2SStefano Zampini #undef __FUNCT__ 782b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit" 7835a95e1ceSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap) 784b1b3d7a2SStefano Zampini { 7859bb4a8caSStefano Zampini IS *faces,*edges,*all_cc,vertices; 7865a95e1ceSStefano Zampini PetscInt i,n_faces,n_edges,n_all_cc; 787b1b3d7a2SStefano Zampini PetscBool is_sorted; 788b1b3d7a2SStefano Zampini PetscErrorCode ierr; 789b1b3d7a2SStefano Zampini 790b1b3d7a2SStefano Zampini PetscFunctionBegin; 791b1b3d7a2SStefano Zampini ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 792b1b3d7a2SStefano Zampini if (!is_sorted) { 793b1b3d7a2SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 794b1b3d7a2SStefano Zampini } 795b1b3d7a2SStefano Zampini ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 796b1b3d7a2SStefano Zampini if (!is_sorted) { 797b1b3d7a2SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 798b1b3d7a2SStefano Zampini } 799b1b3d7a2SStefano Zampini 800b1b3d7a2SStefano Zampini /* reset any previous data */ 801b1b3d7a2SStefano Zampini ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 802b1b3d7a2SStefano Zampini 8035a95e1ceSStefano Zampini /* get index sets for faces and edges (already sorted by global ordering) */ 8049bb4a8caSStefano Zampini ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 805b1b3d7a2SStefano Zampini n_all_cc = n_faces+n_edges; 80608122e43SStefano Zampini ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 80708122e43SStefano Zampini ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 808b1b3d7a2SStefano Zampini ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 809b1b3d7a2SStefano Zampini for (i=0;i<n_faces;i++) { 810b1b3d7a2SStefano Zampini all_cc[i] = faces[i]; 811b1b3d7a2SStefano Zampini } 812b1b3d7a2SStefano Zampini for (i=0;i<n_edges;i++) { 813b1b3d7a2SStefano Zampini all_cc[n_faces+i] = edges[i]; 81408122e43SStefano Zampini ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 815b1b3d7a2SStefano Zampini } 816b1b3d7a2SStefano Zampini ierr = PetscFree(faces);CHKERRQ(ierr); 817b1b3d7a2SStefano Zampini ierr = PetscFree(edges);CHKERRQ(ierr); 818b1b3d7a2SStefano Zampini 8195a95e1ceSStefano Zampini /* Determine if MUMPS cane be used */ 820883469d8SStefano Zampini sub_schurs->use_mumps = PETSC_FALSE; 821883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 8224c6709b3SStefano Zampini sub_schurs->use_mumps = (PetscBool)(!!A); 823883469d8SStefano Zampini #endif 824b1b3d7a2SStefano Zampini 825b1b3d7a2SStefano Zampini /* update info in sub_schurs */ 826aa83b6aeSStefano Zampini if (A) { 8279ab7bb16SStefano Zampini PetscBool isseqaij; 8289ab7bb16SStefano Zampini 8299ab7bb16SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 8309ab7bb16SStefano Zampini if (isseqaij) { 8311e9c79c2SStefano Zampini ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8321e9c79c2SStefano Zampini sub_schurs->A = A; 8339ab7bb16SStefano Zampini } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */ 8349ab7bb16SStefano Zampini ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr); 8359ab7bb16SStefano Zampini } 8361e9c79c2SStefano Zampini } 837b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr); 838b1b3d7a2SStefano Zampini sub_schurs->S = S; 839b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 840b1b3d7a2SStefano Zampini sub_schurs->is_I = is_I; 841b1b3d7a2SStefano Zampini ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 842b1b3d7a2SStefano Zampini sub_schurs->is_B = is_B; 8435db18549SStefano Zampini ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 8445db18549SStefano Zampini sub_schurs->l2gmap = graph->l2gmap; 8455db18549SStefano Zampini ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 8465db18549SStefano Zampini sub_schurs->BtoNmap = BtoNmap; 8475a95e1ceSStefano Zampini sub_schurs->n_subs = n_all_cc; 848b1b3d7a2SStefano Zampini sub_schurs->is_subs = all_cc; 8499bb4a8caSStefano Zampini if (!sub_schurs->use_mumps) { /* for adaptive selection */ 850b96c3477SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 851b96c3477SStefano Zampini ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 852b96c3477SStefano Zampini } 8539bb4a8caSStefano Zampini } 8549bb4a8caSStefano Zampini sub_schurs->is_Ej_com = vertices; 855b1b3d7a2SStefano Zampini 856b96c3477SStefano Zampini 857b96c3477SStefano Zampini /* allocate space for schur complements */ 858b96c3477SStefano Zampini sub_schurs->S_Ej_all = NULL; 859b96c3477SStefano Zampini sub_schurs->sum_S_Ej_all = NULL; 86008122e43SStefano Zampini sub_schurs->sum_S_Ej_inv_all = NULL; 861b96c3477SStefano Zampini sub_schurs->sum_S_Ej_tilda_all = NULL; 862b96c3477SStefano Zampini sub_schurs->is_Ej_all = NULL; 863b1b3d7a2SStefano Zampini PetscFunctionReturn(0); 864b1b3d7a2SStefano Zampini } 865b1b3d7a2SStefano Zampini 866b1b3d7a2SStefano Zampini #undef __FUNCT__ 86734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate" 86834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 86934a97f8cSStefano Zampini { 87034a97f8cSStefano Zampini PCBDDCSubSchurs schurs_ctx; 87134a97f8cSStefano Zampini PetscErrorCode ierr; 87234a97f8cSStefano Zampini 87334a97f8cSStefano Zampini PetscFunctionBegin; 87434a97f8cSStefano Zampini ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 8755ff63025SStefano Zampini schurs_ctx->n_subs = 0; 87634a97f8cSStefano Zampini *sub_schurs = schurs_ctx; 87734a97f8cSStefano Zampini PetscFunctionReturn(0); 87834a97f8cSStefano Zampini } 87934a97f8cSStefano Zampini 88034a97f8cSStefano Zampini #undef __FUNCT__ 88134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy" 88234a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 88334a97f8cSStefano Zampini { 88434a97f8cSStefano Zampini PetscErrorCode ierr; 88534a97f8cSStefano Zampini 88634a97f8cSStefano Zampini PetscFunctionBegin; 88734a97f8cSStefano Zampini ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 88834a97f8cSStefano Zampini ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 88934a97f8cSStefano Zampini PetscFunctionReturn(0); 89034a97f8cSStefano Zampini } 89134a97f8cSStefano Zampini 89234a97f8cSStefano Zampini #undef __FUNCT__ 89334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset" 89434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 89534a97f8cSStefano Zampini { 89634a97f8cSStefano Zampini PetscInt i; 89734a97f8cSStefano Zampini PetscErrorCode ierr; 89834a97f8cSStefano Zampini 89934a97f8cSStefano Zampini PetscFunctionBegin; 9001e9c79c2SStefano Zampini ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 901b1b3d7a2SStefano Zampini ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 902b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 903b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 9045db18549SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 9055db18549SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 90641c3ba1bSStefano Zampini ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 90741c3ba1bSStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 90808122e43SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 909a1337663SStefano Zampini ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 9105db18549SStefano Zampini ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 9119bb4a8caSStefano Zampini ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr); 91208122e43SStefano Zampini ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 91308122e43SStefano Zampini ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 91434a97f8cSStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 915b1b3d7a2SStefano Zampini ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 91634a97f8cSStefano Zampini } 91768270318SStefano Zampini ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 9185ff63025SStefano Zampini if (sub_schurs->n_subs) { 919b1b3d7a2SStefano Zampini ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 9203dc780c3SStefano Zampini } 92134a97f8cSStefano Zampini sub_schurs->n_subs = 0; 92234a97f8cSStefano Zampini PetscFunctionReturn(0); 92334a97f8cSStefano Zampini } 92434a97f8cSStefano Zampini 92534a97f8cSStefano Zampini #undef __FUNCT__ 92634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 9272a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 92834a97f8cSStefano Zampini { 92934a97f8cSStefano Zampini PetscInt i,j,n; 93034a97f8cSStefano Zampini PetscErrorCode ierr; 93134a97f8cSStefano Zampini 93234a97f8cSStefano Zampini PetscFunctionBegin; 93334a97f8cSStefano Zampini n = 0; 93434a97f8cSStefano Zampini for (i=-n_prev;i<0;i++) { 93534a97f8cSStefano Zampini PetscInt start_dof = queue_tip[i]; 93634a97f8cSStefano Zampini for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 93734a97f8cSStefano Zampini PetscInt dof = adjncy[j]; 93834a97f8cSStefano Zampini if (!PetscBTLookup(touched,dof)) { 93934a97f8cSStefano Zampini ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 94034a97f8cSStefano Zampini queue_tip[n] = dof; 94134a97f8cSStefano Zampini n++; 94234a97f8cSStefano Zampini } 94334a97f8cSStefano Zampini } 94434a97f8cSStefano Zampini } 94534a97f8cSStefano Zampini *n_added = n; 94634a97f8cSStefano Zampini PetscFunctionReturn(0); 94734a97f8cSStefano Zampini } 948