Skip to content

Commit

Permalink
rename the variable degree more specificly as polynomial_degree
Browse files Browse the repository at this point in the history
  • Loading branch information
yinghe616 committed Jul 27, 2016
1 parent e3e94aa commit 1566b6c
Show file tree
Hide file tree
Showing 5 changed files with 28 additions and 21 deletions.
13 changes: 10 additions & 3 deletions include/aspect/introspection.h
Original file line number Diff line number Diff line change
Expand Up @@ -188,17 +188,24 @@ namespace aspect
*/
const BaseElements base_elements;

struct Degrees
/**
* A structure that contains the polynomial degree of the finite element
* that correspond to each of the variables in this problem.
*
* If there are compositional fields, they are all discretized with the
* same polynomial degree and, consequently, we only need a single integer.
*/
struct PolynomialDegree
{
unsigned int velocities;
unsigned int temperature;
unsigned int compositional_fields;
};
/**
* A variable that enumerates the base elements of the finite element
* A variable that enumerates the polynomial degree of the finite element
* that correspond to each of the variables in this problem.
*/
const Degrees degrees;
const PolynomialDegree polynomial_degree;

/**
* A structure that contains component masks for each of the variables
Expand Down
4 changes: 2 additions & 2 deletions include/aspect/simulator.h
Original file line number Diff line number Diff line change
Expand Up @@ -266,9 +266,9 @@ namespace aspect

/**
* Look up the polynomial degree order for this temperature or compositional
* field. See Introspection::degree for more information.
* field. See Introspection::polynomial_degree for more information.
*/
unsigned int degree(const Introspection<dim> &introspection) const;
unsigned int polynomial_degree(const Introspection<dim> &introspection) const;
};


Expand Down
4 changes: 2 additions & 2 deletions source/simulator/assembly.cc
Original file line number Diff line number Diff line change
Expand Up @@ -697,7 +697,7 @@ namespace aspect
AdvectionSystem<dim> scratch (finite_element,
finite_element.base_element(advection_field.base_element(introspection)),
*mapping,
QGauss<dim>(advection_field.degree(introspection)
QGauss<dim>(advection_field.polynomial_degree(introspection)
+
(parameters.stokes_velocity_degree+1)/2),
/* Because we can only get here in the continuous case, which never requires
Expand Down Expand Up @@ -3066,7 +3066,7 @@ namespace aspect
// (Note: All compositional fields have the same base element and therefore
// the same composition_degree. Thus, we do not need to find out the degree
// of the current field, but use the global instead)
const unsigned int advection_quadrature_degree = advection_field.degree(introspection)
const unsigned int advection_quadrature_degree = advection_field.polynomial_degree(introspection)
+
(parameters.stokes_velocity_degree+1)/2;

Expand Down
12 changes: 6 additions & 6 deletions source/simulator/helper_functions.cc
Original file line number Diff line number Diff line change
Expand Up @@ -142,12 +142,12 @@ namespace aspect

template <int dim>
unsigned int
Simulator<dim>::AdvectionField::degree(const Introspection<dim> &introspection) const
Simulator<dim>::AdvectionField::polynomial_degree(const Introspection<dim> &introspection) const
{
if (this->is_temperature())
return introspection.degrees.temperature;
return introspection.polynomial_degree.temperature;
else
return introspection.degrees.compositional_fields;
return introspection.polynomial_degree.compositional_fields;
}


Expand Down Expand Up @@ -970,8 +970,8 @@ namespace aspect
* in 2D: the combinations are 21, 12
* in 3D: the combinations are 211, 121, 112
*/
const QGauss<1> quadrature_formula_1 (advection_field.degree(introspection)+1);
const QGaussLobatto<1> quadrature_formula_2 (advection_field.degree(introspection)+1);
const QGauss<1> quadrature_formula_1 (advection_field.polynomial_degree(introspection)+1);
const QGaussLobatto<1> quadrature_formula_2 (advection_field.polynomial_degree(introspection)+1);

const unsigned int n_q_points_1 = quadrature_formula_1.size();
const unsigned int n_q_points_2 = quadrature_formula_2.size();
Expand Down Expand Up @@ -1066,7 +1066,7 @@ namespace aspect
Quadrature<dim> quadrature_formula(quadrature_points);

// Quadrature rules only used for the numerical integration for better accuracy
const QGauss<dim> quadrature_formula_0 (advection_field.degree(introspection)+1);
const QGauss<dim> quadrature_formula_0 (advection_field.polynomial_degree(introspection)+1);

const unsigned int n_q_points_0 = quadrature_formula_0.size();

Expand Down
16 changes: 8 additions & 8 deletions source/simulator/introspection.cc
Original file line number Diff line number Diff line change
Expand Up @@ -86,16 +86,16 @@ namespace aspect
}

template <int dim>
typename Introspection<dim>::Degrees
setup_degrees (const Parameters<dim> &parameters)
typename Introspection<dim>::PolynomialDegree
setup_polynomial_degree (const Parameters<dim> &parameters)
{
typename Introspection<dim>::Degrees degrees;
typename Introspection<dim>::PolynomialDegree polynomial_degree;

degrees.velocities = parameters.stokes_velocity_degree;
degrees.temperature = parameters.temperature_degree;
degrees.compositional_fields = parameters.composition_degree;
polynomial_degree.velocities = parameters.stokes_velocity_degree;
polynomial_degree.temperature = parameters.temperature_degree;
polynomial_degree.compositional_fields = parameters.composition_degree;

return degrees;
return polynomial_degree;
}

template <int dim>
Expand Down Expand Up @@ -180,7 +180,7 @@ namespace aspect
block_indices (internal::setup_blocks<dim>(*this)),
extractors (component_indices),
base_elements (internal::setup_base_elements<dim>(*this)),
degrees (internal::setup_degrees<dim>(parameters)),
polynomial_degree (internal::setup_polynomial_degree<dim>(parameters)),
component_masks (*this),
system_dofs_per_block (n_blocks),
composition_names(parameters.names_of_compositional_fields)
Expand Down

0 comments on commit 1566b6c

Please sign in to comment.