PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::lagrangianSchurApply(Vec X, Vec Y)
{
    const double half_time = d_integrator_time + 0.5 * d_current_dt;

    // The Schur complement is: I-dt*J*inv(L)*S*A/4
    d_ib_implicit_ops->computeLinearizedLagrangianForce(X, half_time);
    d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0);
    d_u_phys_bdry_op->setPatchDataIndex(d_f_idx);
    d_ib_implicit_ops->spreadLinearizedForce(d_f_idx,
                                             d_u_phys_bdry_op,
                                             getProlongRefineSchedules(d_object_name + "::f"),
                                             half_time);
    d_u_scratch_vec->setToScalar(0.0);
    d_hier_velocity_data_ops->copyData(d_f_scratch_vec->getComponentDescriptorIndex(0),
                                       d_f_idx);
#if 0
    d_stokes_solver->setHomogeneousBc(true);
    d_stokes_solver->solveSystem(*d_u_scratch_vec, *d_f_scratch_vec);
#else
    d_u_scratch_vec->scale(
        1.0 / d_stokes_op->getVelocityPoissonSpecifications().getCConstant(), d_f_scratch_vec);
#endif
    d_hier_velocity_data_ops->scale(
        d_u_idx, 0.5, d_u_scratch_vec->getComponentDescriptorIndex(0));
    d_ib_implicit_ops->interpolateLinearizedVelocity(
        d_u_idx,
        getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
        getGhostfillRefineSchedules(d_object_name + "::u"),
        half_time);
    d_ib_implicit_ops->computeLinearizedResidual(X, Y);
    return 0;
} // lagrangianSchurApply
PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::compositeIBJacobianApply(Vec x, Vec f)
{
    PetscErrorCode ierr;
    const double half_time = d_integrator_time + 0.5 * d_current_dt;

    Vec* component_sol_vecs;
    Vec* component_rhs_vecs;
    ierr = VecMultiVecGetSubVecs(x, &component_sol_vecs);
    IBTK_CHKERRQ(ierr);
    ierr = VecMultiVecGetSubVecs(f, &component_rhs_vecs);
    IBTK_CHKERRQ(ierr);

    Pointer<SAMRAIVectorReal<NDIM, double> > u =
        PETScSAMRAIVectorReal::getSAMRAIVector(component_sol_vecs[0]);
    Pointer<SAMRAIVectorReal<NDIM, double> > f_u =
        PETScSAMRAIVectorReal::getSAMRAIVector(component_rhs_vecs[0]);

    Pointer<Variable<NDIM> > u_var = d_ins_hier_integrator->getVelocityVariable();
    const int u_idx = u->getComponentDescriptorIndex(0);
    const int f_u_idx = f_u->getComponentDescriptorIndex(0);

    Vec X = component_sol_vecs[1];
    Vec R = component_rhs_vecs[1];

    // Evaluate the Eulerian terms.
    d_stokes_op->setHomogeneousBc(true);
    d_stokes_op->apply(*u, *f_u);

    d_ib_implicit_ops->computeLinearizedLagrangianForce(X, half_time);
    if (d_enable_logging)
        plog << d_object_name
             << "::integrateHierarchy(): spreading Lagrangian force to the Eulerian grid\n";
    d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0);
    d_u_phys_bdry_op->setPatchDataIndex(d_f_idx);
    d_ib_implicit_ops->spreadLinearizedForce(d_f_idx,
                                             d_u_phys_bdry_op,
                                             getProlongRefineSchedules(d_object_name + "::f"),
                                             half_time);
    d_hier_velocity_data_ops->subtract(f_u_idx, f_u_idx, d_f_idx);
    ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(component_rhs_vecs[0]));
    IBTK_CHKERRQ(ierr);

    // Evaluate the Lagrangian terms.
    d_hier_velocity_data_ops->scale(d_u_idx, 0.5, u_idx);
    d_u_phys_bdry_op->setPatchDataIndex(d_u_idx);
    d_ib_implicit_ops->interpolateLinearizedVelocity(
        d_u_idx,
        getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
        getGhostfillRefineSchedules(d_object_name + "::u"),
        half_time);
    d_ib_implicit_ops->computeLinearizedResidual(X, R);

    // Ensure that PETSc sees that the state of the RHS vector has changed.
    // This is a nasty hack.
    ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(f));
    IBTK_CHKERRQ(ierr);
    return ierr;
} // compositeIBJacobianApply
void
IBHierarchyIntegrator::initializePatchHierarchy(Pointer<PatchHierarchy<NDIM> > hierarchy,
                                                Pointer<GriddingAlgorithm<NDIM> > gridding_alg)
{
    if (d_hierarchy_is_initialized) return;

    // Initialize Eulerian data.
    HierarchyIntegrator::initializePatchHierarchy(hierarchy, gridding_alg);

    // Begin Lagrangian data movement.
    d_ib_method_ops->beginDataRedistribution(hierarchy, gridding_alg);

    // Finish Lagrangian data movement.
    d_ib_method_ops->endDataRedistribution(hierarchy, gridding_alg);

    // Initialize Lagrangian data on the patch hierarchy.
    const int coarsest_ln = 0;
    const int finest_ln = hierarchy->getFinestLevelNumber();
    for (int ln = coarsest_ln; ln <= finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        level->allocatePatchData(d_u_idx, d_integrator_time);
        level->allocatePatchData(d_scratch_data, d_integrator_time);
    }
    VariableDatabase<NDIM>* var_db = VariableDatabase<NDIM>::getDatabase();
    const int u_current_idx =
        var_db->mapVariableAndContextToIndex(d_u_var, getCurrentContext());
    d_hier_velocity_data_ops->copyData(d_u_idx, u_current_idx);
    const bool initial_time = MathUtilities<double>::equalEps(d_integrator_time, d_start_time);
    d_u_phys_bdry_op->setPatchDataIndex(d_u_idx);
    d_ib_method_ops->initializePatchHierarchy(
        hierarchy,
        gridding_alg,
        d_u_idx,
        getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
        getGhostfillRefineSchedules(d_object_name + "::u"),
        d_integrator_step,
        d_integrator_time,
        initial_time);
    for (int ln = coarsest_ln; ln <= finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        level->deallocatePatchData(d_u_idx);
        level->deallocatePatchData(d_scratch_data);
    }

    // Indicate that the hierarchy is initialized.
    d_hierarchy_is_initialized = true;
    return;
} // initializePatchHierarchy
void
IBSimpleHierarchyIntegrator::integrateHierarchy(const double current_time, const double new_time, const int cycle_num)
{
    IBHierarchyIntegrator::integrateHierarchy(current_time, new_time, cycle_num);

    const int finest_level_num = d_hierarchy->getFinestLevelNumber();
    PetscErrorCode ierr;
    const double dt = new_time - current_time;
    Pointer<IBMethod> p_ib_method_ops = d_ib_method_ops;
    LDataManager* l_data_manager = p_ib_method_ops->getLDataManager();

    // Here we implement a simple time integration scheme:
    //
    // (1) Compute Lagrangian forces and spread those forces to the grid: f =
    //     S[X^{n}] F^{n}.
    //
    // (2) Solve the fluid equations using the fluid solver registered with this
    //     class using the forcing computed in (1).
    //
    // (3) Interpolate the Eulerian velocity and update the positions of the
    //     Lagrangian points: X^{n+1} = X^{n} + dt*S^{*}[X^{n}] u^{n+1}.
    //
    // NOTE: We assume here that all IB data are assigned to the finest level of
    // the AMR patch hierarchy.

    ///////////////////////////////////////////////////////////////////////////
    // (1) Compute Lagrangian forces and spread those forces to the grid: f =
    //     S[X^{n}] F^{n}..

    // For simplicity, set the Lagrangian force density to equal zero.
    ierr = VecZeroEntries(d_F_data->getVec());
    IBTK_CHKERRQ(ierr);

    // Spread the forces to the grid.  We use the "current" Lagrangian position
    // data to define the locations from where the forces are spread.
    d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0);
    l_data_manager->spread(d_f_idx,
                           d_F_data,
                           d_X_current_data,
                           d_u_phys_bdry_op,
                           finest_level_num,
                           getProlongRefineSchedules(d_object_name + "::f"),
                           /*F_needs_ghost_fill*/ true,
                           /*X_needs_ghost_fill*/ true);

    // NOTE: Any additional Eulerian forcing should be computed here and added
    // to the data associated with d_f_idx.

    ///////////////////////////////////////////////////////////////////////////
    // (2) Solve the fluid equations using the fluid solver registered with this
    // class using the forcing computed in (1).
    const int ins_num_cycles = d_ins_hier_integrator->getNumberOfCycles();
    for (int ins_cycle_num = 0; ins_cycle_num < ins_num_cycles; ++ins_cycle_num)
    {
        d_ins_hier_integrator->integrateHierarchy(current_time, new_time, ins_cycle_num);
    }

    ///////////////////////////////////////////////////////////////////////////
    // (3) Interpolate the Eulerian velocity and update the positions of the
    // Lagrangian points: X^{n+1} = X^{n} + dt*S^{*}[X^{n}] u^{n+1}.
    //
    // NOTE: We use the "new" velocity data (defined at time level n+1) to
    // determine the velocity of the Lagrangian nodes.  We could also use the
    // "current" data (defined at time level n) or some other velocity field
    // here.  We use the "current" Lagrangian position data to define the
    // locations to where the velocities are interpolated.
    VariableDatabase<NDIM>* var_db = VariableDatabase<NDIM>::getDatabase();
    const int u_new_idx = var_db->mapVariableAndContextToIndex(d_ins_hier_integrator->getVelocityVariable(),
                                                               d_ins_hier_integrator->getNewContext());
    d_hier_velocity_data_ops->copyData(d_u_idx, u_new_idx);
    l_data_manager->interp(d_u_idx,
                           d_U_data,
                           d_X_current_data,
                           finest_level_num,
                           getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
                           getGhostfillRefineSchedules(d_object_name + "::u"),
                           current_time);
    ierr = VecWAXPY(d_X_new_data->getVec(), dt, d_U_data->getVec(), d_X_current_data->getVec());
    IBTK_CHKERRQ(ierr);
    return;
} // integrateHierarchy
PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::compositeIBPCApply(Vec x, Vec y)
{
    PetscErrorCode ierr;
    const double half_time = d_integrator_time + 0.5 * d_current_dt;

    Vec* component_x_vecs;
    Vec* component_y_vecs;
    ierr = VecMultiVecGetSubVecs(x, &component_x_vecs);
    IBTK_CHKERRQ(ierr);
    ierr = VecMultiVecGetSubVecs(y, &component_y_vecs);
    IBTK_CHKERRQ(ierr);

    Pointer<SAMRAIVectorReal<NDIM, double> > eul_x =
        PETScSAMRAIVectorReal::getSAMRAIVector(component_x_vecs[0]);
    Pointer<SAMRAIVectorReal<NDIM, double> > eul_y =
        PETScSAMRAIVectorReal::getSAMRAIVector(component_y_vecs[0]);

    Vec lag_x = component_x_vecs[1];
    Vec lag_y = component_y_vecs[1];

    // The full (nonlinear) system is:
    //
    //   L u(n+1) = S*F[X(n+1/2)] + f
    //   X(n+1) - X(n) = dt*U(n+1/2)
    //
    // where:
    //
    //   L = Eulerian operator (i.e. Stokes)
    //   F = Lagrangian force operator (potentially nonlinear)
    //   S = spreading operator
    //   J = interpolation operator = S^*
    //   X(n+1/2) = (X(n+1)+X(n))/2
    //   f = explicit right-hand side term
    //
    // For simplicity, only "lagged" S and J are considered, i.e., S and J are
    // not functions of the unknown X(n+1/2), but rather of some lagged
    // approximation to X(n+1/2).  This does not affect the stability of the
    // time stepping scheme, and the lagged values can be chosen so that the
    // overall scheme is second-order accurate.
    //
    // The linearized system is:
    //
    //   [L         -S*A/2] [u]
    //   [-dt*J/2   I     ] [X]
    //
    // where:
    //
    //   L = Eulerian operator
    //   A = dF/dX = Lagrangian operator
    //   S = spreading operator
    //   J = interpolation operator = S^*
    //
    // The Lagrangian Schur complement preconditioner is P = (4)*(3)*(2)*(1),
    // which is the inverse of the linearized system.
    //
    // (1) = [inv(L)  0]  ==>  [I         -inv(L)*S*A/2]
    //       [0       I]       [-dt*J/2   I            ]
    //
    // (2) = [I        0]  ==>  [I   -inv(L)*S*A/2      ]
    //       [dt*J/2   I]       [0   I-dt*J*inv(L)*S*A/4]
    //
    // Sc = Schur complement = I-dt*J*inv(L)*S*A/4
    //
    // (3) = [I   0      ]  ==>  [I   -inv(L)*S*A/2]
    //       [0   inv(Sc)]       [0   I            ]
    //
    // (4) = [I   inv(L)*S*A/2]  ==>  [I   0]
    //       [0   I           ]       [0   I]

    // Step 1: eul_y := inv(L)*eul_x
    eul_y->setToScalar(0.0);
    d_stokes_solver->setHomogeneousBc(true);
    d_stokes_solver->solveSystem(*eul_y, *eul_x);

    // Step 2: lag_y := lag_x + dt*J*eul_y/2
    d_hier_velocity_data_ops->scale(d_u_idx, -0.5, eul_y->getComponentDescriptorIndex(0));
    d_ib_implicit_ops->interpolateLinearizedVelocity(
        d_u_idx,
        getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
        getGhostfillRefineSchedules(d_object_name + "::u"),
        half_time);
    d_ib_implicit_ops->computeLinearizedResidual(lag_x, lag_y);

    // Step 3: lag_y := inv(Sc)*lag_y
    ierr = KSPSolve(d_schur_solver, lag_y, lag_y);
    IBTK_CHKERRQ(ierr);

    // Step 4: eul_y := eul_y + inv(L)*S*A*lag_y/2
    d_ib_implicit_ops->computeLinearizedLagrangianForce(lag_y, half_time);
    d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0);
    d_u_phys_bdry_op->setPatchDataIndex(d_f_idx);
    d_ib_implicit_ops->spreadLinearizedForce(d_f_idx,
                                             d_u_phys_bdry_op,
                                             getProlongRefineSchedules(d_object_name + "::f"),
                                             half_time);
    d_u_scratch_vec->setToScalar(0.0);
    d_f_scratch_vec->setToScalar(0.0);
    d_hier_velocity_data_ops->copyData(d_f_scratch_vec->getComponentDescriptorIndex(0),
                                       d_f_idx);
    d_stokes_solver->setHomogeneousBc(true);
    d_stokes_solver->solveSystem(*d_u_scratch_vec, *d_f_scratch_vec);
    eul_y->add(eul_y, d_u_scratch_vec);

    // Ensure that PETSc sees that the state of the RHS vector has changed.
    // This is a nasty hack.
    ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(component_y_vecs[0]));
    IBTK_CHKERRQ(ierr);
    ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(y));
    IBTK_CHKERRQ(ierr);

    return ierr;
} // compositeIBPCApply
void IBImplicitStaggeredHierarchyIntegrator::postprocessIntegrateHierarchy(
    const double current_time,
    const double new_time,
    const bool skip_synchronize_new_state_data,
    const int num_cycles)
{
    IBHierarchyIntegrator::postprocessIntegrateHierarchy(
        current_time, new_time, skip_synchronize_new_state_data, num_cycles);

    const int coarsest_ln = 0;
    const int finest_ln = d_hierarchy->getFinestLevelNumber();
    const double dt = new_time - current_time;
    VariableDatabase<NDIM>* var_db = VariableDatabase<NDIM>::getDatabase();
    const int u_new_idx = var_db->mapVariableAndContextToIndex(
        d_ins_hier_integrator->getVelocityVariable(), d_ins_hier_integrator->getNewContext());

    // Interpolate the Eulerian velocity to the curvilinear mesh.
    d_hier_velocity_data_ops->copyData(d_u_idx, u_new_idx);
    if (d_enable_logging)
        plog << d_object_name << "::postprocessIntegrateHierarchy(): interpolating Eulerian "
                                 "velocity to the Lagrangian mesh\n";
    d_ib_implicit_ops->interpolateVelocity(
        d_u_idx,
        getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"),
        getGhostfillRefineSchedules(d_object_name + "::u"),
        new_time);

    // Synchronize new state data.
    if (!skip_synchronize_new_state_data)
    {
        if (d_enable_logging)
            plog << d_object_name
                 << "::postprocessIntegrateHierarchy(): synchronizing updated data\n";
        synchronizeHierarchyData(NEW_DATA);
    }

    // Determine the CFL number.
    double cfl_max = 0.0;
    PatchCellDataOpsReal<NDIM, double> patch_cc_ops;
    PatchSideDataOpsReal<NDIM, double> patch_sc_ops;
    for (int ln = coarsest_ln; ln <= finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        for (PatchLevel<NDIM>::Iterator p(level); p; p++)
        {
            Pointer<Patch<NDIM> > patch = level->getPatch(p());
            const Box<NDIM>& patch_box = patch->getBox();
            const Pointer<CartesianPatchGeometry<NDIM> > pgeom = patch->getPatchGeometry();
            const double* const dx = pgeom->getDx();
            const double dx_min = *(std::min_element(dx, dx + NDIM));
            Pointer<CellData<NDIM, double> > u_cc_new_data = patch->getPatchData(u_new_idx);
            Pointer<SideData<NDIM, double> > u_sc_new_data = patch->getPatchData(u_new_idx);
            double u_max = 0.0;
            if (u_cc_new_data) u_max = patch_cc_ops.maxNorm(u_cc_new_data, patch_box);
            if (u_sc_new_data) u_max = patch_sc_ops.maxNorm(u_sc_new_data, patch_box);
            cfl_max = std::max(cfl_max, u_max * dt / dx_min);
        }
    }
    cfl_max = SAMRAI_MPI::maxReduction(cfl_max);
    d_regrid_cfl_estimate += cfl_max;
    if (d_enable_logging)
        plog << d_object_name << "::postprocessIntegrateHierarchy(): CFL number = " << cfl_max
             << "\n";
    if (d_enable_logging)
        plog << d_object_name
             << "::postprocessIntegrateHierarchy(): estimated upper bound on IB "
                "point displacement since last regrid = " << d_regrid_cfl_estimate << "\n";

    // Deallocate the fluid solver.
    const int ins_num_cycles = d_ins_hier_integrator->getNumberOfCycles();
    d_ins_hier_integrator->postprocessIntegrateHierarchy(
        current_time, new_time, skip_synchronize_new_state_data, ins_num_cycles);

    // Deallocate IB data.
    d_ib_implicit_ops->postprocessIntegrateData(current_time, new_time, num_cycles);

    // Deallocate Eulerian scratch data.
    for (int ln = coarsest_ln; ln <= finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        level->deallocatePatchData(d_u_idx);
        level->deallocatePatchData(d_f_idx);
        level->deallocatePatchData(d_scratch_data);
        level->deallocatePatchData(d_new_data);
    }

    // Execute any registered callbacks.
    executePostprocessIntegrateHierarchyCallbackFcns(
        current_time, new_time, skip_synchronize_new_state_data, num_cycles);
    return;
} // postprocessIntegrateHierarchy
void IBImplicitStaggeredHierarchyIntegrator::integrateHierarchy(const double current_time,
                                                                const double new_time,
                                                                const int cycle_num)
{
    IBHierarchyIntegrator::integrateHierarchy(current_time, new_time, cycle_num);

    Pointer<INSStaggeredHierarchyIntegrator> ins_hier_integrator = d_ins_hier_integrator;
    TBOX_ASSERT(ins_hier_integrator);

    PetscErrorCode ierr;
    int n_local;

    const int coarsest_ln = 0;
    const int finest_ln = d_hierarchy->getFinestLevelNumber();
    //  const double half_time = current_time+0.5*(new_time-current_time);

    VariableDatabase<NDIM>* var_db = VariableDatabase<NDIM>::getDatabase();
    Pointer<VariableContext> current_ctx = ins_hier_integrator->getCurrentContext();
    Pointer<VariableContext> scratch_ctx = ins_hier_integrator->getScratchContext();
    Pointer<VariableContext> new_ctx = ins_hier_integrator->getNewContext();

    const int wgt_cc_idx = d_hier_math_ops->getCellWeightPatchDescriptorIndex();
    const int wgt_sc_idx = d_hier_math_ops->getSideWeightPatchDescriptorIndex();

    Pointer<Variable<NDIM> > u_var = ins_hier_integrator->getVelocityVariable();
    //  const int u_current_idx = var_db->mapVariableAndContextToIndex(u_var, current_ctx);
    const int u_scratch_idx = var_db->mapVariableAndContextToIndex(u_var, scratch_ctx);
    //  const int u_new_idx     = var_db->mapVariableAndContextToIndex(u_var, new_ctx);

    Pointer<Variable<NDIM> > p_var = ins_hier_integrator->getPressureVariable();
    //  const int p_current_idx = var_db->mapVariableAndContextToIndex(p_var, current_ctx);
    const int p_scratch_idx = var_db->mapVariableAndContextToIndex(p_var, scratch_ctx);
    //  const int p_new_idx     = var_db->mapVariableAndContextToIndex(p_var, new_ctx);

    // Skip all cycles in the INS solver --- we advance the state data here.
    ins_hier_integrator->skipCycle(current_time, new_time, cycle_num);

    // Setup Eulerian vectors used in solving the implicit IB equations.
    Pointer<SAMRAIVectorReal<NDIM, double> > eul_sol_vec = new SAMRAIVectorReal<NDIM, double>(
        d_object_name + "::eulerian_sol_vec", d_hierarchy, coarsest_ln, finest_ln);
    eul_sol_vec->addComponent(u_var, u_scratch_idx, wgt_sc_idx, d_hier_velocity_data_ops);
    eul_sol_vec->addComponent(p_var, p_scratch_idx, wgt_cc_idx, d_hier_pressure_data_ops);

    Pointer<SAMRAIVectorReal<NDIM, double> > eul_rhs_vec =
        eul_sol_vec->cloneVector(d_object_name + "::eulerian_rhs_vec");
    eul_rhs_vec->allocateVectorData(current_time);

    d_u_scratch_vec = eul_sol_vec->cloneVector(d_object_name + "::u_scratch_vec");
    d_f_scratch_vec = eul_rhs_vec->cloneVector(d_object_name + "::f_scratch_vec");
    d_u_scratch_vec->allocateVectorData(current_time);
    d_f_scratch_vec->allocateVectorData(current_time);

    ins_hier_integrator->setupSolverVectors(
        eul_sol_vec, eul_rhs_vec, current_time, new_time, cycle_num);

    d_stokes_solver = ins_hier_integrator->getStokesSolver();
    Pointer<KrylovLinearSolver> p_stokes_solver = d_stokes_solver;
    TBOX_ASSERT(p_stokes_solver);
    d_stokes_op = p_stokes_solver->getOperator();
    TBOX_ASSERT(d_stokes_op);

    // Setup Lagrangian vectors used in solving the implicit IB equations.
    Vec lag_sol_petsc_vec, lag_rhs_petsc_vec;
    d_ib_implicit_ops->createSolverVecs(lag_sol_petsc_vec, lag_rhs_petsc_vec);
    d_ib_implicit_ops->setupSolverVecs(lag_sol_petsc_vec, lag_rhs_petsc_vec);

    // Indicate that the current approximation to position of the structure
    // should be used for Lagrangian-Eulerian coupling.
    d_ib_implicit_ops->updateFixedLEOperators();

    // Setup multi-vec objects to store the composite solution and
    // right-hand-side vectors.
    Vec eul_sol_petsc_vec =
        PETScSAMRAIVectorReal::createPETScVector(eul_sol_vec, PETSC_COMM_WORLD);
    Vec eul_rhs_petsc_vec =
        PETScSAMRAIVectorReal::createPETScVector(eul_rhs_vec, PETSC_COMM_WORLD);

    Vec sol_petsc_vecs[] = { eul_sol_petsc_vec, lag_sol_petsc_vec };
    Vec rhs_petsc_vecs[] = { eul_rhs_petsc_vec, lag_rhs_petsc_vec };

    Vec composite_sol_petsc_vec, composite_rhs_petsc_vec, composite_res_petsc_vec;
    ierr = VecCreateMultiVec(PETSC_COMM_WORLD, 2, sol_petsc_vecs, &composite_sol_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = VecCreateMultiVec(PETSC_COMM_WORLD, 2, rhs_petsc_vecs, &composite_rhs_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = VecDuplicate(composite_rhs_petsc_vec, &composite_res_petsc_vec);
    IBTK_CHKERRQ(ierr);

    // Solve the implicit IB equations.
    d_ib_implicit_ops->preprocessSolveFluidEquations(current_time, new_time, cycle_num);

    SNES snes;
    ierr = SNESCreate(PETSC_COMM_WORLD, &snes);
    IBTK_CHKERRQ(ierr);
    ierr = SNESSetFunction(snes, composite_res_petsc_vec, compositeIBFunction_SAMRAI, this);
    IBTK_CHKERRQ(ierr);
    ierr = SNESSetOptionsPrefix(snes, "ib_");
    IBTK_CHKERRQ(ierr);

    Mat jac;
    ierr = VecGetLocalSize(composite_sol_petsc_vec, &n_local);
    IBTK_CHKERRQ(ierr);
    ierr = MatCreateShell(
        PETSC_COMM_WORLD, n_local, n_local, PETSC_DETERMINE, PETSC_DETERMINE, this, &jac);
    IBTK_CHKERRQ(ierr);
    ierr = MatShellSetOperation(
        jac, MATOP_MULT, reinterpret_cast<void (*)(void)>(compositeIBJacobianApply_SAMRAI));
    IBTK_CHKERRQ(ierr);
    ierr = SNESSetJacobian(snes, jac, jac, compositeIBJacobianSetup_SAMRAI, this);
    IBTK_CHKERRQ(ierr);

    Mat schur;
    ierr = VecGetLocalSize(lag_sol_petsc_vec, &n_local);
    IBTK_CHKERRQ(ierr);
    ierr = MatCreateShell(
        PETSC_COMM_WORLD, n_local, n_local, PETSC_DETERMINE, PETSC_DETERMINE, this, &schur);
    IBTK_CHKERRQ(ierr);
    ierr = MatShellSetOperation(
        schur, MATOP_MULT, reinterpret_cast<void (*)(void)>(lagrangianSchurApply_SAMRAI));
    IBTK_CHKERRQ(ierr);
    ierr = KSPCreate(PETSC_COMM_WORLD, &d_schur_solver);
    IBTK_CHKERRQ(ierr);
    ierr = KSPSetOptionsPrefix(d_schur_solver, "ib_schur_");
    IBTK_CHKERRQ(ierr);
    ierr = KSPSetOperators(d_schur_solver, schur, schur, SAME_PRECONDITIONER);
    IBTK_CHKERRQ(ierr);
    PC schur_pc;
    ierr = KSPGetPC(d_schur_solver, &schur_pc);
    IBTK_CHKERRQ(ierr);
    ierr = PCSetType(schur_pc, PCNONE);
    IBTK_CHKERRQ(ierr);
    ierr = KSPSetFromOptions(d_schur_solver);
    IBTK_CHKERRQ(ierr);

    KSP snes_ksp;
    ierr = SNESGetKSP(snes, &snes_ksp);
    IBTK_CHKERRQ(ierr);
    ierr = KSPSetType(snes_ksp, KSPFGMRES);
    IBTK_CHKERRQ(ierr);
    PC snes_pc;
    ierr = KSPGetPC(snes_ksp, &snes_pc);
    IBTK_CHKERRQ(ierr);
    ierr = PCSetType(snes_pc, PCSHELL);
    IBTK_CHKERRQ(ierr);
    ierr = PCShellSetContext(snes_pc, this);
    IBTK_CHKERRQ(ierr);
    ierr = PCShellSetApply(snes_pc, compositeIBPCApply_SAMRAI);
    IBTK_CHKERRQ(ierr);

    ierr = SNESSetFromOptions(snes);
    IBTK_CHKERRQ(ierr);
    ierr = SNESSolve(snes, composite_rhs_petsc_vec, composite_sol_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = SNESDestroy(&snes);
    IBTK_CHKERRQ(ierr);
    ierr = MatDestroy(&jac);
    IBTK_CHKERRQ(ierr);
    ierr = MatDestroy(&schur);
    IBTK_CHKERRQ(ierr);
    ierr = KSPDestroy(&d_schur_solver);
    IBTK_CHKERRQ(ierr);

    d_ib_implicit_ops->postprocessSolveFluidEquations(current_time, new_time, cycle_num);

    // Reset Eulerian solver vectors and Eulerian state data.
    ins_hier_integrator->resetSolverVectors(
        eul_sol_vec, eul_rhs_vec, current_time, new_time, cycle_num);

    // Interpolate the Eulerian velocity to the curvilinear mesh.
    d_ib_implicit_ops->setUpdatedPosition(lag_sol_petsc_vec);
#if 0
    d_hier_velocity_data_ops->linearSum(d_u_idx, 0.5, u_current_idx, 0.5, u_new_idx);
    if (d_enable_logging) plog << d_object_name << "::integrateHierarchy(): interpolating Eulerian velocity to the Lagrangian mesh\n";
    d_ib_implicit_ops->interpolateVelocity(d_u_idx, getCoarsenSchedules(d_object_name+"::u::CONSERVATIVE_COARSEN"), getGhostfillRefineSchedules(d_object_name+"::u"), half_time);

    // Compute the final value of the updated positions of the Lagrangian
    // structure.
    d_ib_implicit_ops->midpointStep(current_time, new_time);
#endif

    // Deallocate temporary data.
    ierr = VecDestroy(&composite_sol_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = VecDestroy(&composite_rhs_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = VecDestroy(&composite_res_petsc_vec);
    IBTK_CHKERRQ(ierr);
    PETScSAMRAIVectorReal::destroyPETScVector(eul_sol_petsc_vec);
    PETScSAMRAIVectorReal::destroyPETScVector(eul_rhs_petsc_vec);
    eul_rhs_vec->freeVectorComponents();
    d_u_scratch_vec->freeVectorComponents();
    d_f_scratch_vec->freeVectorComponents();
    ierr = VecDestroy(&lag_sol_petsc_vec);
    IBTK_CHKERRQ(ierr);
    ierr = VecDestroy(&lag_rhs_petsc_vec);
    IBTK_CHKERRQ(ierr);

    // Execute any registered callbacks.
    executeIntegrateHierarchyCallbackFcns(current_time, new_time, cycle_num);
    return;
} // integrateHierarchy