From: Wolfgang Bangerth Date: Fri, 1 May 2020 03:20:31 +0000 (-0600) Subject: Convert steps 4,5,6 to FEValues::quadrature_point_indices(). X-Git-Tag: v9.2.0-rc1~152^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=49de7b35f9b0cc2b7f66d114947531775e6775b5;p=dealii.git Convert steps 4,5,6 to FEValues::quadrature_point_indices(). --- diff --git a/examples/step-3/step-3.cc b/examples/step-3/step-3.cc index 71472edbcf..91d6337d92 100644 --- a/examples/step-3/step-3.cc +++ b/examples/step-3/step-3.cc @@ -339,8 +339,8 @@ void Step3::assemble_system() // but maybe not so in higher level languages like C++, but serve // the current purpose quite well. - // For use further down below, we define two shortcuts for values that will - // be used very frequently. First, an abbreviation for the number of degrees + // For use further down below, we define a shortcut for a value that will + // be used very frequently. Namely, an abbreviation for the number of degrees // of freedom on each cell (since we are in 2D and degrees of freedom are // associated with vertices only, this number is four, but we rather want to // write the definition of this variable in a way that does not preclude us @@ -348,19 +348,23 @@ void Step3::assemble_system() // number of degrees of freedom per cell, or work in a different space // dimension). // - // Secondly, we also define an abbreviation for the number of quadrature - // points (here that should be four). In general, it is a good idea to use - // their symbolic names instead of hard-coding these numbers even if you know - // them, since you may want to change the quadrature formula and/or finite - // element at some time; the program will just work with these changes, - // without the need to change anything in this function. + // In general, it is a good idea to use a symbolic name instead of + // hard-coding these numbers even if you know them, since for example, + // you may want to change the finite element at some time. Changing the + // element would have to be done in a different function and it is easy + // to forget to make a corresponding change in another part of the program. + // It is better to not rely on your own calculations, but instead ask + // the right object for the information: Here, we ask the finite element + // to tell us about the number of degrees of freedom per cell and we + // will get the correct number regardless of the space dimension or + // polynomial degree we may have chosen elsewhere in the program. // - // The shortcuts, finally, are only defined to make the following loops a - // bit more readable. You will see them in many places in larger programs, - // and `dofs_per_cell` and `n_q_points` are more or less by convention the - // standard names for these purposes: + // The shortcut here, defined primarily to discuss the basic concept + // and not because it saves a lot of typing, will then make the following + // loops a bit more readable. You will see such shortcuts in many places in + // larger programs, and `dofs_per_cell` is one that is more or less the + // conventional name for this kind of object. const unsigned int dofs_per_cell = fe.dofs_per_cell; - const unsigned int n_q_points = quadrature_formula.size(); // Now, we said that we wanted to assemble the global matrix and vector // cell-by-cell. We could write the results directly into the global matrix, @@ -414,7 +418,7 @@ void Step3::assemble_system() // Now it is time to start integration over the cell, which we // do by looping over all quadrature points, which we will // number by q_index. - for (unsigned int q_index = 0; q_index < n_q_points; ++q_index) + for (const unsigned int q_index : fe_values.quadrature_point_indices()) { // First assemble the matrix: For the Laplace problem, the // matrix on each cell is the integral over the gradients of diff --git a/examples/step-4/step-4.cc b/examples/step-4/step-4.cc index c68a664082..324a55cf83 100644 --- a/examples/step-4/step-4.cc +++ b/examples/step-4/step-4.cc @@ -337,12 +337,12 @@ void Step4::assemble_system() update_values | update_gradients | update_quadrature_points | update_JxW_values); - // We then again define a few abbreviations. The values of these variables - // of course depend on the dimension which we are presently using. However, - // the FE and Quadrature classes do all the necessary work for you and you - // don't have to care about the dimension dependent parts: + // We then again define the same abbreviation as in the previous program. + // The value of this variable of course depends on the dimension which we + // are presently using, but the FiniteElement class does all the necessary + // work for you and you don't have to care about the dimension dependent + // parts: const unsigned int dofs_per_cell = fe.dofs_per_cell; - const unsigned int n_q_points = quadrature_formula.size(); FullMatrix cell_matrix(dofs_per_cell, dofs_per_cell); Vector cell_rhs(dofs_per_cell); @@ -372,7 +372,7 @@ void Step4::assemble_system() // difference to how we did things in step-3: Instead of using a // constant right hand side with value 1, we use the object representing // the right hand side and evaluate it at the quadrature points: - for (unsigned int q_index = 0; q_index < n_q_points; ++q_index) + for (const unsigned int q_index : fe_values.quadrature_point_indices()) for (unsigned int i = 0; i < dofs_per_cell; ++i) { for (unsigned int j = 0; j < dofs_per_cell; ++j) diff --git a/examples/step-5/step-5.cc b/examples/step-5/step-5.cc index 4da722c5b1..e3f6f9a0c1 100644 --- a/examples/step-5/step-5.cc +++ b/examples/step-5/step-5.cc @@ -170,7 +170,6 @@ void Step5::assemble_system() update_quadrature_points | update_JxW_values); const unsigned int dofs_per_cell = fe.dofs_per_cell; - const unsigned int n_q_points = quadrature_formula.size(); FullMatrix cell_matrix(dofs_per_cell, dofs_per_cell); Vector cell_rhs(dofs_per_cell); @@ -189,7 +188,7 @@ void Step5::assemble_system() fe_values.reinit(cell); - for (unsigned int q_index = 0; q_index < n_q_points; ++q_index) + for (const unsigned int q_index : fe_values.quadrature_point_indices()) { const double current_coefficient = coefficient(fe_values.quadrature_point(q_index));