/* ---------------------------------------------------------------------
*
- * Copyright (C) 2010 - 2019 by the deal.II authors and
+ * Copyright (C) 2010 - 2020 by the deal.II authors and
* & Jean-Paul Pelteret and Andrew McBride
*
* This file is part of the deal.II library.
scratch.fe_values[J_fe].get_function_values(
scratch.solution_total, scratch.solution_values_J_total);
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point :
+ scratch.fe_values.quadrature_point_indices())
lqph[q_point]->update_values(scratch.solution_grads_u_total[q_point],
scratch.solution_values_p_total[q_point],
scratch.solution_values_J_total[q_point]);
quadrature_point_history.get_data(cell);
Assert(lqph.size() == n_q_points, ExcInternalError());
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point : fe_values.quadrature_point_indices())
{
const double det_F_qp = lqph[q_point]->get_det_F();
const double JxW = fe_values.JxW(q_point);
quadrature_point_history.get_data(cell);
Assert(lqph.size() == n_q_points, ExcInternalError());
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point : fe_values.quadrature_point_indices())
{
const double det_F_qp = lqph[q_point]->get_det_F();
const double J_tilde_qp = lqph[q_point]->get_J_tilde();
quadrature_point_history.get_data(cell);
Assert(lqph.size() == n_q_points, ExcInternalError());
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point :
+ scratch.fe_values.quadrature_point_indices())
{
const Tensor<2, dim> F_inv = lqph[q_point]->get_F_inv();
- for (unsigned int k = 0; k < dofs_per_cell; ++k)
+ for (const unsigned int k : scratch.fe_values.dof_indices())
{
const unsigned int k_group = fe.system_to_base_index(k).first.first;
//
// In doing so, we first extract some configuration dependent variables
// from our quadrature history objects for the current quadrature point.
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point :
+ scratch.fe_values.quadrature_point_indices())
{
const Tensor<2, dim> tau = lqph[q_point]->get_tau();
const SymmetricTensor<4, dim> Jc = lqph[q_point]->get_Jc();
const std::vector<Tensor<2, dim>> &grad_Nx = scratch.grad_Nx[q_point];
const double JxW = scratch.fe_values.JxW(q_point);
- for (unsigned int i = 0; i < dofs_per_cell; ++i)
+ for (const unsigned int i : scratch.fe_values.dof_indices())
{
const unsigned int component_i =
fe.system_to_component_index(i).first;
const unsigned int i_group = fe.system_to_base_index(i).first.first;
- for (unsigned int j = 0; j <= i; ++j)
+ for (const unsigned int j :
+ scratch.fe_values.dof_indices_ending_at(i))
{
const unsigned int component_j =
fe.system_to_component_index(j).first;
// Finally, we need to copy the lower half of the local matrix into the
// upper half:
- for (unsigned int i = 0; i < dofs_per_cell; ++i)
- for (unsigned int j = i + 1; j < dofs_per_cell; ++j)
+ for (const unsigned int i : scratch.fe_values.dof_indices())
+ for (const unsigned int j :
+ scratch.fe_values.dof_indices_starting_at(i + 1))
data.cell_matrix(i, j) = data.cell_matrix(j, i);
}
quadrature_point_history.get_data(cell);
Assert(lqph.size() == n_q_points, ExcInternalError());
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point :
+ scratch.fe_values.quadrature_point_indices())
{
const Tensor<2, dim> F_inv = lqph[q_point]->get_F_inv();
- for (unsigned int k = 0; k < dofs_per_cell; ++k)
+ for (const unsigned int k : scratch.fe_values.dof_indices())
{
const unsigned int k_group = fe.system_to_base_index(k).first.first;
}
}
- for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+ for (const unsigned int q_point :
+ scratch.fe_values.quadrature_point_indices())
{
const SymmetricTensor<2, dim> tau = lqph[q_point]->get_tau();
const double det_F = lqph[q_point]->get_det_F();
// definition of the rhs as the negative
// of the residual, these contributions
// are subtracted.
- for (unsigned int i = 0; i < dofs_per_cell; ++i)
+ for (const unsigned int i : scratch.fe_values.dof_indices())
{
const unsigned int i_group = fe.system_to_base_index(i).first.first;
{
scratch.fe_face_values.reinit(cell, face);
- for (unsigned int f_q_point = 0; f_q_point < n_q_points_f;
- ++f_q_point)
+ for (const unsigned int f_q_point :
+ scratch.fe_face_values.quadrature_point_indices())
{
const Tensor<1, dim> &N =
scratch.fe_face_values.normal_vector(f_q_point);
const double pressure = p0 * parameters.p_p0 * time_ramp;
const Tensor<1, dim> traction = pressure * N;
- for (unsigned int i = 0; i < dofs_per_cell; ++i)
+ for (const unsigned int i : scratch.fe_values.dof_indices())
{
const unsigned int i_group =
fe.system_to_base_index(i).first.first;