void
INSCollocatedPPMConvectiveOperator::applyConvectiveOperator(
    const int U_idx,
    const int N_idx)
{
    IBAMR_TIMER_START(t_apply_convective_operator);
#ifdef DEBUG_CHECK_ASSERTIONS
    if (!d_is_initialized)
    {
        TBOX_ERROR("INSCollocatedPPMConvectiveOperator::applyConvectiveOperator():\n"
                   << "  operator must be initialized prior to call to applyConvectiveOperator\n");
    }
#endif

    // Setup communications algorithm.
    Pointer<CartesianGridGeometry<NDIM> > grid_geom = d_hierarchy->getGridGeometry();
    Pointer<RefineAlgorithm<NDIM> > refine_alg = new RefineAlgorithm<NDIM>();
    Pointer<RefineOperator<NDIM> > refine_op = grid_geom->lookupRefineOperator(d_U_var, "CONSERVATIVE_LINEAR_REFINE");
    refine_alg->registerRefine(d_U_scratch_idx, U_idx, d_U_scratch_idx, refine_op);

    // Extrapolate from cell centers to cell faces.
    for (int ln = d_coarsest_ln; ln <= d_finest_ln; ++ln)
    {
        refine_alg->resetSchedule(d_ghostfill_scheds[ln]);
        d_ghostfill_scheds[ln]->fillData(d_solution_time);
        d_ghostfill_alg->resetSchedule(d_ghostfill_scheds[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 IntVector<NDIM>& patch_lower = patch_box.lower();
            const IntVector<NDIM>& patch_upper = patch_box.upper();

            Pointer<CellData<NDIM,double> > U_data = patch->getPatchData(d_U_scratch_idx);
            const IntVector<NDIM>& U_data_gcw = U_data->getGhostCellWidth();
#ifdef DEBUG_CHECK_ASSERTIONS
            TBOX_ASSERT(U_data_gcw.min() == U_data_gcw.max());
#endif
            Pointer<FaceData<NDIM,double> > u_ADV_data = patch->getPatchData(d_u_idx);
            const IntVector<NDIM>& u_ADV_data_gcw = u_ADV_data->getGhostCellWidth();
#ifdef DEBUG_CHECK_ASSERTIONS
            TBOX_ASSERT(u_ADV_data_gcw.min() == u_ADV_data_gcw.max());
#endif
            Pointer<FaceData<NDIM,double> > u_extrap_data = patch->getPatchData(d_u_extrap_idx);
            const IntVector<NDIM>& u_extrap_data_gcw = u_extrap_data->getGhostCellWidth();
#ifdef DEBUG_CHECK_ASSERTIONS
            TBOX_ASSERT(u_extrap_data_gcw.min() == u_extrap_data_gcw.max());
#endif
            CellData<NDIM,double>& U0_data = *U_data;
            CellData<NDIM,double>  U1_data(patch_box, 1, U_data_gcw);
#if (NDIM == 3)
            CellData<NDIM,double>  U2_data(patch_box, 1, U_data_gcw);
#endif
            CellData<NDIM,double>  dU_data(patch_box, 1, U_data_gcw);
            CellData<NDIM,double> U_L_data(patch_box, 1, U_data_gcw);
            CellData<NDIM,double> U_R_data(patch_box, 1, U_data_gcw);

            // Extrapolate from cell centers to cell faces.
            for (unsigned int axis = 0; axis < NDIM; ++axis)
            {
                GODUNOV_EXTRAPOLATE_FC(
#if (NDIM == 2)
                    patch_lower(0), patch_upper(0),
                    patch_lower(1), patch_upper(1),
                    U_data_gcw(0), U_data_gcw(1),
                    U0_data.getPointer(axis), U1_data.getPointer(),
                    dU_data.getPointer(), U_L_data.getPointer(), U_R_data.getPointer(),
                    u_ADV_data_gcw   (0), u_ADV_data_gcw   (1),
                    u_extrap_data_gcw(0), u_extrap_data_gcw(1),
                    u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),
                    u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis)
#endif
#if (NDIM == 3)
                    patch_lower(0), patch_upper(0),
                    patch_lower(1), patch_upper(1),
                    patch_lower(2), patch_upper(2),
                    U_data_gcw(0), U_data_gcw(1), U_data_gcw(2),
                    U0_data.getPointer(axis), U1_data.getPointer(), U2_data.getPointer(),
                    dU_data.getPointer(), U_L_data.getPointer(), U_R_data.getPointer(),
                    u_ADV_data_gcw   (0), u_ADV_data_gcw   (1), u_ADV_data_gcw   (2),
                    u_extrap_data_gcw(0), u_extrap_data_gcw(1), u_extrap_data_gcw(2),
                    u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),      u_ADV_data   ->getPointer(2),
                    u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis), u_extrap_data->getPointer(2,axis)
#endif
                                       );
            }

            // If we are using conservative or skew-symmetric differencing,
            // compute the advective fluxes.  These need to be synchronized on
            // the patch hierarchy.
            if (d_difference_form == CONSERVATIVE || d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM,double> > u_ADV_data = patch->getPatchData(d_u_idx);
                const IntVector<NDIM>& u_ADV_data_gcw = u_ADV_data->getGhostCellWidth();
                Pointer<FaceData<NDIM,double> > u_flux_data = patch->getPatchData(d_u_flux_idx);
                const IntVector<NDIM>& u_flux_data_gcw = u_flux_data->getGhostCellWidth();
                for (unsigned int axis = 0; axis < NDIM; ++axis)
                {
                    static const double dt = 1.0;
                    ADVECT_FLUX_FC(
                        dt,
#if (NDIM == 2)
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
//                      u_extrap_data_gcw(0), u_extrap_data_gcw(1),
                        u_ADV_data_gcw   (0), u_ADV_data_gcw   (1),
                        u_extrap_data_gcw(0), u_extrap_data_gcw(1),
                        u_flux_data_gcw  (0), u_flux_data_gcw  (1),
//                      u_extrap_data->getPointer(0,0),    u_extrap_data->getPointer(1,1),
                        u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),
                        u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis),
                        u_flux_data  ->getPointer(0,axis), u_flux_data  ->getPointer(1,axis)
#endif
#if (NDIM == 3)
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
                        patch_lower(2), patch_upper(2),
//                      u_extrap_data_gcw(0), u_extrap_data_gcw(1), u_extrap_data_gcw(2),
                        u_ADV_data_gcw   (0), u_ADV_data_gcw   (1), u_ADV_data_gcw   (2),
                        u_extrap_data_gcw(0), u_extrap_data_gcw(1), u_extrap_data_gcw(2),
                        u_flux_data_gcw  (0), u_flux_data_gcw  (1), u_flux_data_gcw  (2),
//                      u_extrap_data->getPointer(0,0),    u_extrap_data->getPointer(1,1),    u_extrap_data->getPointer(2,2),
                        u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),      u_ADV_data   ->getPointer(2),
                        u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis), u_extrap_data->getPointer(2,axis),
                        u_flux_data  ->getPointer(0,axis), u_flux_data  ->getPointer(1,axis), u_flux_data  ->getPointer(2,axis)
#endif
                                   );
                }
            }
        }
    }

    // Synchronize data on the patch hierarchy.
    for (int ln = d_finest_ln; ln > d_coarsest_ln; --ln)
    {
        d_coarsen_scheds[ln]->coarsenData();
    }

    // Difference values on the patches.
    for (int ln = d_coarsest_ln; ln <= d_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 IntVector<NDIM>& patch_lower = patch_box.lower();
            const IntVector<NDIM>& patch_upper = patch_box.upper();

            const Pointer<CartesianPatchGeometry<NDIM> > patch_geom = patch->getPatchGeometry();
            const double* const dx = patch_geom->getDx();

            Pointer<CellData<NDIM,double> > N_data = patch->getPatchData(N_idx);
            const IntVector<NDIM>& N_data_gcw = N_data->getGhostCellWidth();

            if (d_difference_form == ADVECTIVE || d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM,double> > u_ADV_data = patch->getPatchData(d_u_idx);
                const IntVector<NDIM>& u_ADV_data_gcw = u_ADV_data->getGhostCellWidth();
                Pointer<FaceData<NDIM,double> > u_extrap_data = patch->getPatchData(d_u_extrap_idx);
                const IntVector<NDIM>& u_extrap_data_gcw = u_extrap_data->getGhostCellWidth();
                for (unsigned int axis = 0; axis < NDIM; ++axis)
                {
                    ADVECT_DERIVATIVE_FC(
                        dx,
#if (NDIM == 2)
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
//                      u_extrap_data_gcw(0), u_extrap_data_gcw(1),
                        u_ADV_data_gcw   (0), u_ADV_data_gcw   (1),
                        u_extrap_data_gcw(0), u_extrap_data_gcw(1),
//                      u_extrap_data->getPointer(0,0),    u_extrap_data->getPointer(1,1),
                        u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),
                        u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis),
                        N_data_gcw(0), N_data_gcw(1),
#endif
#if (NDIM == 3)
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
                        patch_lower(2), patch_upper(2),
//                      u_extrap_data_gcw(0), u_extrap_data_gcw(1), u_extrap_data_gcw(2),
                        u_ADV_data_gcw   (0), u_ADV_data_gcw   (1), u_ADV_data_gcw   (2),
                        u_extrap_data_gcw(0), u_extrap_data_gcw(1), u_extrap_data_gcw(2),
//                      u_extrap_data->getPointer(0,0),    u_extrap_data->getPointer(1,1),    u_extrap_data->getPointer(2,2),
                        u_ADV_data   ->getPointer(0),      u_ADV_data   ->getPointer(1),      u_ADV_data   ->getPointer(2),
                        u_extrap_data->getPointer(0,axis), u_extrap_data->getPointer(1,axis), u_extrap_data->getPointer(2,axis),
                        N_data_gcw(0), N_data_gcw(1), N_data_gcw(2),
#endif
                        N_data->getPointer(axis));
                }
            }

            if (d_difference_form == CONSERVATIVE)
            {
                Pointer<FaceData<NDIM,double> > u_flux_data = patch->getPatchData(d_u_flux_idx);
                const IntVector<NDIM>& u_flux_data_gcw = u_flux_data->getGhostCellWidth();
                for (unsigned int axis = 0; axis < NDIM; ++axis)
                {
                    static const double alpha = 1.0;
                    F_TO_C_DIV_FC(
                        N_data->getPointer(axis), N_data_gcw.min(),
                        alpha,
#if (NDIM == 2)
                        u_flux_data->getPointer(0,axis), u_flux_data->getPointer(1,axis), u_flux_data_gcw.min(),
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
#endif
#if (NDIM == 3)
                        u_flux_data->getPointer(0,axis), u_flux_data->getPointer(1,axis), u_flux_data->getPointer(2,axis), u_flux_data_gcw.min(),
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
                        patch_lower(2), patch_upper(2),
#endif
                        dx);
                }
            }

            if (d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM,double> > u_flux_data = patch->getPatchData(d_u_flux_idx);
                const IntVector<NDIM>& u_flux_data_gcw = u_flux_data->getGhostCellWidth();
                for (unsigned int axis = 0; axis < NDIM; ++axis)
                {
                    static const double alpha = 0.5;
                    static const double beta  = 0.5;
                    F_TO_C_DIV_ADD_FC(
                        N_data->getPointer(axis), N_data_gcw.min(),
                        alpha,
#if (NDIM == 2)
                        u_flux_data->getPointer(0,axis), u_flux_data->getPointer(1,axis), u_flux_data_gcw.min(),
                        beta,
                        N_data->getPointer(axis), N_data_gcw.min(),
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
#endif
#if (NDIM == 3)
                        u_flux_data->getPointer(0,axis), u_flux_data->getPointer(1,axis), u_flux_data->getPointer(2,axis), u_flux_data_gcw.min(),
                        beta,
                        N_data->getPointer(axis), N_data_gcw.min(),
                        patch_lower(0), patch_upper(0),
                        patch_lower(1), patch_upper(1),
                        patch_lower(2), patch_upper(2),
#endif
                        dx);
                }
            }
        }
    }
    IBAMR_TIMER_STOP(t_apply_convective_operator);
    return;
}// applyConvectiveOperator
void
AdvDiffCenteredConvectiveOperator::applyConvectiveOperator(const int Q_idx, const int N_idx)
{
    IBAMR_TIMER_START(t_apply_convective_operator);
#if !defined(NDEBUG)
    if (!d_is_initialized)
    {
        TBOX_ERROR("AdvDiffCenteredConvectiveOperator::applyConvectiveOperator():\n"
                   << "  operator must be initialized prior to call to applyConvectiveOperator\n");
    }
#endif

    // Allocate scratch data.
    for (int ln = d_coarsest_ln; ln <= d_finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        level->allocatePatchData(d_Q_scratch_idx);
        level->allocatePatchData(d_q_extrap_idx);
        if (d_difference_form == CONSERVATIVE || d_difference_form == SKEW_SYMMETRIC)
            level->allocatePatchData(d_q_flux_idx);
    }

    // Setup communications algorithm.
    Pointer<CartesianGridGeometry<NDIM> > grid_geom = d_hierarchy->getGridGeometry();
    Pointer<RefineAlgorithm<NDIM> > refine_alg = new RefineAlgorithm<NDIM>();
    Pointer<RefineOperator<NDIM> > refine_op = grid_geom->lookupRefineOperator(d_Q_var, "CONSERVATIVE_LINEAR_REFINE");
    refine_alg->registerRefine(d_Q_scratch_idx, Q_idx, d_Q_scratch_idx, refine_op);

    // Extrapolate from cell centers to cell faces.
    for (int ln = d_coarsest_ln; ln <= d_finest_ln; ++ln)
    {
        refine_alg->resetSchedule(d_ghostfill_scheds[ln]);
        d_ghostfill_scheds[ln]->fillData(d_solution_time);
        d_ghostfill_alg->resetSchedule(d_ghostfill_scheds[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 IntVector<NDIM>& patch_lower = patch_box.lower();
            const IntVector<NDIM>& patch_upper = patch_box.upper();

            Pointer<CellData<NDIM, double> > Q_data = patch->getPatchData(d_Q_scratch_idx);
            const IntVector<NDIM>& Q_data_gcw = Q_data->getGhostCellWidth();
#if !defined(NDEBUG)
            TBOX_ASSERT(Q_data_gcw.min() == Q_data_gcw.max());
#endif
            Pointer<FaceData<NDIM, double> > u_ADV_data = patch->getPatchData(d_u_idx);
            const IntVector<NDIM>& u_ADV_data_gcw = u_ADV_data->getGhostCellWidth();
#if !defined(NDEBUG)
            TBOX_ASSERT(u_ADV_data_gcw.min() == u_ADV_data_gcw.max());
#endif
            Pointer<FaceData<NDIM, double> > q_extrap_data = patch->getPatchData(d_q_extrap_idx);
            const IntVector<NDIM>& q_extrap_data_gcw = q_extrap_data->getGhostCellWidth();
#if !defined(NDEBUG)
            TBOX_ASSERT(q_extrap_data_gcw.min() == q_extrap_data_gcw.max());
#endif
            // Enforce physical boundary conditions at inflow boundaries.
            AdvDiffPhysicalBoundaryUtilities::setPhysicalBoundaryConditions(
                Q_data,
                u_ADV_data,
                patch,
                d_bc_coefs,
                d_solution_time,
                /*inflow_boundary_only*/ d_outflow_bdry_extrap_type != "NONE",
                d_homogeneous_bc);

            // Interpolate from cell centers to cell faces.
            for (unsigned int d = 0; d < d_Q_data_depth; ++d)
            {
                C_TO_F_CWISE_INTERP_2ND_FC(
#if (NDIM == 2)
                    q_extrap_data->getPointer(0, d),
                    q_extrap_data->getPointer(1, d),
                    q_extrap_data_gcw.min(),
                    Q_data->getPointer(d),
                    Q_data_gcw.min(),
                    patch_lower(0),
                    patch_upper(0),
                    patch_lower(1),
                    patch_upper(1)
#endif
#if (NDIM == 3)
                        q_extrap_data->getPointer(0, d),
                    q_extrap_data->getPointer(1, d),
                    q_extrap_data->getPointer(2, d),
                    q_extrap_data_gcw.min(),
                    Q_data->getPointer(d),
                    Q_data_gcw.min(),
                    patch_lower(0),
                    patch_upper(0),
                    patch_lower(1),
                    patch_upper(1),
                    patch_lower(2),
                    patch_upper(2)
#endif
                        );
            }

            // If we are using conservative or skew-symmetric differencing,
            // compute the advective fluxes.  These need to be synchronized on
            // the patch hierarchy.
            if (d_difference_form == CONSERVATIVE || d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM, double> > q_flux_data = patch->getPatchData(d_q_flux_idx);
                const IntVector<NDIM>& q_flux_data_gcw = q_flux_data->getGhostCellWidth();
                for (unsigned int d = 0; d < d_Q_data_depth; ++d)
                {
                    static const double dt = 1.0;
                    ADVECT_FLUX_FC(dt,
#if (NDIM == 2)
                                   patch_lower(0),
                                   patch_upper(0),
                                   patch_lower(1),
                                   patch_upper(1),
                                   u_ADV_data_gcw(0),
                                   u_ADV_data_gcw(1),
                                   q_extrap_data_gcw(0),
                                   q_extrap_data_gcw(1),
                                   q_flux_data_gcw(0),
                                   q_flux_data_gcw(1),
                                   u_ADV_data->getPointer(0),
                                   u_ADV_data->getPointer(1),
                                   q_extrap_data->getPointer(0, d),
                                   q_extrap_data->getPointer(1, d),
                                   q_flux_data->getPointer(0, d),
                                   q_flux_data->getPointer(1, d)
#endif
#if (NDIM == 3)
                                       patch_lower(0),
                                   patch_upper(0),
                                   patch_lower(1),
                                   patch_upper(1),
                                   patch_lower(2),
                                   patch_upper(2),
                                   u_ADV_data_gcw(0),
                                   u_ADV_data_gcw(1),
                                   u_ADV_data_gcw(2),
                                   q_extrap_data_gcw(0),
                                   q_extrap_data_gcw(1),
                                   q_extrap_data_gcw(2),
                                   q_flux_data_gcw(0),
                                   q_flux_data_gcw(1),
                                   q_flux_data_gcw(2),
                                   u_ADV_data->getPointer(0),
                                   u_ADV_data->getPointer(1),
                                   u_ADV_data->getPointer(2),
                                   q_extrap_data->getPointer(0, d),
                                   q_extrap_data->getPointer(1, d),
                                   q_extrap_data->getPointer(2, d),
                                   q_flux_data->getPointer(0, d),
                                   q_flux_data->getPointer(1, d),
                                   q_flux_data->getPointer(2, d)
#endif
                                       );
                }
            }
        }
    }

    // Synchronize data on the patch hierarchy.
    for (int ln = d_finest_ln; ln > d_coarsest_ln; --ln)
    {
        d_coarsen_scheds[ln]->coarsenData();
    }

    // Difference values on the patches.
    for (int ln = d_coarsest_ln; ln <= d_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 IntVector<NDIM>& patch_lower = patch_box.lower();
            const IntVector<NDIM>& patch_upper = patch_box.upper();

            const Pointer<CartesianPatchGeometry<NDIM> > patch_geom = patch->getPatchGeometry();
            const double* const dx = patch_geom->getDx();

            Pointer<CellData<NDIM, double> > N_data = patch->getPatchData(N_idx);
            const IntVector<NDIM>& N_data_gcw = N_data->getGhostCellWidth();

            if (d_difference_form == ADVECTIVE || d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM, double> > u_ADV_data = patch->getPatchData(d_u_idx);
                const IntVector<NDIM>& u_ADV_data_gcw = u_ADV_data->getGhostCellWidth();
                Pointer<FaceData<NDIM, double> > q_extrap_data = patch->getPatchData(d_q_extrap_idx);
                const IntVector<NDIM>& q_extrap_data_gcw = q_extrap_data->getGhostCellWidth();
                for (unsigned int d = 0; d < d_Q_data_depth; ++d)
                {
                    ADVECT_DERIVATIVE_FC(dx,
#if (NDIM == 2)
                                         patch_lower(0),
                                         patch_upper(0),
                                         patch_lower(1),
                                         patch_upper(1),
                                         u_ADV_data_gcw(0),
                                         u_ADV_data_gcw(1),
                                         q_extrap_data_gcw(0),
                                         q_extrap_data_gcw(1),
                                         u_ADV_data->getPointer(0),
                                         u_ADV_data->getPointer(1),
                                         q_extrap_data->getPointer(0, d),
                                         q_extrap_data->getPointer(1, d),
                                         N_data_gcw(0),
                                         N_data_gcw(1),
#endif
#if (NDIM == 3)
                                         patch_lower(0),
                                         patch_upper(0),
                                         patch_lower(1),
                                         patch_upper(1),
                                         patch_lower(2),
                                         patch_upper(2),
                                         u_ADV_data_gcw(0),
                                         u_ADV_data_gcw(1),
                                         u_ADV_data_gcw(2),
                                         q_extrap_data_gcw(0),
                                         q_extrap_data_gcw(1),
                                         q_extrap_data_gcw(2),
                                         u_ADV_data->getPointer(0),
                                         u_ADV_data->getPointer(1),
                                         u_ADV_data->getPointer(2),
                                         q_extrap_data->getPointer(0, d),
                                         q_extrap_data->getPointer(1, d),
                                         q_extrap_data->getPointer(2, d),
                                         N_data_gcw(0),
                                         N_data_gcw(1),
                                         N_data_gcw(2),
#endif
                                         N_data->getPointer(d));
                }
            }

            if (d_difference_form == CONSERVATIVE)
            {
                Pointer<FaceData<NDIM, double> > q_flux_data = patch->getPatchData(d_q_flux_idx);
                const IntVector<NDIM>& q_flux_data_gcw = q_flux_data->getGhostCellWidth();
                for (unsigned int d = 0; d < d_Q_data_depth; ++d)
                {
                    static const double alpha = 1.0;
                    F_TO_C_DIV_FC(N_data->getPointer(d),
                                  N_data_gcw.min(),
                                  alpha,
#if (NDIM == 2)
                                  q_flux_data->getPointer(0, d),
                                  q_flux_data->getPointer(1, d),
                                  q_flux_data_gcw.min(),
                                  patch_lower(0),
                                  patch_upper(0),
                                  patch_lower(1),
                                  patch_upper(1),
#endif
#if (NDIM == 3)
                                  q_flux_data->getPointer(0, d),
                                  q_flux_data->getPointer(1, d),
                                  q_flux_data->getPointer(2, d),
                                  q_flux_data_gcw.min(),
                                  patch_lower(0),
                                  patch_upper(0),
                                  patch_lower(1),
                                  patch_upper(1),
                                  patch_lower(2),
                                  patch_upper(2),
#endif
                                  dx);
                }
            }

            if (d_difference_form == SKEW_SYMMETRIC)
            {
                Pointer<FaceData<NDIM, double> > q_flux_data = patch->getPatchData(d_q_flux_idx);
                const IntVector<NDIM>& q_flux_data_gcw = q_flux_data->getGhostCellWidth();
                for (unsigned int d = 0; d < d_Q_data_depth; ++d)
                {
                    static const double alpha = 0.5;
                    static const double beta = 0.5;
                    F_TO_C_DIV_ADD_FC(N_data->getPointer(d),
                                      N_data_gcw.min(),
                                      alpha,
#if (NDIM == 2)
                                      q_flux_data->getPointer(0, d),
                                      q_flux_data->getPointer(1, d),
                                      q_flux_data_gcw.min(),
                                      beta,
                                      N_data->getPointer(d),
                                      N_data_gcw.min(),
                                      patch_lower(0),
                                      patch_upper(0),
                                      patch_lower(1),
                                      patch_upper(1),
#endif
#if (NDIM == 3)
                                      q_flux_data->getPointer(0, d),
                                      q_flux_data->getPointer(1, d),
                                      q_flux_data->getPointer(2, d),
                                      q_flux_data_gcw.min(),
                                      beta,
                                      N_data->getPointer(d),
                                      N_data_gcw.min(),
                                      patch_lower(0),
                                      patch_upper(0),
                                      patch_lower(1),
                                      patch_upper(1),
                                      patch_lower(2),
                                      patch_upper(2),
#endif
                                      dx);
                }
            }
        }
    }

    // Deallocate scratch data.
    for (int ln = d_coarsest_ln; ln <= d_finest_ln; ++ln)
    {
        Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln);
        level->deallocatePatchData(d_Q_scratch_idx);
        level->deallocatePatchData(d_q_extrap_idx);
        if (d_difference_form == CONSERVATIVE || d_difference_form == SKEW_SYMMETRIC)
            level->deallocatePatchData(d_q_flux_idx);
    }

    IBAMR_TIMER_STOP(t_apply_convective_operator);
    return;
} // applyConvectiveOperator