xref: /petsc/src/ksp/pc/impls/is/nn/nn.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
144b9ad928SBarry Smith    Notes:
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 */
186849ba73SBarry Smith static PetscErrorCode PCSetUp_NN(PC pc)
194b9ad928SBarry Smith {
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   }
274b9ad928SBarry Smith   PetscFunctionReturn(0);
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:
354b9ad928SBarry Smith .  pc - the preconditioner context
364b9ad928SBarry 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  */
436849ba73SBarry Smith static PetscErrorCode PCApply_NN(PC pc,Vec r,Vec z)
444b9ad928SBarry Smith {
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   */
72*d0609cedSBarry 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));
934b9ad928SBarry Smith   PetscFunctionReturn(0);
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 */
1066849ba73SBarry Smith static PetscErrorCode PCDestroy_NN(PC pc)
1074b9ad928SBarry Smith {
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));
1264b9ad928SBarry Smith   PetscFunctionReturn(0);
1274b9ad928SBarry Smith }
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
130dfc58137SBarry Smith /*MC
13137a17b4dSBarry Smith    PCNN - Balancing Neumann-Neumann for scalar elliptic PDEs.
1324b9ad928SBarry Smith 
13337a17b4dSBarry Smith    Options Database Keys:
13437a17b4dSBarry Smith +    -pc_nn_turn_off_first_balancing - do not balance the residual before solving the local Neumann problems
13537a17b4dSBarry Smith                                        (this skips the first coarse grid solve in the preconditioner)
1365c9740d6SBarry Smith .    -pc_nn_turn_off_second_balancing - do not balance the solution solving the local Neumann problems
13737a17b4dSBarry Smith                                        (this skips the second coarse grid solve in the preconditioner)
1385c9740d6SBarry Smith .    -pc_is_damp_fixed <fact> -
1395c9740d6SBarry Smith .    -pc_is_remove_nullspace_fixed -
1405c9740d6SBarry Smith .    -pc_is_set_damping_factor_floating <fact> -
1415c9740d6SBarry Smith .    -pc_is_not_damp_floating -
142a2b725a8SWilliam Gropp -    -pc_is_not_remove_nullspace_floating -
1434b9ad928SBarry Smith 
1448eb2139fSSatish Balay    Level: intermediate
14537a17b4dSBarry Smith 
14695452b02SPatrick Sanan    Notes:
14795452b02SPatrick Sanan     The matrix used with this preconditioner must be of type MATIS
14837a17b4dSBarry Smith 
1491b99c236SBarry Smith           Unlike more 'conventional' Neumann-Neumann preconditioners this iterates over ALL the
1501b99c236SBarry Smith           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
1511b99c236SBarry Smith           on the subdomains; though in our experience using approximate solvers is slower.).
1521b99c236SBarry Smith 
1535c9740d6SBarry Smith           Options for the coarse grid preconditioner can be set with -nn_coarse_pc_xxx
1545c9740d6SBarry Smith           Options for the Dirichlet subproblem preconditioner can be set with -is_localD_pc_xxx
1555c9740d6SBarry Smith           Options for the Neumann subproblem preconditioner can be set with -is_localN_pc_xxx
15637a17b4dSBarry Smith 
15737a17b4dSBarry Smith    Contributed by Paulo Goldfeld
15837a17b4dSBarry Smith 
15978392ef1SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
16037a17b4dSBarry Smith M*/
161b2573a8aSBarry Smith 
1628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_NN(PC pc)
1634b9ad928SBarry Smith {
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   */
1719566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&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;
1964b9ad928SBarry Smith   PetscFunctionReturn(0);
1974b9ad928SBarry Smith }
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2004b9ad928SBarry Smith /*
2014b9ad928SBarry Smith    PCNNCreateCoarseMatrix -
2024b9ad928SBarry Smith */
203dfbe8321SBarry Smith PetscErrorCode PCNNCreateCoarseMatrix(PC pc)
2044b9ad928SBarry Smith {
2054b9ad928SBarry Smith   MPI_Request    *send_request, *recv_request;
20613f74950SBarry Smith   PetscInt       i, j, k;
2074b9ad928SBarry Smith   PetscScalar    *mat;     /* Sub-matrix with this subdomain's contribution to the coarse matrix             */
2084b9ad928SBarry Smith   PetscScalar    **DZ_OUT; /* proc[k].DZ_OUT[i][] = bit of vector to be sent from processor k to processor i */
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith   /* aliasing some names */
2114b9ad928SBarry Smith   PC_IS       *pcis     = (PC_IS*)(pc->data);
2124b9ad928SBarry Smith   PC_NN       *pcnn     = (PC_NN*)pc->data;
21313f74950SBarry Smith   PetscInt    n_neigh   = pcis->n_neigh;
21413f74950SBarry Smith   PetscInt    *neigh    = pcis->neigh;
21513f74950SBarry Smith   PetscInt    *n_shared = pcis->n_shared;
21613f74950SBarry Smith   PetscInt    **shared  = pcis->shared;
2174b9ad928SBarry Smith   PetscScalar **DZ_IN;   /* Must be initialized after memory allocation. */
2184b9ad928SBarry Smith 
2194b9ad928SBarry Smith   PetscFunctionBegin;
2204b9ad928SBarry Smith   /* Allocate memory for mat (the +1 is to handle the case n_neigh equal to zero) */
2219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_neigh*n_neigh+1,&mat));
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith   /* Allocate memory for DZ */
2244b9ad928SBarry Smith   /* Notice that DZ_OUT[0] is allocated some space that is never used. */
2254b9ad928SBarry Smith   /* This is just in order to DZ_OUT and DZ_IN to have exactly the same form. */
2264b9ad928SBarry Smith   {
22713f74950SBarry Smith     PetscInt size_of_Z = 0;
2289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&pcnn->DZ_IN));
2294b9ad928SBarry Smith     DZ_IN = pcnn->DZ_IN;
2309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&DZ_OUT));
2312fa5cd67SKarl Rupp     for (i=0; i<n_neigh; i++) size_of_Z += n_shared[i];
2329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_IN[0]));
2339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_OUT[0]));
2344b9ad928SBarry Smith   }
2354b9ad928SBarry Smith   for (i=1; i<n_neigh; i++) {
2364b9ad928SBarry Smith     DZ_IN[i]  = DZ_IN [i-1] + n_shared[i-1];
2374b9ad928SBarry Smith     DZ_OUT[i] = DZ_OUT[i-1] + n_shared[i-1];
2384b9ad928SBarry Smith   }
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   /* Set the values of DZ_OUT, in order to send this info to the neighbours */
2414b9ad928SBarry Smith   /* First, set the auxiliary array pcis->work_N. */
2429566063dSJacob Faibussowitsch   PetscCall(PCISScatterArrayNToVecB(pcis->work_N,pcis->D,INSERT_VALUES,SCATTER_REVERSE,pc));
2434b9ad928SBarry Smith   for (i=1; i<n_neigh; i++) {
2444b9ad928SBarry Smith     for (j=0; j<n_shared[i]; j++) {
2454b9ad928SBarry Smith       DZ_OUT[i][j] = pcis->work_N[shared[i][j]];
2464b9ad928SBarry Smith     }
2474b9ad928SBarry Smith   }
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith   /* Non-blocking send/receive the common-interface chunks of scaled nullspaces */
2504b9ad928SBarry Smith   /* Notice that send_request[] and recv_request[] could have one less element. */
2514b9ad928SBarry Smith   /* We make them longer to have request[i] corresponding to neigh[i].          */
2524b9ad928SBarry Smith   {
25313f74950SBarry Smith     PetscMPIInt tag;
2549566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)pc,&tag));
2559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n_neigh+1,&send_request,n_neigh+1,&recv_request));
2564b9ad928SBarry Smith     for (i=1; i<n_neigh; i++) {
2579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend((void*)(DZ_OUT[i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,PetscObjectComm((PetscObject)pc),&(send_request[i])));
2589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv((void*)(DZ_IN [i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,PetscObjectComm((PetscObject)pc),&(recv_request[i])));
2594b9ad928SBarry Smith     }
2604b9ad928SBarry Smith   }
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith   /* Set DZ_IN[0][] (recall that neigh[0]==rank, always) */
2632fa5cd67SKarl Rupp   for (j=0; j<n_shared[0]; j++) DZ_IN[0][j] = pcis->work_N[shared[0][j]];
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith   /* Start computing with local D*Z while communication goes on.    */
2664b9ad928SBarry Smith   /* Apply Schur complement. The result is "stored" in vec (more    */
2674b9ad928SBarry Smith   /* precisely, vec points to the result, stored in pc_nn->vec1_B)  */
2684b9ad928SBarry Smith   /* and also scattered to pcnn->work_N.                            */
2699566063dSJacob 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));
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith   /* Compute the first column, while completing the receiving. */
2724b9ad928SBarry Smith   for (i=0; i<n_neigh; i++) {
2734b9ad928SBarry Smith     MPI_Status  stat;
27413f74950SBarry Smith     PetscMPIInt ind=0;
2759566063dSJacob Faibussowitsch     if (i>0) {PetscCallMPI(MPI_Waitany(n_neigh-1,recv_request+1,&ind,&stat)); ind++;}
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. */
3794b9ad928SBarry Smith   PetscFunctionReturn(0);
3804b9ad928SBarry Smith }
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
3834b9ad928SBarry Smith /*
3844b9ad928SBarry Smith    PCNNApplySchurToChunk -
3854b9ad928SBarry Smith 
3864b9ad928SBarry Smith    Input parameters:
3874b9ad928SBarry Smith .  pcnn
3884b9ad928SBarry Smith .  n - size of chunk
3894b9ad928SBarry Smith .  idx - indices of chunk
3904b9ad928SBarry Smith .  chunk - values
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Output parameters:
3934b9ad928SBarry Smith .  array_N - result of Schur complement applied to chunk, scattered to big array
3944b9ad928SBarry Smith .  vec1_B  - result of Schur complement applied to chunk
3954b9ad928SBarry Smith .  vec2_B  - garbage (used as work space)
3964b9ad928SBarry Smith .  vec1_D  - garbage (used as work space)
3974b9ad928SBarry Smith .  vec2_D  - garbage (used as work space)
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith */
40013f74950SBarry Smith 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)
4014b9ad928SBarry Smith {
40213f74950SBarry Smith   PetscInt       i;
4034b9ad928SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
4044b9ad928SBarry Smith 
4054b9ad928SBarry Smith   PetscFunctionBegin;
4069566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(array_N, pcis->n));
4072fa5cd67SKarl Rupp   for (i=0; i<n; i++) array_N[idx[i]] = chunk[i];
4089566063dSJacob Faibussowitsch   PetscCall(PCISScatterArrayNToVecB(array_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc));
4099566063dSJacob Faibussowitsch   PetscCall(PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D));
4109566063dSJacob Faibussowitsch   PetscCall(PCISScatterArrayNToVecB(array_N,vec1_B,INSERT_VALUES,SCATTER_REVERSE,pc));
4114b9ad928SBarry Smith   PetscFunctionReturn(0);
4124b9ad928SBarry Smith }
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
4154b9ad928SBarry Smith /*
4164b9ad928SBarry Smith    PCNNApplyInterfacePreconditioner - Apply the interface preconditioner, i.e.,
4174b9ad928SBarry Smith                                       the preconditioner for the Schur complement.
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith    Input parameter:
4204b9ad928SBarry Smith .  r - global vector of interior and interface nodes. The values on the interior nodes are NOT used.
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith    Output parameters:
4234b9ad928SBarry Smith .  z - global vector of interior and interface nodes. The values on the interface are the result of
4244b9ad928SBarry Smith        the application of the interface preconditioner to the interface part of r. The values on the
4254b9ad928SBarry Smith        interior nodes are garbage.
4264b9ad928SBarry Smith .  work_N - array of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
4274b9ad928SBarry Smith .  vec1_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
4284b9ad928SBarry Smith .  vec2_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
4294b9ad928SBarry Smith .  vec3_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
4304b9ad928SBarry Smith .  vec1_D - vector of local interior nodes; returns garbage (used as work space)
4314b9ad928SBarry Smith .  vec2_D - vector of local interior nodes; returns garbage (used as work space)
4324b9ad928SBarry Smith .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
4334b9ad928SBarry Smith .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith */
4362fa5cd67SKarl Rupp 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)
4374b9ad928SBarry Smith {
4384b9ad928SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
4394b9ad928SBarry Smith 
4404b9ad928SBarry Smith   PetscFunctionBegin;
4414b9ad928SBarry Smith   /*
4424b9ad928SBarry Smith     First balancing step.
4434b9ad928SBarry Smith   */
4444b9ad928SBarry Smith   {
445ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
4469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc_nn_turn_off_first_balancing",&flg,NULL));
4474b9ad928SBarry Smith     if (!flg) {
4489566063dSJacob Faibussowitsch       PetscCall(PCNNBalancing(pc,r,(Vec)0,z,vec1_B,vec2_B,(Vec)0,vec1_D,vec2_D,work_N));
4494b9ad928SBarry Smith     } else {
4509566063dSJacob Faibussowitsch       PetscCall(VecCopy(r,z));
4514b9ad928SBarry Smith     }
4524b9ad928SBarry Smith   }
4534b9ad928SBarry Smith 
4544b9ad928SBarry Smith   /*
4554b9ad928SBarry Smith     Extract the local interface part of z and scale it by D
4564b9ad928SBarry Smith   */
4579566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,z,vec1_B,INSERT_VALUES,SCATTER_FORWARD));
4589566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->global_to_B,z,vec1_B,INSERT_VALUES,SCATTER_FORWARD));
4599566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(vec2_B,pcis->D,vec1_B));
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith   /* Neumann Solver */
4629566063dSJacob Faibussowitsch   PetscCall(PCISApplyInvSchur(pc,vec2_B,vec1_B,vec1_N,vec2_N));
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith   /*
4654b9ad928SBarry Smith     Second balancing step.
4664b9ad928SBarry Smith   */
4674b9ad928SBarry Smith   {
468ace3abfcSBarry Smith     PetscBool flg = PETSC_FALSE;
4699566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc_turn_off_second_balancing",&flg,NULL));
4704b9ad928SBarry Smith     if (!flg) {
4719566063dSJacob Faibussowitsch       PetscCall(PCNNBalancing(pc,r,vec1_B,z,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D,work_N));
4724b9ad928SBarry Smith     } else {
4739566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(vec2_B,pcis->D,vec1_B));
4749566063dSJacob Faibussowitsch       PetscCall(VecSet(z,0.0));
4759566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE));
4769566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd  (pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE));
4774b9ad928SBarry Smith     }
4784b9ad928SBarry Smith   }
4794b9ad928SBarry Smith   PetscFunctionReturn(0);
4804b9ad928SBarry Smith }
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
4834b9ad928SBarry Smith /*
4844b9ad928SBarry Smith    PCNNBalancing - Computes z, as given in equations (15) and (16) (if the
4854b9ad928SBarry Smith                    input argument u is provided), or s, as given in equations
4864b9ad928SBarry Smith                    (12) and (13), if the input argument u is a null vector.
4874b9ad928SBarry Smith                    Notice that the input argument u plays the role of u_i in
4884b9ad928SBarry Smith                    equation (14). The equation numbers refer to [Man93].
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith    Input Parameters:
4914b9ad928SBarry Smith .  pcnn - NN preconditioner context.
4924b9ad928SBarry Smith .  r - MPI vector of all nodes (interior and interface). It's preserved.
4934b9ad928SBarry Smith .  u - (Optional) sequential vector of local interface nodes. It's preserved UNLESS vec3_B is null.
4944b9ad928SBarry Smith 
4954b9ad928SBarry Smith    Output Parameters:
4964b9ad928SBarry Smith .  z - MPI vector of interior and interface nodes. Returns s or z (see description above).
4974b9ad928SBarry Smith .  vec1_B - Sequential vector of local interface nodes. Workspace.
4984b9ad928SBarry Smith .  vec2_B - Sequential vector of local interface nodes. Workspace.
4994b9ad928SBarry Smith .  vec3_B - (Optional) sequential vector of local interface nodes. Workspace.
5004b9ad928SBarry Smith .  vec1_D - Sequential vector of local interior nodes. Workspace.
5014b9ad928SBarry Smith .  vec2_D - Sequential vector of local interior nodes. Workspace.
5024b9ad928SBarry Smith .  work_N - Array of all local nodes (interior and interface). Workspace.
5034b9ad928SBarry Smith 
5044b9ad928SBarry Smith */
5052fa5cd67SKarl Rupp 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)
5064b9ad928SBarry Smith {
50713f74950SBarry Smith   PetscInt       k;
5084b9ad928SBarry Smith   PetscScalar    value;
5094b9ad928SBarry Smith   PetscScalar    *lambda;
5104b9ad928SBarry Smith   PC_NN          *pcnn = (PC_NN*)(pc->data);
5114b9ad928SBarry Smith   PC_IS          *pcis = (PC_IS*)(pc->data);
5124b9ad928SBarry Smith 
5134b9ad928SBarry Smith   PetscFunctionBegin;
5149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyCoarse,pc,0,0,0));
5154b9ad928SBarry Smith   if (u) {
5162fa5cd67SKarl Rupp     if (!vec3_B) vec3_B = u;
5179566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(vec1_B,pcis->D,u));
5189566063dSJacob Faibussowitsch     PetscCall(VecSet(z,0.0));
5199566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE));
5209566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd  (pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE));
5219566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5229566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5239566063dSJacob Faibussowitsch     PetscCall(PCISApplySchur(pc,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D));
5249566063dSJacob Faibussowitsch     PetscCall(VecScale(vec3_B,-1.0));
5259566063dSJacob Faibussowitsch     PetscCall(VecCopy(r,z));
5269566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,vec3_B,z,ADD_VALUES,SCATTER_REVERSE));
5279566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd  (pcis->global_to_B,vec3_B,z,ADD_VALUES,SCATTER_REVERSE));
5284b9ad928SBarry Smith   } else {
5299566063dSJacob Faibussowitsch     PetscCall(VecCopy(r,z));
5304b9ad928SBarry Smith   }
5319566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5329566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5339566063dSJacob Faibussowitsch   PetscCall(PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_REVERSE,pc));
5342fa5cd67SKarl 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]];
5354b9ad928SBarry Smith   value *= pcnn->factor_coarse_rhs;  /* This factor is set in CreateCoarseMatrix(). */
5364b9ad928SBarry Smith   {
53713f74950SBarry Smith     PetscMPIInt rank;
5389566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
5399566063dSJacob Faibussowitsch     PetscCall(VecSetValue(pcnn->coarse_b,rank,value,INSERT_VALUES));
5404b9ad928SBarry Smith     /*
5414b9ad928SBarry Smith        Since we are only inserting local values (one value actually) we don't need to do the
5424b9ad928SBarry Smith        reduction that tells us there is no data that needs to be moved. Hence we comment out these
5439566063dSJacob Faibussowitsch        PetscCall(VecAssemblyBegin(pcnn->coarse_b));
5449566063dSJacob Faibussowitsch        PetscCall(VecAssemblyEnd  (pcnn->coarse_b));
5454b9ad928SBarry Smith     */
5464b9ad928SBarry Smith   }
5479566063dSJacob Faibussowitsch   PetscCall(KSPSolve(pcnn->ksp_coarse,pcnn->coarse_b,pcnn->coarse_x));
5489566063dSJacob Faibussowitsch   if (!u) PetscCall(VecScale(pcnn->coarse_x,-1.0));
5499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pcnn->coarse_x,&lambda));
5502fa5cd67SKarl Rupp   for (k=0; k<pcis->n_shared[0]; k++) work_N[pcis->shared[0][k]] = *lambda * pcnn->DZ_IN[0][k];
5519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pcnn->coarse_x,&lambda));
5529566063dSJacob Faibussowitsch   PetscCall(PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc));
5539566063dSJacob Faibussowitsch   PetscCall(VecSet(z,0.0));
5549566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE));
5559566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE));
5564b9ad928SBarry Smith   if (!u) {
5579566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5589566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD));
5599566063dSJacob Faibussowitsch     PetscCall(PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D));
5609566063dSJacob Faibussowitsch     PetscCall(VecCopy(r,z));
5614b9ad928SBarry Smith   }
5629566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE));
5639566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd  (pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE));
5649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyCoarse,pc,0,0,0));
5654b9ad928SBarry Smith   PetscFunctionReturn(0);
5664b9ad928SBarry Smith }
5674b9ad928SBarry Smith 
5684b9ad928SBarry Smith /*  -------   E N D   O F   T H E   C O D E   -------  */
5694b9ad928SBarry Smith /*                                                     */
5704b9ad928SBarry Smith /*  From now on, "footnotes" (or "historical notes").  */
5714b9ad928SBarry Smith /*                                                     */
5724b9ad928SBarry Smith /*  -------------------------------------------------  */
5734b9ad928SBarry Smith 
5744b9ad928SBarry Smith /* --------------------------------------------------------------------------
5754b9ad928SBarry Smith    Historical note 01
5764b9ad928SBarry Smith    -------------------------------------------------------------------------- */
5774b9ad928SBarry Smith /*
5784b9ad928SBarry Smith    We considered the possibility of an alternative D_i that would still
5794b9ad928SBarry Smith    provide a partition of unity (i.e., $ \sum_i  N_i D_i N_i^T = I $).
5804b9ad928SBarry Smith    The basic principle was still the pseudo-inverse of the counting
5814b9ad928SBarry Smith    function; the difference was that we would not count subdomains
5824b9ad928SBarry Smith    that do not contribute to the coarse space (i.e., not pure-Neumann
5834b9ad928SBarry Smith    subdomains).
5844b9ad928SBarry Smith 
5854b9ad928SBarry Smith    This turned out to be a bad idea:  we would solve trivial Neumann
5864b9ad928SBarry Smith    problems in the not pure-Neumann subdomains, since we would be scaling
5874b9ad928SBarry Smith    the balanced residual by zero.
5884b9ad928SBarry Smith */
5894b9ad928SBarry Smith 
5904b9ad928SBarry Smith /* --------------------------------------------------------------------------
5914b9ad928SBarry Smith    Historical note 02
5924b9ad928SBarry Smith    -------------------------------------------------------------------------- */
5934b9ad928SBarry Smith /*
5944b9ad928SBarry Smith    We tried an alternative coarse problem, that would eliminate exactly a
5954b9ad928SBarry Smith    constant error. Turned out not to improve the overall convergence.
5964b9ad928SBarry Smith */
597