14b9ad928SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/is/nn/nn.h> 34b9ad928SBarry Smith 44b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 54b9ad928SBarry Smith /* 64b9ad928SBarry Smith PCSetUp_NN - Prepares for the use of the NN preconditioner 74b9ad928SBarry Smith by setting data structures and options. 84b9ad928SBarry Smith 94b9ad928SBarry Smith Input Parameter: 104b9ad928SBarry Smith . pc - the preconditioner context 114b9ad928SBarry Smith 124b9ad928SBarry Smith Application Interface Routine: PCSetUp() 134b9ad928SBarry Smith 14f1580f4eSBarry Smith Note: 154b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 164b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 174b9ad928SBarry Smith */ 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_NN(PC pc) 19d71ae5a4SJacob Faibussowitsch { 204b9ad928SBarry Smith PetscFunctionBegin; 214b9ad928SBarry Smith if (!pc->setupcalled) { 224b9ad928SBarry Smith /* Set up all the "iterative substructuring" common block */ 239566063dSJacob Faibussowitsch PetscCall(PCISSetUp(pc, PETSC_TRUE, PETSC_TRUE)); 244b9ad928SBarry Smith /* Create the coarse matrix. */ 259566063dSJacob Faibussowitsch PetscCall(PCNNCreateCoarseMatrix(pc)); 264b9ad928SBarry Smith } 273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284b9ad928SBarry Smith } 294b9ad928SBarry Smith 304b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 314b9ad928SBarry Smith /* 324b9ad928SBarry Smith PCApply_NN - Applies the NN preconditioner to a vector. 334b9ad928SBarry Smith 344b9ad928SBarry Smith Input Parameters: 35*2fe279fdSBarry Smith + pc - the preconditioner context 36*2fe279fdSBarry Smith - r - input vector (global) 374b9ad928SBarry Smith 384b9ad928SBarry Smith Output Parameter: 394b9ad928SBarry Smith . z - output vector (global) 404b9ad928SBarry Smith 414b9ad928SBarry Smith Application Interface Routine: PCApply() 424b9ad928SBarry Smith */ 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_NN(PC pc, Vec r, Vec z) 44d71ae5a4SJacob Faibussowitsch { 454b9ad928SBarry Smith PC_IS *pcis = (PC_IS *)(pc->data); 464b9ad928SBarry Smith PetscScalar m_one = -1.0; 474b9ad928SBarry Smith Vec w = pcis->vec1_global; 484b9ad928SBarry Smith 494b9ad928SBarry Smith PetscFunctionBegin; 504b9ad928SBarry Smith /* 514b9ad928SBarry Smith Dirichlet solvers. 524b9ad928SBarry Smith Solving $ B_I^{(i)}r_I^{(i)} $ at each processor. 534b9ad928SBarry Smith Storing the local results at vec2_D 544b9ad928SBarry Smith */ 559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 569566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 579566063dSJacob Faibussowitsch PetscCall(KSPSolve(pcis->ksp_D, pcis->vec1_D, pcis->vec2_D)); 584b9ad928SBarry Smith 594b9ad928SBarry Smith /* 604b9ad928SBarry Smith Computing $ r_B - \sum_j \tilde R_j^T A_{BI}^{(j)} (B_I^{(j)}r_I^{(j)}) $ . 614b9ad928SBarry Smith Storing the result in the interface portion of the global vector w. 624b9ad928SBarry Smith */ 639566063dSJacob Faibussowitsch PetscCall(MatMult(pcis->A_BI, pcis->vec2_D, pcis->vec1_B)); 649566063dSJacob Faibussowitsch PetscCall(VecScale(pcis->vec1_B, m_one)); 659566063dSJacob Faibussowitsch PetscCall(VecCopy(r, w)); 669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, w, ADD_VALUES, SCATTER_REVERSE)); 679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, w, ADD_VALUES, SCATTER_REVERSE)); 684b9ad928SBarry Smith 694b9ad928SBarry Smith /* 704b9ad928SBarry Smith Apply the interface preconditioner 714b9ad928SBarry Smith */ 72d0609cedSBarry Smith PetscCall(PCNNApplyInterfacePreconditioner(pc, w, z, pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec3_B, pcis->vec1_D, pcis->vec3_D, pcis->vec1_N, pcis->vec2_N)); 734b9ad928SBarry Smith 744b9ad928SBarry Smith /* 754b9ad928SBarry Smith Computing $ t_I^{(i)} = A_{IB}^{(i)} \tilde R_i z_B $ 764b9ad928SBarry Smith The result is stored in vec1_D. 774b9ad928SBarry Smith */ 789566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 799566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 809566063dSJacob Faibussowitsch PetscCall(MatMult(pcis->A_IB, pcis->vec1_B, pcis->vec1_D)); 814b9ad928SBarry Smith 824b9ad928SBarry Smith /* 834b9ad928SBarry Smith Dirichlet solvers. 844b9ad928SBarry Smith Computing $ B_I^{(i)}t_I^{(i)} $ and sticking into the global vector the blocks 854b9ad928SBarry Smith $ B_I^{(i)}r_I^{(i)} - B_I^{(i)}t_I^{(i)} $. 864b9ad928SBarry Smith */ 879566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 889566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 899566063dSJacob Faibussowitsch PetscCall(KSPSolve(pcis->ksp_D, pcis->vec1_D, pcis->vec2_D)); 909566063dSJacob Faibussowitsch PetscCall(VecScale(pcis->vec2_D, m_one)); 919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, ADD_VALUES, SCATTER_REVERSE)); 929566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, ADD_VALUES, SCATTER_REVERSE)); 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 944b9ad928SBarry Smith } 954b9ad928SBarry Smith 964b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 974b9ad928SBarry Smith /* 984b9ad928SBarry Smith PCDestroy_NN - Destroys the private context for the NN preconditioner 994b9ad928SBarry Smith that was created with PCCreate_NN(). 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith Input Parameter: 1024b9ad928SBarry Smith . pc - the preconditioner context 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith Application Interface Routine: PCDestroy() 1054b9ad928SBarry Smith */ 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_NN(PC pc) 107d71ae5a4SJacob Faibussowitsch { 1084b9ad928SBarry Smith PC_NN *pcnn = (PC_NN *)pc->data; 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith PetscFunctionBegin; 1119566063dSJacob Faibussowitsch PetscCall(PCISDestroy(pc)); 1124b9ad928SBarry Smith 1139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pcnn->coarse_mat)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pcnn->coarse_x)); 1159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pcnn->coarse_b)); 1169566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&pcnn->ksp_coarse)); 1174b9ad928SBarry Smith if (pcnn->DZ_IN) { 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(pcnn->DZ_IN[0])); 1199566063dSJacob Faibussowitsch PetscCall(PetscFree(pcnn->DZ_IN)); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith /* 1234b9ad928SBarry Smith Free the private data structure that was hanging off the PC 1244b9ad928SBarry Smith */ 1259566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1274b9ad928SBarry Smith } 1284b9ad928SBarry Smith 129dfc58137SBarry Smith /*MC 13037a17b4dSBarry Smith PCNN - Balancing Neumann-Neumann for scalar elliptic PDEs. 1314b9ad928SBarry Smith 13237a17b4dSBarry Smith Options Database Keys: 13337a17b4dSBarry Smith + -pc_nn_turn_off_first_balancing - do not balance the residual before solving the local Neumann problems 13437a17b4dSBarry Smith (this skips the first coarse grid solve in the preconditioner) 1355c9740d6SBarry Smith . -pc_nn_turn_off_second_balancing - do not balance the solution solving the local Neumann problems 13637a17b4dSBarry Smith (this skips the second coarse grid solve in the preconditioner) 1375c9740d6SBarry Smith . -pc_is_damp_fixed <fact> - 1385c9740d6SBarry Smith . -pc_is_remove_nullspace_fixed - 1395c9740d6SBarry Smith . -pc_is_set_damping_factor_floating <fact> - 1405c9740d6SBarry Smith . -pc_is_not_damp_floating - 141a2b725a8SWilliam Gropp - -pc_is_not_remove_nullspace_floating - 1424b9ad928SBarry Smith 143f1580f4eSBarry Smith Options Database prefixes for the subsolvers this preconditioner uses: 144f1580f4eSBarry Smith + -nn_coarse_pc_ - for the coarse grid preconditioner 145f1580f4eSBarry Smith . -is_localD_pc_ - for the Dirichlet subproblem preconditioner 146f1580f4eSBarry Smith - -is_localN_pc_ - for the Neumann subproblem preconditioner 147f1580f4eSBarry Smith 1488eb2139fSSatish Balay Level: intermediate 14937a17b4dSBarry Smith 15095452b02SPatrick Sanan Notes: 151f1580f4eSBarry Smith The matrix used with this preconditioner must be of type `MATIS` 15237a17b4dSBarry Smith 1531b99c236SBarry Smith Unlike more 'conventional' Neumann-Neumann preconditioners this iterates over ALL the 1541b99c236SBarry Smith degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers 1551b99c236SBarry Smith on the subdomains; though in our experience using approximate solvers is slower.). 1561b99c236SBarry Smith 15737a17b4dSBarry Smith Contributed by Paulo Goldfeld 15837a17b4dSBarry Smith 159f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MATIS`, `PCBDDC` 16037a17b4dSBarry Smith M*/ 161b2573a8aSBarry Smith 162d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_NN(PC pc) 163d71ae5a4SJacob Faibussowitsch { 1644b9ad928SBarry Smith PC_NN *pcnn; 1654b9ad928SBarry Smith 1664b9ad928SBarry Smith PetscFunctionBegin; 1674b9ad928SBarry Smith /* 1684b9ad928SBarry Smith Creates the private data structure for this preconditioner and 1694b9ad928SBarry Smith attach it to the PC object. 1704b9ad928SBarry Smith */ 1714dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pcnn)); 1724b9ad928SBarry Smith pc->data = (void *)pcnn; 1734b9ad928SBarry Smith 1749566063dSJacob Faibussowitsch PetscCall(PCISCreate(pc)); 1750a545947SLisandro Dalcin pcnn->coarse_mat = NULL; 1760a545947SLisandro Dalcin pcnn->coarse_x = NULL; 1770a545947SLisandro Dalcin pcnn->coarse_b = NULL; 1780a545947SLisandro Dalcin pcnn->ksp_coarse = NULL; 1790a545947SLisandro Dalcin pcnn->DZ_IN = NULL; 1804b9ad928SBarry Smith 1814b9ad928SBarry Smith /* 1824b9ad928SBarry Smith Set the pointers for the functions that are provided above. 1834b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 1844b9ad928SBarry Smith are called, they will automatically call these functions. Note we 1854b9ad928SBarry Smith choose not to provide a couple of these functions since they are 1864b9ad928SBarry Smith not needed. 1874b9ad928SBarry Smith */ 1884b9ad928SBarry Smith pc->ops->apply = PCApply_NN; 1890a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 1904b9ad928SBarry Smith pc->ops->setup = PCSetUp_NN; 1914b9ad928SBarry Smith pc->ops->destroy = PCDestroy_NN; 1920a545947SLisandro Dalcin pc->ops->view = NULL; 1930a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1940a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 1950a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1974b9ad928SBarry Smith } 1984b9ad928SBarry Smith 1994b9ad928SBarry Smith /* 2004b9ad928SBarry Smith PCNNCreateCoarseMatrix - 2014b9ad928SBarry Smith */ 202d71ae5a4SJacob Faibussowitsch PetscErrorCode PCNNCreateCoarseMatrix(PC pc) 203d71ae5a4SJacob Faibussowitsch { 2044b9ad928SBarry Smith MPI_Request *send_request, *recv_request; 20513f74950SBarry Smith PetscInt i, j, k; 2064b9ad928SBarry Smith PetscScalar *mat; /* Sub-matrix with this subdomain's contribution to the coarse matrix */ 2074b9ad928SBarry Smith PetscScalar **DZ_OUT; /* proc[k].DZ_OUT[i][] = bit of vector to be sent from processor k to processor i */ 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith /* aliasing some names */ 2104b9ad928SBarry Smith PC_IS *pcis = (PC_IS *)(pc->data); 2114b9ad928SBarry Smith PC_NN *pcnn = (PC_NN *)pc->data; 21213f74950SBarry Smith PetscInt n_neigh = pcis->n_neigh; 21313f74950SBarry Smith PetscInt *neigh = pcis->neigh; 21413f74950SBarry Smith PetscInt *n_shared = pcis->n_shared; 21513f74950SBarry Smith PetscInt **shared = pcis->shared; 2164b9ad928SBarry Smith PetscScalar **DZ_IN; /* Must be initialized after memory allocation. */ 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith PetscFunctionBegin; 2194b9ad928SBarry Smith /* Allocate memory for mat (the +1 is to handle the case n_neigh equal to zero) */ 2209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_neigh * n_neigh + 1, &mat)); 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith /* Allocate memory for DZ */ 2234b9ad928SBarry Smith /* Notice that DZ_OUT[0] is allocated some space that is never used. */ 2244b9ad928SBarry Smith /* This is just in order to DZ_OUT and DZ_IN to have exactly the same form. */ 2254b9ad928SBarry Smith { 22613f74950SBarry Smith PetscInt size_of_Z = 0; 2279566063dSJacob Faibussowitsch PetscCall(PetscMalloc((n_neigh + 1) * sizeof(PetscScalar *), &pcnn->DZ_IN)); 2284b9ad928SBarry Smith DZ_IN = pcnn->DZ_IN; 2299566063dSJacob Faibussowitsch PetscCall(PetscMalloc((n_neigh + 1) * sizeof(PetscScalar *), &DZ_OUT)); 2302fa5cd67SKarl Rupp for (i = 0; i < n_neigh; i++) size_of_Z += n_shared[i]; 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc((size_of_Z + 1) * sizeof(PetscScalar), &DZ_IN[0])); 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc((size_of_Z + 1) * sizeof(PetscScalar), &DZ_OUT[0])); 2334b9ad928SBarry Smith } 2344b9ad928SBarry Smith for (i = 1; i < n_neigh; i++) { 2354b9ad928SBarry Smith DZ_IN[i] = DZ_IN[i - 1] + n_shared[i - 1]; 2364b9ad928SBarry Smith DZ_OUT[i] = DZ_OUT[i - 1] + n_shared[i - 1]; 2374b9ad928SBarry Smith } 2384b9ad928SBarry Smith 2394b9ad928SBarry Smith /* Set the values of DZ_OUT, in order to send this info to the neighbours */ 2404b9ad928SBarry Smith /* First, set the auxiliary array pcis->work_N. */ 2419566063dSJacob Faibussowitsch PetscCall(PCISScatterArrayNToVecB(pcis->work_N, pcis->D, INSERT_VALUES, SCATTER_REVERSE, pc)); 2424b9ad928SBarry Smith for (i = 1; i < n_neigh; i++) { 243ad540459SPierre Jolivet for (j = 0; j < n_shared[i]; j++) DZ_OUT[i][j] = pcis->work_N[shared[i][j]]; 2444b9ad928SBarry Smith } 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith /* Non-blocking send/receive the common-interface chunks of scaled nullspaces */ 2474b9ad928SBarry Smith /* Notice that send_request[] and recv_request[] could have one less element. */ 2484b9ad928SBarry Smith /* We make them longer to have request[i] corresponding to neigh[i]. */ 2494b9ad928SBarry Smith { 25013f74950SBarry Smith PetscMPIInt tag; 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)pc, &tag)); 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n_neigh + 1, &send_request, n_neigh + 1, &recv_request)); 2534b9ad928SBarry Smith for (i = 1; i < n_neigh; i++) { 2549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend((void *)(DZ_OUT[i]), n_shared[i], MPIU_SCALAR, neigh[i], tag, PetscObjectComm((PetscObject)pc), &(send_request[i]))); 2559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv((void *)(DZ_IN[i]), n_shared[i], MPIU_SCALAR, neigh[i], tag, PetscObjectComm((PetscObject)pc), &(recv_request[i]))); 2564b9ad928SBarry Smith } 2574b9ad928SBarry Smith } 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith /* Set DZ_IN[0][] (recall that neigh[0]==rank, always) */ 2602fa5cd67SKarl Rupp for (j = 0; j < n_shared[0]; j++) DZ_IN[0][j] = pcis->work_N[shared[0][j]]; 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith /* Start computing with local D*Z while communication goes on. */ 2634b9ad928SBarry Smith /* Apply Schur complement. The result is "stored" in vec (more */ 2644b9ad928SBarry Smith /* precisely, vec points to the result, stored in pc_nn->vec1_B) */ 2654b9ad928SBarry Smith /* and also scattered to pcnn->work_N. */ 2669566063dSJacob Faibussowitsch PetscCall(PCNNApplySchurToChunk(pc, n_shared[0], shared[0], DZ_IN[0], pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec1_D, pcis->vec2_D)); 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith /* Compute the first column, while completing the receiving. */ 2694b9ad928SBarry Smith for (i = 0; i < n_neigh; i++) { 2704b9ad928SBarry Smith MPI_Status stat; 27113f74950SBarry Smith PetscMPIInt ind = 0; 2729371c9d4SSatish Balay if (i > 0) { 2739371c9d4SSatish Balay PetscCallMPI(MPI_Waitany(n_neigh - 1, recv_request + 1, &ind, &stat)); 2749371c9d4SSatish Balay ind++; 2759371c9d4SSatish Balay } 2764b9ad928SBarry Smith mat[ind * n_neigh + 0] = 0.0; 2772fa5cd67SKarl Rupp for (k = 0; k < n_shared[ind]; k++) mat[ind * n_neigh + 0] += DZ_IN[ind][k] * pcis->work_N[shared[ind][k]]; 2784b9ad928SBarry Smith } 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith /* Compute the remaining of the columns */ 2814b9ad928SBarry Smith for (j = 1; j < n_neigh; j++) { 2829566063dSJacob Faibussowitsch PetscCall(PCNNApplySchurToChunk(pc, n_shared[j], shared[j], DZ_IN[j], pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec1_D, pcis->vec2_D)); 2834b9ad928SBarry Smith for (i = 0; i < n_neigh; i++) { 2844b9ad928SBarry Smith mat[i * n_neigh + j] = 0.0; 2852fa5cd67SKarl Rupp for (k = 0; k < n_shared[i]; k++) mat[i * n_neigh + j] += DZ_IN[i][k] * pcis->work_N[shared[i][k]]; 2864b9ad928SBarry Smith } 2874b9ad928SBarry Smith } 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith /* Complete the sending. */ 2904b9ad928SBarry Smith if (n_neigh > 1) { 2914b9ad928SBarry Smith MPI_Status *stat; 2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_neigh - 1, &stat)); 2939566063dSJacob Faibussowitsch if (n_neigh - 1) PetscCallMPI(MPI_Waitall(n_neigh - 1, &(send_request[1]), stat)); 2949566063dSJacob Faibussowitsch PetscCall(PetscFree(stat)); 2954b9ad928SBarry Smith } 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith /* Free the memory for the MPI requests */ 2989566063dSJacob Faibussowitsch PetscCall(PetscFree2(send_request, recv_request)); 2994b9ad928SBarry Smith 3004b9ad928SBarry Smith /* Free the memory for DZ_OUT */ 3014b9ad928SBarry Smith if (DZ_OUT) { 3029566063dSJacob Faibussowitsch PetscCall(PetscFree(DZ_OUT[0])); 3039566063dSJacob Faibussowitsch PetscCall(PetscFree(DZ_OUT)); 3044b9ad928SBarry Smith } 3054b9ad928SBarry Smith 3064b9ad928SBarry Smith { 30713f74950SBarry Smith PetscMPIInt size; 3089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 3094b9ad928SBarry Smith /* Create the global coarse vectors (rhs and solution). */ 3109566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc), 1, size, &(pcnn->coarse_b))); 3119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pcnn->coarse_b, &(pcnn->coarse_x))); 312f204ca49SKris Buschelman /* Create and set the global coarse AIJ matrix. */ 3139566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &(pcnn->coarse_mat))); 3149566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pcnn->coarse_mat, 1, 1, size, size)); 3159566063dSJacob Faibussowitsch PetscCall(MatSetType(pcnn->coarse_mat, MATAIJ)); 3169566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pcnn->coarse_mat, 1, NULL)); 3179566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pcnn->coarse_mat, 1, NULL, n_neigh, NULL)); 3189566063dSJacob Faibussowitsch PetscCall(MatSetOption(pcnn->coarse_mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 3199566063dSJacob Faibussowitsch PetscCall(MatSetOption(pcnn->coarse_mat, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE)); 3209566063dSJacob Faibussowitsch PetscCall(MatSetValues(pcnn->coarse_mat, n_neigh, neigh, n_neigh, neigh, mat, ADD_VALUES)); 3219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pcnn->coarse_mat, MAT_FINAL_ASSEMBLY)); 3229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pcnn->coarse_mat, MAT_FINAL_ASSEMBLY)); 3234b9ad928SBarry Smith } 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith { 32613f74950SBarry Smith PetscMPIInt rank; 3274b9ad928SBarry Smith PetscScalar one = 1.0; 3289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 3294b9ad928SBarry Smith /* "Zero out" rows of not-purely-Neumann subdomains */ 3304b9ad928SBarry Smith if (pcis->pure_neumann) { /* does NOT zero the row; create an empty index set. The reason is that MatZeroRows() is collective. */ 3319566063dSJacob Faibussowitsch PetscCall(MatZeroRows(pcnn->coarse_mat, 0, NULL, one, NULL, NULL)); 3324b9ad928SBarry Smith } else { /* here it DOES zero the row, since it's not a floating subdomain. */ 3331c82e0cfSMatthew Knepley PetscInt row = (PetscInt)rank; 3349566063dSJacob Faibussowitsch PetscCall(MatZeroRows(pcnn->coarse_mat, 1, &row, one, NULL, NULL)); 3354b9ad928SBarry Smith } 3364b9ad928SBarry Smith } 3374b9ad928SBarry Smith 3384b9ad928SBarry Smith /* Create the coarse linear solver context */ 3394b9ad928SBarry Smith { 3404b9ad928SBarry Smith PC pc_ctx, inner_pc; 34183ab6a24SBarry Smith KSP inner_ksp; 34283ab6a24SBarry Smith 3439566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &pcnn->ksp_coarse)); 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcnn->ksp_coarse, (PetscObject)pc, 2)); 3459566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(pcnn->ksp_coarse, pcnn->coarse_mat, pcnn->coarse_mat)); 3469566063dSJacob Faibussowitsch PetscCall(KSPGetPC(pcnn->ksp_coarse, &pc_ctx)); 3479566063dSJacob Faibussowitsch PetscCall(PCSetType(pc_ctx, PCREDUNDANT)); 3489566063dSJacob Faibussowitsch PetscCall(KSPSetType(pcnn->ksp_coarse, KSPPREONLY)); 3499566063dSJacob Faibussowitsch PetscCall(PCRedundantGetKSP(pc_ctx, &inner_ksp)); 3509566063dSJacob Faibussowitsch PetscCall(KSPGetPC(inner_ksp, &inner_pc)); 3519566063dSJacob Faibussowitsch PetscCall(PCSetType(inner_pc, PCLU)); 3529566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(pcnn->ksp_coarse, "nn_coarse_")); 3539566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(pcnn->ksp_coarse)); 3544b9ad928SBarry Smith /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 3559566063dSJacob Faibussowitsch PetscCall(KSPSetUp(pcnn->ksp_coarse)); 3564b9ad928SBarry Smith } 3574b9ad928SBarry Smith 3584b9ad928SBarry Smith /* Free the memory for mat */ 3599566063dSJacob Faibussowitsch PetscCall(PetscFree(mat)); 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith /* for DEBUGGING, save the coarse matrix to a file. */ 3624b9ad928SBarry Smith { 363ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 3649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_nn_save_coarse_matrix", &flg, NULL)); 3654b9ad928SBarry Smith if (flg) { 3664b9ad928SBarry Smith PetscViewer viewer; 3679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "coarse.m", &viewer)); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); 3699566063dSJacob Faibussowitsch PetscCall(MatView(pcnn->coarse_mat, viewer)); 3709566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 3719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 3724b9ad928SBarry Smith } 3734b9ad928SBarry Smith } 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith /* Set the variable pcnn->factor_coarse_rhs. */ 3764b9ad928SBarry Smith pcnn->factor_coarse_rhs = (pcis->pure_neumann) ? 1.0 : 0.0; 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith /* See historical note 02, at the bottom of this file. */ 3793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3804b9ad928SBarry Smith } 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith /* 3834b9ad928SBarry Smith PCNNApplySchurToChunk - 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Input parameters: 3864b9ad928SBarry Smith . pcnn 3874b9ad928SBarry Smith . n - size of chunk 3884b9ad928SBarry Smith . idx - indices of chunk 3894b9ad928SBarry Smith . chunk - values 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith Output parameters: 3924b9ad928SBarry Smith . array_N - result of Schur complement applied to chunk, scattered to big array 3934b9ad928SBarry Smith . vec1_B - result of Schur complement applied to chunk 3944b9ad928SBarry Smith . vec2_B - garbage (used as work space) 3954b9ad928SBarry Smith . vec1_D - garbage (used as work space) 3964b9ad928SBarry Smith . vec2_D - garbage (used as work space) 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith */ 399d71ae5a4SJacob Faibussowitsch PetscErrorCode PCNNApplySchurToChunk(PC pc, PetscInt n, PetscInt *idx, PetscScalar *chunk, PetscScalar *array_N, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 400d71ae5a4SJacob Faibussowitsch { 40113f74950SBarry Smith PetscInt i; 4024b9ad928SBarry Smith PC_IS *pcis = (PC_IS *)(pc->data); 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith PetscFunctionBegin; 4059566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(array_N, pcis->n)); 4062fa5cd67SKarl Rupp for (i = 0; i < n; i++) array_N[idx[i]] = chunk[i]; 4079566063dSJacob Faibussowitsch PetscCall(PCISScatterArrayNToVecB(array_N, vec2_B, INSERT_VALUES, SCATTER_FORWARD, pc)); 4089566063dSJacob Faibussowitsch PetscCall(PCISApplySchur(pc, vec2_B, vec1_B, (Vec)0, vec1_D, vec2_D)); 4099566063dSJacob Faibussowitsch PetscCall(PCISScatterArrayNToVecB(array_N, vec1_B, INSERT_VALUES, SCATTER_REVERSE, pc)); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4114b9ad928SBarry Smith } 4124b9ad928SBarry Smith 4134b9ad928SBarry Smith /* 4144b9ad928SBarry Smith PCNNApplyInterfacePreconditioner - Apply the interface preconditioner, i.e., 4154b9ad928SBarry Smith the preconditioner for the Schur complement. 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith Input parameter: 4184b9ad928SBarry Smith . r - global vector of interior and interface nodes. The values on the interior nodes are NOT used. 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith Output parameters: 4214b9ad928SBarry Smith . z - global vector of interior and interface nodes. The values on the interface are the result of 4224b9ad928SBarry Smith the application of the interface preconditioner to the interface part of r. The values on the 4234b9ad928SBarry Smith interior nodes are garbage. 4244b9ad928SBarry Smith . work_N - array of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 4254b9ad928SBarry Smith . vec1_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 4264b9ad928SBarry Smith . vec2_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 4274b9ad928SBarry Smith . vec3_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 4284b9ad928SBarry Smith . vec1_D - vector of local interior nodes; returns garbage (used as work space) 4294b9ad928SBarry Smith . vec2_D - vector of local interior nodes; returns garbage (used as work space) 4304b9ad928SBarry Smith . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 4314b9ad928SBarry Smith . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith */ 434d71ae5a4SJacob Faibussowitsch PetscErrorCode PCNNApplyInterfacePreconditioner(PC pc, Vec r, Vec z, PetscScalar *work_N, Vec vec1_B, Vec vec2_B, Vec vec3_B, Vec vec1_D, Vec vec2_D, Vec vec1_N, Vec vec2_N) 435d71ae5a4SJacob Faibussowitsch { 4364b9ad928SBarry Smith PC_IS *pcis = (PC_IS *)(pc->data); 4374b9ad928SBarry Smith 4384b9ad928SBarry Smith PetscFunctionBegin; 4394b9ad928SBarry Smith /* 4404b9ad928SBarry Smith First balancing step. 4414b9ad928SBarry Smith */ 4424b9ad928SBarry Smith { 443ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 4449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_nn_turn_off_first_balancing", &flg, NULL)); 4454b9ad928SBarry Smith if (!flg) { 4469566063dSJacob Faibussowitsch PetscCall(PCNNBalancing(pc, r, (Vec)0, z, vec1_B, vec2_B, (Vec)0, vec1_D, vec2_D, work_N)); 4474b9ad928SBarry Smith } else { 4489566063dSJacob Faibussowitsch PetscCall(VecCopy(r, z)); 4494b9ad928SBarry Smith } 4504b9ad928SBarry Smith } 4514b9ad928SBarry Smith 4524b9ad928SBarry Smith /* 4534b9ad928SBarry Smith Extract the local interface part of z and scale it by D 4544b9ad928SBarry Smith */ 4559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, z, vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 4569566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, z, vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 4579566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(vec2_B, pcis->D, vec1_B)); 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith /* Neumann Solver */ 4609566063dSJacob Faibussowitsch PetscCall(PCISApplyInvSchur(pc, vec2_B, vec1_B, vec1_N, vec2_N)); 4614b9ad928SBarry Smith 4624b9ad928SBarry Smith /* 4634b9ad928SBarry Smith Second balancing step. 4644b9ad928SBarry Smith */ 4654b9ad928SBarry Smith { 466ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 4679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_turn_off_second_balancing", &flg, NULL)); 4684b9ad928SBarry Smith if (!flg) { 4699566063dSJacob Faibussowitsch PetscCall(PCNNBalancing(pc, r, vec1_B, z, vec2_B, vec3_B, (Vec)0, vec1_D, vec2_D, work_N)); 4704b9ad928SBarry Smith } else { 4719566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(vec2_B, pcis->D, vec1_B)); 4729566063dSJacob Faibussowitsch PetscCall(VecSet(z, 0.0)); 4739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE)); 4749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE)); 4754b9ad928SBarry Smith } 4764b9ad928SBarry Smith } 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4784b9ad928SBarry Smith } 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith /* 4814b9ad928SBarry Smith PCNNBalancing - Computes z, as given in equations (15) and (16) (if the 4824b9ad928SBarry Smith input argument u is provided), or s, as given in equations 4834b9ad928SBarry Smith (12) and (13), if the input argument u is a null vector. 4844b9ad928SBarry Smith Notice that the input argument u plays the role of u_i in 4854b9ad928SBarry Smith equation (14). The equation numbers refer to [Man93]. 4864b9ad928SBarry Smith 4874b9ad928SBarry Smith Input Parameters: 488*2fe279fdSBarry Smith + pcnn - NN preconditioner context. 4894b9ad928SBarry Smith . r - MPI vector of all nodes (interior and interface). It's preserved. 490*2fe279fdSBarry Smith - u - (Optional) sequential vector of local interface nodes. It's preserved UNLESS vec3_B is null. 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith Output Parameters: 493*2fe279fdSBarry Smith + z - MPI vector of interior and interface nodes. Returns s or z (see description above). 4944b9ad928SBarry Smith . vec1_B - Sequential vector of local interface nodes. Workspace. 4954b9ad928SBarry Smith . vec2_B - Sequential vector of local interface nodes. Workspace. 4964b9ad928SBarry Smith . vec3_B - (Optional) sequential vector of local interface nodes. Workspace. 4974b9ad928SBarry Smith . vec1_D - Sequential vector of local interior nodes. Workspace. 4984b9ad928SBarry Smith . vec2_D - Sequential vector of local interior nodes. Workspace. 499*2fe279fdSBarry Smith - work_N - Array of all local nodes (interior and interface). Workspace. 5004b9ad928SBarry Smith 5014b9ad928SBarry Smith */ 502d71ae5a4SJacob Faibussowitsch PetscErrorCode PCNNBalancing(PC pc, Vec r, Vec u, Vec z, Vec vec1_B, Vec vec2_B, Vec vec3_B, Vec vec1_D, Vec vec2_D, PetscScalar *work_N) 503d71ae5a4SJacob Faibussowitsch { 50413f74950SBarry Smith PetscInt k; 5054b9ad928SBarry Smith PetscScalar value; 5064b9ad928SBarry Smith PetscScalar *lambda; 5074b9ad928SBarry Smith PC_NN *pcnn = (PC_NN *)(pc->data); 5084b9ad928SBarry Smith PC_IS *pcis = (PC_IS *)(pc->data); 5094b9ad928SBarry Smith 5104b9ad928SBarry Smith PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyCoarse, pc, 0, 0, 0)); 5124b9ad928SBarry Smith if (u) { 5132fa5cd67SKarl Rupp if (!vec3_B) vec3_B = u; 5149566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(vec1_B, pcis->D, u)); 5159566063dSJacob Faibussowitsch PetscCall(VecSet(z, 0.0)); 5169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 5179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 5189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5209566063dSJacob Faibussowitsch PetscCall(PCISApplySchur(pc, vec2_B, vec3_B, (Vec)0, vec1_D, vec2_D)); 5219566063dSJacob Faibussowitsch PetscCall(VecScale(vec3_B, -1.0)); 5229566063dSJacob Faibussowitsch PetscCall(VecCopy(r, z)); 5239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, vec3_B, z, ADD_VALUES, SCATTER_REVERSE)); 5249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, vec3_B, z, ADD_VALUES, SCATTER_REVERSE)); 5254b9ad928SBarry Smith } else { 5269566063dSJacob Faibussowitsch PetscCall(VecCopy(r, z)); 5274b9ad928SBarry Smith } 5289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5309566063dSJacob Faibussowitsch PetscCall(PCISScatterArrayNToVecB(work_N, vec2_B, INSERT_VALUES, SCATTER_REVERSE, pc)); 5312fa5cd67SKarl Rupp for (k = 0, value = 0.0; k < pcis->n_shared[0]; k++) value += pcnn->DZ_IN[0][k] * work_N[pcis->shared[0][k]]; 5324b9ad928SBarry Smith value *= pcnn->factor_coarse_rhs; /* This factor is set in CreateCoarseMatrix(). */ 5334b9ad928SBarry Smith { 53413f74950SBarry Smith PetscMPIInt rank; 5359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 5369566063dSJacob Faibussowitsch PetscCall(VecSetValue(pcnn->coarse_b, rank, value, INSERT_VALUES)); 5374b9ad928SBarry Smith /* 5384b9ad928SBarry Smith Since we are only inserting local values (one value actually) we don't need to do the 5394b9ad928SBarry Smith reduction that tells us there is no data that needs to be moved. Hence we comment out these 5409566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(pcnn->coarse_b)); 5419566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd (pcnn->coarse_b)); 5424b9ad928SBarry Smith */ 5434b9ad928SBarry Smith } 5449566063dSJacob Faibussowitsch PetscCall(KSPSolve(pcnn->ksp_coarse, pcnn->coarse_b, pcnn->coarse_x)); 5459566063dSJacob Faibussowitsch if (!u) PetscCall(VecScale(pcnn->coarse_x, -1.0)); 5469566063dSJacob Faibussowitsch PetscCall(VecGetArray(pcnn->coarse_x, &lambda)); 5472fa5cd67SKarl Rupp for (k = 0; k < pcis->n_shared[0]; k++) work_N[pcis->shared[0][k]] = *lambda * pcnn->DZ_IN[0][k]; 5489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pcnn->coarse_x, &lambda)); 5499566063dSJacob Faibussowitsch PetscCall(PCISScatterArrayNToVecB(work_N, vec2_B, INSERT_VALUES, SCATTER_FORWARD, pc)); 5509566063dSJacob Faibussowitsch PetscCall(VecSet(z, 0.0)); 5519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE)); 5529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE)); 5534b9ad928SBarry Smith if (!u) { 5549566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5559566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD)); 5569566063dSJacob Faibussowitsch PetscCall(PCISApplySchur(pc, vec2_B, vec1_B, (Vec)0, vec1_D, vec2_D)); 5579566063dSJacob Faibussowitsch PetscCall(VecCopy(r, z)); 5584b9ad928SBarry Smith } 5599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 5609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 5619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyCoarse, pc, 0, 0, 0)); 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5634b9ad928SBarry Smith } 5644b9ad928SBarry Smith 5654b9ad928SBarry Smith /* */ 5664b9ad928SBarry Smith /* From now on, "footnotes" (or "historical notes"). */ 5674b9ad928SBarry Smith /* */ 5684b9ad928SBarry Smith /* 569f1580f4eSBarry Smith Historical note 01 570f1580f4eSBarry Smith 5714b9ad928SBarry Smith We considered the possibility of an alternative D_i that would still 5724b9ad928SBarry Smith provide a partition of unity (i.e., $ \sum_i N_i D_i N_i^T = I $). 5734b9ad928SBarry Smith The basic principle was still the pseudo-inverse of the counting 5744b9ad928SBarry Smith function; the difference was that we would not count subdomains 5754b9ad928SBarry Smith that do not contribute to the coarse space (i.e., not pure-Neumann 5764b9ad928SBarry Smith subdomains). 5774b9ad928SBarry Smith 5784b9ad928SBarry Smith This turned out to be a bad idea: we would solve trivial Neumann 5794b9ad928SBarry Smith problems in the not pure-Neumann subdomains, since we would be scaling 5804b9ad928SBarry Smith the balanced residual by zero. 5814b9ad928SBarry Smith */ 5824b9ad928SBarry Smith 5834b9ad928SBarry Smith /* 584f1580f4eSBarry Smith Historical note 02 585f1580f4eSBarry Smith 5864b9ad928SBarry Smith We tried an alternative coarse problem, that would eliminate exactly a 5874b9ad928SBarry Smith constant error. Turned out not to improve the overall convergence. 5884b9ad928SBarry Smith */ 589