From 9e1296734ad0061c6b71b178d1838ccbed0d25d0 Mon Sep 17 00:00:00 2001 From: Markus Buerg Date: Fri, 14 May 2010 19:48:43 +0000 Subject: [PATCH] Documentation and style changes. git-svn-id: https://svn.dealii.org/trunk@21130 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/deal.II/include/fe/fe_values.h | 87 ++++++++++++++++---------- deal.II/deal.II/source/fe/fe_values.cc | 63 ++++++++++++------- 2 files changed, 95 insertions(+), 55 deletions(-) diff --git a/deal.II/deal.II/include/fe/fe_values.h b/deal.II/deal.II/include/fe/fe_values.h index 5eec718392..4f6ddd7b96 100644 --- a/deal.II/deal.II/include/fe/fe_values.h +++ b/deal.II/deal.II/include/fe/fe_values.h @@ -501,15 +501,14 @@ namespace FEValuesViews typedef double divergence_type; /** - * A typedef for the type of the curl - * of the view this class + * A typedef for the type of the + * curl of the view this class * represents. Here, for a set of - * dim components of the - * finite element, the curl is of type - * Tensor@<1, 1@> (i.e. a - * scalar) in 2d, and of type - * Tensor@<1, 3@> (i.e. a - * vector) in 3d. + * spacedim=2 components + * of the finite element, the curl is + * a Tensor@<1, 1@>. For + * spacedim=3 it is a + * Tensor@<1, dim@>. */ typedef Tensor<1, (spacedim == 3)? 3 : 1> curl_type; @@ -616,11 +615,28 @@ namespace FEValuesViews * the vector components selected by * this view, for the shape function * and quadrature point selected by the - * arguments. + * arguments. For 1d this function does + * not make any sense. Thus it is not + * implemented for spacedim=1. + * In 2d the curl is defined as + * \begin{equation*} + * \operatorname{curl}(u):=\frac{du_2}{dx} + * -\frac{du_1}{dy}, + * \end{equation*} + * whereas in 3d it is given by + * \begin{equation*} + * \operatorname{curl}(u):=\left( + * \begin{array}{c} + * \frac{du_3}{dy}-\frac{du_2}{dz}\\ + * \frac{du_1}{dz}-\frac{du_3}{dx}\\ + * \frac{du_2}{dx}-\frac{du_1}{dy} + * \end{array} + * \right). + * \end{equation*} */ curl_type curl (const unsigned int shape_function, - const unsigned int q_point) const; + const unsigned int q_point) const; /** * Return the Hessian (the tensor of @@ -713,7 +729,7 @@ namespace FEValuesViews * * There is no equivalent function such * as - * FEValuesBase::get_function_divergences + * FEValuesBase::get_function_gradients * in the FEValues classes but the * information can be obtained from * FEValuesBase::get_function_gradients, @@ -736,7 +752,7 @@ namespace FEValuesViews * * There is no equivalent function such * as - * FEValuesBase::get_function_curls + * FEValuesBase::get_function_gradients * in the FEValues classes but the * information can be obtained from * FEValuesBase::get_function_gradients, @@ -744,7 +760,7 @@ namespace FEValuesViews */ template void get_function_curls (const InputVector& fe_function, - std::vector& curls) const; + std::vector& curls) const; /** * Return the Hessians of the selected @@ -3778,19 +3794,15 @@ namespace FEValuesViews template inline typename Vector::curl_type - Vector::curl (const unsigned int shape_function, - const unsigned int q_point) const - { - // this function works like in the case - // above + Vector::curl (const unsigned int shape_function, const unsigned int q_point) const { + // this function works like in the case above typedef FEValuesBase FVB; Assert (shape_function < fe_values.fe->dofs_per_cell, - ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell)); + ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell)); Assert (fe_values.update_flags & update_gradients, - typename FVB::ExcAccessToUninitializedField()); - // same as for the scalar case except - // that we have one more index + typename FVB::ExcAccessToUninitializedField()); + // same as for the scalar case except that we have one more index const int snc = shape_function_data[shape_function].single_nonzero_component; if (snc == -2) @@ -3798,8 +3810,10 @@ namespace FEValuesViews else switch (dim) { - case 1: + case 1: { + Assert (false, ExcMessage("Computing the curl in 1d is not a useful operation")); return curl_type (); + } case 2: { if (snc != -1) { @@ -3824,10 +3838,12 @@ namespace FEValuesViews return_value[0] = 0.0; if (shape_function_data[shape_function].is_nonzero_shape_function_component[0]) - return_value[0] -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1]; + return_value[0] + -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1]; if (shape_function_data[shape_function].is_nonzero_shape_function_component[1]) - return_value[0] += fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0]; + return_value[0] + += fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0]; return return_value; } @@ -3868,27 +3884,32 @@ namespace FEValuesViews return_value[i] = 0.0; if (shape_function_data[shape_function].is_nonzero_shape_function_component[0]) { - return_value[1] += fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2]; - return_value[2] -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1]; + return_value[1] + += fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2]; + return_value[2] + -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1]; } if (shape_function_data[shape_function].is_nonzero_shape_function_component[1]) { - return_value[0] -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2]; - return_value[2] += fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0]; + return_value[0] + -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2]; + return_value[2] + += fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0]; } if (shape_function_data[shape_function].is_nonzero_shape_function_component[2]) { - return_value[0] += fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1]; - return_value[1] -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0]; + return_value[0] + += fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1]; + return_value[1] + -= fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0]; } return return_value; } } } - } +} - template inline typename Vector::hessian_type diff --git a/deal.II/deal.II/source/fe/fe_values.cc b/deal.II/deal.II/source/fe/fe_values.cc index ff01035778..f1dc1bc73a 100644 --- a/deal.II/deal.II/source/fe/fe_values.cc +++ b/deal.II/deal.II/source/fe/fe_values.cc @@ -758,28 +758,34 @@ namespace FEValuesViews template template void - Vector::get_function_curls (const InputVector &fe_function, - std::vector &curls) const { + Vector:: + get_function_curls (const InputVector &fe_function, + std::vector &curls) const { typedef FEValuesBase FVB; Assert (fe_values.update_flags & update_gradients, - typename FVB::ExcAccessToUninitializedField()); + typename FVB::ExcAccessToUninitializedField()); Assert (curls.size() == fe_values.n_quadrature_points, - ExcDimensionMismatch (curls.size(), fe_values.n_quadrature_points)); + ExcDimensionMismatch (curls.size(), fe_values.n_quadrature_points)); Assert (fe_values.present_cell.get () != 0, - ExcMessage ("FEValues object is not reinit'ed to any cell")); + ExcMessage ("FEValues object is not reinited to any cell")); Assert (fe_function.size () == fe_values.present_cell->n_dofs_for_dof_handler (), - ExcDimensionMismatch (fe_function.size (), fe_values.present_cell->n_dofs_for_dof_handler ())); - // get function values of dofs on this - // cell + ExcDimensionMismatch (fe_function.size (), fe_values.present_cell->n_dofs_for_dof_handler ())); + // get function values of dofs on this cell dealii::Vector dof_values (fe_values.dofs_per_cell); fe_values.present_cell->get_interpolated_dof_values (fe_function, dof_values); std::fill (curls.begin (), curls.end (), curl_type ()); switch (dim) { + case 1: { + Assert (false, ExcMessage("Computing the curl in 1d is not a useful operation")); + break; + } + case 2: { - for (unsigned int shape_function = 0; shape_function < fe_values.fe->dofs_per_cell; ++shape_function) { + for (unsigned int shape_function = 0; + shape_function < fe_values.fe->dofs_per_cell; ++shape_function) { const int snc = shape_function_data[shape_function].single_nonzero_component; if (snc == -2) @@ -792,19 +798,22 @@ namespace FEValuesViews continue; if (snc != -1) { - const unsigned int comp = shape_function_data[shape_function].single_nonzero_component_index; + const unsigned int comp = + shape_function_data[shape_function].single_nonzero_component_index; const Tensor<1, spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[snc][0]; switch (shape_function_data[shape_function].single_nonzero_component_index) { case 0: { - for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) + for (unsigned int q_point = 0; + q_point < fe_values.n_quadrature_points; ++q_point) curls[q_point][0] = -1.0 * value * (*shape_gradient_ptr++)[1]; break; } default: - for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) + for (unsigned int q_point = 0; + q_point < fe_values.n_quadrature_points; ++q_point) curls[q_point][0] = value * (*shape_gradient_ptr)[0]; } } @@ -814,14 +823,16 @@ namespace FEValuesViews curls[q_point][0] = 0; if (shape_function_data[shape_function].is_nonzero_shape_function_component[0]) { - const Tensor<1,spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][0]; + const Tensor<1,spacedim> *shape_gradient_ptr = + &fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][0]; for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) curls[q_point][0] -= value * (*shape_gradient_ptr++)[1]; } if (shape_function_data[shape_function].is_nonzero_shape_function_component[1]) { - const Tensor<1,spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][0]; + const Tensor<1,spacedim> *shape_gradient_ptr = + &fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][0]; for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) curls[q_point][0] += value * (*shape_gradient_ptr++)[0]; @@ -832,7 +843,8 @@ namespace FEValuesViews } case 3: { - for (unsigned int shape_function = 0; shape_function < fe_values.fe->dofs_per_cell; ++shape_function) { + for (unsigned int shape_function = 0; + shape_function < fe_values.fe->dofs_per_cell; ++shape_function) { const int snc = shape_function_data[shape_function].single_nonzero_component; if (snc == -2) @@ -845,12 +857,14 @@ namespace FEValuesViews continue; if (snc != -1) { - const unsigned int comp = shape_function_data[shape_function].single_nonzero_component_index; + const unsigned int comp = + shape_function_data[shape_function].single_nonzero_component_index; const Tensor<1, spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[snc][0]; switch (shape_function_data[shape_function].single_nonzero_component_index) { case 0: { - for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { + for (unsigned int q_point = 0; + q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][0] = 0; curls[q_point][1] = value * (*shape_gradient_ptr)[2]; curls[q_point][2] = -1.0 * value * (*shape_gradient_ptr++)[1]; @@ -860,7 +874,8 @@ namespace FEValuesViews } case 1: { - for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { + for (unsigned int q_point = 0; + q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][0] = -1.0 * value * (*shape_gradient_ptr)[2]; curls[q_point][1] = 0; curls[q_point][2] = value * (*shape_gradient_ptr++)[0]; @@ -870,7 +885,8 @@ namespace FEValuesViews } default: - for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { + for (unsigned int q_point = 0; + q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][0] = value * (*shape_gradient_ptr)[1]; curls[q_point][1] = -1.0 * value * (*shape_gradient_ptr++)[0]; curls[q_point][2] = 0; @@ -884,7 +900,8 @@ namespace FEValuesViews curls[q_point][d] = 0; if (shape_function_data[shape_function].is_nonzero_shape_function_component[0]) { - const Tensor<1,spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][0]; + const Tensor<1,spacedim> *shape_gradient_ptr = + &fe_values.shape_gradients[shape_function_data[shape_function].row_index[0]][0]; for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][1] += value * (*shape_gradient_ptr++)[2]; @@ -893,7 +910,8 @@ namespace FEValuesViews } if (shape_function_data[shape_function].is_nonzero_shape_function_component[1]) { - const Tensor<1,spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][0]; + const Tensor<1,spacedim> *shape_gradient_ptr = + &fe_values.shape_gradients[shape_function_data[shape_function].row_index[1]][0]; for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][0] -= value * (*shape_gradient_ptr++)[2]; @@ -902,7 +920,8 @@ namespace FEValuesViews } if (shape_function_data[shape_function].is_nonzero_shape_function_component[2]) { - const Tensor<1,spacedim> *shape_gradient_ptr = &fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][0]; + const Tensor<1,spacedim> *shape_gradient_ptr = + &fe_values.shape_gradients[shape_function_data[shape_function].row_index[2]][0]; for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points; ++q_point) { curls[q_point][0] += value * (*shape_gradient_ptr++)[1]; -- 2.39.5