From: David Wells Date: Thu, 22 Dec 2016 18:28:08 +0000 (-0500) Subject: Add documentation for enumerations. X-Git-Tag: v8.5.0-rc1~315^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=eebf91c67c392bc30d6c0746d348eaefdbac924c;p=dealii.git Add documentation for enumerations. Individual values in an enumeration will not show up in doxygen unless they have a documentation string, so add a little description for each enumeration value. --- diff --git a/include/deal.II/base/exceptions.h b/include/deal.II/base/exceptions.h index ae3e10f9c2..4846a17ee7 100644 --- a/include/deal.II/base/exceptions.h +++ b/include/deal.II/base/exceptions.h @@ -224,26 +224,40 @@ namespace deal_II_exceptions /** * Conditionally abort the program. * - * Depending on whether disable_abort_on_exception was called, this - * function either aborts the program flow by printing the error message - * provided by @p exc and calling std::abort(), or throws @p exc - * instead (if @p nothrow is set to false). + * Depending on whether deal_II_exceptions::disable_abort_on_exception was + * called, this function either aborts the program flow by printing the + * error message provided by @p exc and calling std::abort(), or + * throws @p exc instead (if @p nothrow is set to false). * - * If the boolean @p nothrow is set to true and disable_abort_on_exception - * was called, the exception type is just printed to deallog and program - * flow continues. This is useful if throwing an exception is prohibited + * If the boolean @p nothrow is set to true and + * deal_II_exceptions::disable_abort_on_exception was called, the + * exception type is just printed to deallog and program flow + * continues. This is useful if throwing an exception is prohibited * (e.g. in a destructor with noexcept(true) or * throw()). */ void abort (const ExceptionBase &exc, bool nothrow = false); /** - * An enum describing how to treat an exception in issue_error + * An enum describing how to treat an exception in issue_error. */ enum ExceptionHandling { + /** + * Abort the program by calling std::abort unless + * deal_II_exceptions::disable_abort_on_exception has been called: in + * that case the program will throw an exception. + */ abort_on_exception, + /** + * Throw the exception normally. + */ throw_on_exception, + /** + * Call std::abort as long as + * deal_II_exceptions::disable_abort_on_exception has not been called: + * if it has, then just print a description of the exception to deallog. + */ abort_nothrow_on_exception }; diff --git a/include/deal.II/base/geometry_info.h b/include/deal.II/base/geometry_info.h index 71f850db45..9c52e86764 100644 --- a/include/deal.II/base/geometry_info.h +++ b/include/deal.II/base/geometry_info.h @@ -52,9 +52,21 @@ public: */ enum Object { + /** + * A vertex. + */ vertex = 0, + /** + * A line. + */ line = 1, + /** + * A quadrilateral. + */ quad = 2, + /** + * A hexahedron. + */ hex = 3 }; @@ -139,8 +151,14 @@ struct RefinementPossibilities */ enum Possibilities { - no_refinement= 0, + /** + * Do not perform refinement. + */ + no_refinement = 0, + /** + * Perform isotropic refinement. + */ isotropic_refinement = static_cast(-1) }; }; @@ -196,9 +214,17 @@ struct RefinementPossibilities<1> */ enum Possibilities { - no_refinement= 0, + /** + * Do not refine. + */ + no_refinement = 0, + /** + * Perform a cut in the x-direction. + */ cut_x = 1, - + /** + * Perform isotropic refinement. + */ isotropic_refinement = cut_x }; }; @@ -255,11 +281,26 @@ struct RefinementPossibilities<2> */ enum Possibilities { - no_refinement= 0, - cut_x = 1, - cut_y = 2, - cut_xy = cut_x | cut_y, + /** + * Do not refine. + */ + no_refinement = 0, + /** + * Perform a cut in the x-direction. + */ + cut_x = 1, + /** + * Perform a cut in the y-direction. + */ + cut_y = 2, + /** + * Perform cuts in the x- and y-directions. + */ + cut_xy = cut_x | cut_y, + /** + * Perform isotropic refinement. + */ isotropic_refinement = cut_xy }; }; @@ -316,15 +357,42 @@ struct RefinementPossibilities<3> */ enum Possibilities { - no_refinement= 0, - cut_x = 1, - cut_y = 2, - cut_xy = cut_x | cut_y, - cut_z = 4, - cut_xz = cut_x | cut_z, - cut_yz = cut_y | cut_z, - cut_xyz = cut_x | cut_y | cut_z, + /** + * Do not refine. + */ + no_refinement = 0, + /** + * Perform a cut in the x-direction. + */ + cut_x = 1, + /** + * Perform a cut in the y-direction. + */ + cut_y = 2, + /** + * Perform a cut in the x and y-directions. + */ + cut_xy = cut_x | cut_y, + /** + * Perform a cut in the z-direction. + */ + cut_z = 4, + /** + * Perform a cuts in the x- and y-directions. + */ + cut_xz = cut_x | cut_z, + /** + * Perform a cuts in the x- and y-directions. + */ + cut_yz = cut_y | cut_z, + /** + * Perform a cuts in the x-, y-, and z-directions. + */ + cut_xyz = cut_x | cut_y | cut_z, + /** + * Perform isotropic refinement. + */ isotropic_refinement = cut_xyz }; }; @@ -470,8 +538,14 @@ namespace internal */ enum Possibilities { + /** + * Do not refine. + */ case_none = 0, + /** + * Refine isotropically. + */ case_isotropic = static_cast(-1) }; }; @@ -496,8 +570,14 @@ namespace internal */ enum Possibilities { + /** + * Do not refine. + */ case_none = 0, + /** + * Refine isotropically. + */ case_isotropic = case_none }; }; @@ -524,8 +604,14 @@ namespace internal */ enum Possibilities { + /** + * Do not refine. + */ case_none = 0, + /** + * Refine isotropically. + */ case_isotropic = case_none }; }; @@ -554,9 +640,17 @@ namespace internal */ enum Possibilities { + /** + * Do not refine. + */ case_none = 0, + /** + * Cut in the x-direction. + */ case_x = 1, - + /** + * Refine isotropically. + */ case_isotropic = case_x }; }; diff --git a/include/deal.II/base/parameter_handler.h b/include/deal.II/base/parameter_handler.h index f395c5e463..1178195a18 100644 --- a/include/deal.II/base/parameter_handler.h +++ b/include/deal.II/base/parameter_handler.h @@ -791,7 +791,17 @@ namespace Patterns * Files can be used for input or output. This can be specified in the * constructor by choosing the flag type. */ - enum FileType {input = 0, output = 1}; + enum FileType + { + /** + * Open for input. + */ + input = 0, + /** + * Open for output. + */ + output = 1 + }; /** * Constructor. The type of the file can be specified by choosing the @@ -2548,7 +2558,14 @@ private: */ enum MultipleEntryType { - variant, array + /** + * A variant entry. + */ + variant, + /** + * An array entry. + */ + array }; /** diff --git a/include/deal.II/base/quadrature_lib.h b/include/deal.II/base/quadrature_lib.h index 67c3cf0b35..d9c879223d 100644 --- a/include/deal.II/base/quadrature_lib.h +++ b/include/deal.II/base/quadrature_lib.h @@ -583,7 +583,17 @@ public: /* EndPoint is used to specify which of the two endpoints of the unit interval * is used also as quadrature point */ - enum EndPoint { left,right }; + enum EndPoint + { + /** + * Left end point. + */ + left, + /** + * Right end point. + */ + right + }; /// Generate a formula with n quadrature points QGaussRadauChebyshev(const unsigned int n, EndPoint ep=QGaussRadauChebyshev::left); diff --git a/include/deal.II/base/table_handler.h b/include/deal.II/base/table_handler.h index d5fed63aa1..76cbee05d2 100644 --- a/include/deal.II/base/table_handler.h +++ b/include/deal.II/base/table_handler.h @@ -322,9 +322,22 @@ public: */ enum TextOutputFormat { + /** + * Print the table with headers. + */ table_with_headers, + /** + * Print the table with separate lines for each column label. + */ table_with_separate_column_description, + /** + * Like table_with_separate_column_description, but without aligning the + * column containing the column labels. + */ simple_table_with_separate_column_description, + /** + * Print the table in org mode format. + */ org_mode_table }; diff --git a/include/deal.II/base/timer.h b/include/deal.II/base/timer.h index 0692aabf6e..ac2c028b09 100644 --- a/include/deal.II/base/timer.h +++ b/include/deal.II/base/timer.h @@ -431,9 +431,21 @@ public: */ enum OutputFrequency { + /** + * Generate output after every call. + */ every_call, + /** + * Generate output in summary at the end. + */ summary, + /** + * Generate output both after every call and in summary at the end. + */ every_call_and_summary, + /** + * Never generate any output. + */ never }; @@ -443,8 +455,17 @@ public: */ enum OutputType { + /** + * Output CPU times. + */ cpu_times, + /** + * Output wall clock times. + */ wall_times, + /** + * Output both CPU and wall clock times. + */ cpu_and_wall_times }; diff --git a/include/deal.II/fe/fe_base.h b/include/deal.II/fe/fe_base.h index 97dc0f2835..626fb3177c 100644 --- a/include/deal.II/fe/fe_base.h +++ b/include/deal.II/fe/fe_base.h @@ -94,10 +94,25 @@ namespace FiniteElementDomination */ enum Domination { + /** + * The current element dominates. + */ this_element_dominates, + /** + * The other element dominates. + */ other_element_dominates, + /** + * Neither element dominates. + */ neither_element_dominates, + /** + * Either element may dominate. + */ either_element_can_dominate, + /** + * There are no requirements. + */ no_requirements }; diff --git a/include/deal.II/fe/fe_update_flags.h b/include/deal.II/fe/fe_update_flags.h index 5477216581..e094f8618e 100644 --- a/include/deal.II/fe/fe_update_flags.h +++ b/include/deal.II/fe/fe_update_flags.h @@ -355,9 +355,22 @@ namespace CellSimilarity { enum Similarity { + /** + * The cells differ by something besides a translation or inverted + * translations. + */ none, + /** + * The cells differ by a translation. + */ translation, + /** + * The cells differ by an inverted translation. + */ inverted_translation, + /** + * The next cell is not valid. + */ invalid_next_cell }; } diff --git a/include/deal.II/grid/grid_out.h b/include/deal.II/grid/grid_out.h index bb9c89dc32..2cd4a58d81 100644 --- a/include/deal.II/grid/grid_out.h +++ b/include/deal.II/grid/grid_out.h @@ -285,7 +285,14 @@ namespace GridOutFlags */ enum SizeType { - width, height + /** + * Scale with the width. + */ + width, + /** + * Scale with the height. + */ + height }; /** diff --git a/include/deal.II/grid/tria_accessor.h b/include/deal.II/grid/tria_accessor.h index 8bb37682bc..5f2f0a9e63 100644 --- a/include/deal.II/grid/tria_accessor.h +++ b/include/deal.II/grid/tria_accessor.h @@ -1860,8 +1860,17 @@ public: */ enum VertexKind { + /** + * Left vertex. + */ left_vertex, + /** + * Interior vertex. + */ interior_vertex, + /** + * Right vertex. + */ right_vertex }; diff --git a/include/deal.II/lac/arpack_solver.h b/include/deal.II/lac/arpack_solver.h index 887e358d67..4fba4b6ec9 100644 --- a/include/deal.II/lac/arpack_solver.h +++ b/include/deal.II/lac/arpack_solver.h @@ -118,14 +118,44 @@ public: */ enum WhichEigenvalues { + /** + * The algebraically largest eigenvalues. + */ algebraically_largest, + /** + * The algebraically smallest eigenvalues. + */ algebraically_smallest, + /** + * The eigenvalue with the largest magnitudes. + */ largest_magnitude, + /** + * The eigenvalue with the smallest magnitudes. + */ smallest_magnitude, + /** + * The eigenvalues with the largest real parts. + */ largest_real_part, + /** + * The eigenvalues with the smallest real parts. + */ smallest_real_part, + /** + * The eigenvalues with the largest imaginary parts. + */ largest_imaginary_part, + /** + * The eigenvalues with the smallest imaginary parts. + */ smallest_imaginary_part, + /** + * Compute half of the eigenvalues from the high end of the spectrum and + * the other half from the low end. If the number of requested + * eigenvectors is odd, then the extra eigenvector comes from the high end + * of the spectrum. + */ both_ends }; diff --git a/include/deal.II/lac/parpack_solver.h b/include/deal.II/lac/parpack_solver.h index fe49ef86a3..2280f686c8 100644 --- a/include/deal.II/lac/parpack_solver.h +++ b/include/deal.II/lac/parpack_solver.h @@ -160,14 +160,44 @@ public: */ enum WhichEigenvalues { + /** + * The algebraically largest eigenvalues. + */ algebraically_largest, + /** + * The algebraically smallest eigenvalues. + */ algebraically_smallest, + /** + * The eigenvalue with the largest magnitudes. + */ largest_magnitude, + /** + * The eigenvalue with the smallest magnitudes. + */ smallest_magnitude, + /** + * The eigenvalues with the largest real parts. + */ largest_real_part, + /** + * The eigenvalues with the smallest real parts. + */ smallest_real_part, + /** + * The eigenvalues with the largest imaginary parts. + */ largest_imaginary_part, + /** + * The eigenvalues with the smallest imaginary parts. + */ smallest_imaginary_part, + /** + * Compute half of the eigenvalues from the high end of the spectrum and + * the other half from the low end. If the number of requested + * eigenvectors is odd, then the extra eigenvector comes from the high end + * of the spectrum. + */ both_ends }; diff --git a/include/deal.II/lac/sparse_vanka.h b/include/deal.II/lac/sparse_vanka.h index 553bfec9a6..0b4319b2c7 100644 --- a/include/deal.II/lac/sparse_vanka.h +++ b/include/deal.II/lac/sparse_vanka.h @@ -517,7 +517,14 @@ public: */ enum BlockingStrategy { - index_intervals, adaptive + /** + * Block by index intervals. + */ + index_intervals, + /** + * Block with an adaptive strategy. + */ + adaptive }; /** diff --git a/include/deal.II/lac/trilinos_solver.h b/include/deal.II/lac/trilinos_solver.h index dc9902e608..44798f1a98 100644 --- a/include/deal.II/lac/trilinos_solver.h +++ b/include/deal.II/lac/trilinos_solver.h @@ -74,7 +74,29 @@ namespace TrilinosWrappers * one of the specialized derived classes when the solver should be set at * runtime. Currently enabled options are: */ - enum SolverName {cg, cgs, gmres, bicgstab, tfqmr} solver_name; + enum SolverName + { + /** + * Use the conjugate gradient (CG) algorithm. + */ + cg, + /** + * Use the conjugate gradient squared (CGS) algorithm. + */ + cgs, + /** + * Use the generalized minimum residual (GMRES) algorithm. + */ + gmres, + /** + * Use the biconjugate gradient stabilized (BICGStab) algorithm. + */ + bicgstab, + /** + * Use the transpose-free quasi-minimal residual (TFQMR) method. + */ + tfqmr + } solver_name; /** * Standardized data struct to pipe additional data to the solver. diff --git a/include/deal.II/lac/vector.h b/include/deal.II/lac/vector.h index 8a141b455d..c6af741f76 100644 --- a/include/deal.II/lac/vector.h +++ b/include/deal.II/lac/vector.h @@ -88,7 +88,21 @@ namespace parallel */ struct VectorOperation { - enum values { unknown, insert, add }; + enum values + { + /** + * The current operation is unknown. + */ + unknown, + /** + * The current operation is an insertion. + */ + insert, + /** + * The current operation is an addition. + */ + add + }; }; diff --git a/include/deal.II/matrix_free/fe_evaluation.h b/include/deal.II/matrix_free/fe_evaluation.h index cf88f73523..ca1e3af0ae 100644 --- a/include/deal.II/matrix_free/fe_evaluation.h +++ b/include/deal.II/matrix_free/fe_evaluation.h @@ -5013,8 +5013,22 @@ namespace internal */ enum EvaluatorVariant { + /** + * Do not use anything more than the tensor product structure of the + * finite element. + */ evaluate_general, + /** + * Perform evaluation by exploiting symmetry in the finite element: i.e., + * skip some computations by utilizing the symmetry in the shape functions + * and quadrature points. + */ evaluate_symmetric, + /** + * Use symmetry to apply the operator to even and odd parts of the input + * vector separately: see the documentation of the EvaluatorTensorProduct + * specialization for more information. + */ evaluate_evenodd }; diff --git a/include/deal.II/matrix_free/helper_functions.h b/include/deal.II/matrix_free/helper_functions.h index f6aeffbae8..7766dbbba1 100644 --- a/include/deal.II/matrix_free/helper_functions.h +++ b/include/deal.II/matrix_free/helper_functions.h @@ -135,7 +135,26 @@ namespace internal /** * Data type to identify cell type. */ - enum CellType {cartesian=0, affine=1, general=2, undefined=3}; + enum CellType + { + /** + * The cell is Cartesian. + */ + cartesian = 0, + /** + * The cell may be described with an affine mapping. + */ + affine = 1, + /** + * There is no special information available for compressing the + * representation of the cell. + */ + general = 2, + /** + * The cell type is undefined. + */ + undefined = 3 + }; /** diff --git a/include/deal.II/matrix_free/matrix_free.h b/include/deal.II/matrix_free/matrix_free.h index d2e1e9c4a1..b9515d756a 100644 --- a/include/deal.II/matrix_free/matrix_free.h +++ b/include/deal.II/matrix_free/matrix_free.h @@ -138,9 +138,30 @@ public: struct AdditionalData { /** - * Collects options for task parallelism. + * Collects options for task parallelism. See the documentation of the + * member variable MatrixFree::AdditionalData::tasks_parallel_scheme for a + * thorough description. */ - enum TasksParallelScheme {none, partition_partition, partition_color, color}; + enum TasksParallelScheme + { + /** + * Perform application in serial. + */ + none, + /** + * Partition the cells into two levels and afterwards form chunks. + */ + partition_partition, + /** + * Partition on the global level and color cells within the partitions. + */ + partition_color, + /** + * Use the traditional coloring algorithm: this is like + * TasksParallelScheme::partition_color, but only uses one partition. + */ + color + }; /** * Constructor for AdditionalData. @@ -879,7 +900,17 @@ private: std::vector > > dof_handler; std::vector > > hp_dof_handler; - enum ActiveDoFHandler { usual, hp } active_dof_handler; + enum ActiveDoFHandler + { + /** + * Use DoFHandler. + */ + usual, + /** + * Use hp::DoFHandler. + */ + hp + } active_dof_handler; unsigned int n_dof_handlers; unsigned int level; }; diff --git a/include/deal.II/meshworker/loop.h b/include/deal.II/meshworker/loop.h index 7f5a5df30f..f61db468b6 100644 --- a/include/deal.II/meshworker/loop.h +++ b/include/deal.II/meshworker/loop.h @@ -95,10 +95,24 @@ namespace MeshWorker */ bool ghost_cells; + /** + * Enumeration describing when to do assembly on a face: see, e.g., + * MeshWorker::LoopControl::faces_to_ghost for an example of how the value + * of this enumeration is interpreted in a particular circumstance. + */ enum FaceOption { + /** + * Do not perform assembly on a face. + */ never, + /** + * Perform assembly on one face. + */ one, + /** + * Perform assembly on both faces. + */ both }; diff --git a/include/deal.II/numerics/data_out_stack.h b/include/deal.II/numerics/data_out_stack.h index 4a1e071d06..d3fe8e8fb8 100644 --- a/include/deal.II/numerics/data_out_stack.h +++ b/include/deal.II/numerics/data_out_stack.h @@ -113,7 +113,17 @@ public: * Data type declaring the two types of vectors which are used in this * class. */ - enum VectorType { cell_vector, dof_vector }; + enum VectorType + { + /** + * The data describes one value for each cell. + */ + cell_vector, + /** + * The data describes one value for each DoF. + */ + dof_vector + }; /** * Destructor. Only declared to make it @p virtual. diff --git a/include/deal.II/numerics/histogram.h b/include/deal.II/numerics/histogram.h index fb8588633a..344643e18b 100644 --- a/include/deal.II/numerics/histogram.h +++ b/include/deal.II/numerics/histogram.h @@ -74,7 +74,14 @@ public: */ enum IntervalSpacing { - linear, logarithmic + /** + * Space intervals linearly. + */ + linear, + /** + * Space intervals logarithmically. + */ + logarithmic }; diff --git a/include/deal.II/numerics/solution_transfer.h b/include/deal.II/numerics/solution_transfer.h index 524e164148..ea494bad97 100644 --- a/include/deal.II/numerics/solution_transfer.h +++ b/include/deal.II/numerics/solution_transfer.h @@ -433,7 +433,18 @@ private: */ enum PreparationState { - none, pure_refinement, coarsening_and_refinement + /** + * The SolutionTransfer is not yet prepared. + */ + none, + /** + * The SolutionTransfer is prepared for purely refinement. + */ + pure_refinement, + /** + * The SolutionTransfer is prepared for coarsening and refinement. + */ + coarsening_and_refinement }; /** diff --git a/include/deal.II/numerics/time_dependent.h b/include/deal.II/numerics/time_dependent.h index 660e2cede6..912cd19d78 100644 --- a/include/deal.II/numerics/time_dependent.h +++ b/include/deal.II/numerics/time_dependent.h @@ -408,7 +408,14 @@ public: */ enum Direction { - forward, backward + /** + * Go in the forward direction. + */ + forward, + /** + * Go in the backward direction. + */ + backward }; /** @@ -643,8 +650,17 @@ public: */ enum SolutionState { + /** + * Solve the primal problem next. + */ primal_problem = 0x0, + /** + * Solve the dual problem next. + */ dual_problem = 0x1, + /** + * Perform postprocessing next. + */ postprocess = 0x2 }; @@ -1272,6 +1288,9 @@ public: */ enum SolutionState { + /** + * Perform grid refinement next. + */ grid_refinement = 0x1000 };