Reference documentation for deal.II version GIT relicensing-214-g6e74dec06b 2024-03-27 18:10:01+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
PointValueHistory< dim > Class Template Reference

#include <deal.II/numerics/point_value_history.h>

Public Member Functions

 PointValueHistory (const unsigned int n_independent_variables=0)
 
 PointValueHistory (const DoFHandler< dim > &dof_handler, const unsigned int n_independent_variables=0)
 
 PointValueHistory (const PointValueHistory &point_value_history)
 
PointValueHistoryoperator= (const PointValueHistory &point_value_history)
 
 ~PointValueHistory ()
 
void add_point (const Point< dim > &location)
 
void add_points (const std::vector< Point< dim > > &locations)
 
void add_field_name (const std::string &vector_name, const ComponentMask &component_mask={})
 
void add_field_name (const std::string &vector_name, const unsigned int n_components)
 
void add_component_names (const std::string &vector_name, const std::vector< std::string > &component_names)
 
void add_independent_names (const std::vector< std::string > &independent_names)
 
template<typename VectorType >
void evaluate_field (const std::string &name, const VectorType &solution)
 
template<typename VectorType >
void evaluate_field (const std::vector< std::string > &names, const VectorType &solution, const DataPostprocessor< dim > &data_postprocessor, const Quadrature< dim > &quadrature)
 
template<typename VectorType >
void evaluate_field (const std::string &name, const VectorType &solution, const DataPostprocessor< dim > &data_postprocessor, const Quadrature< dim > &quadrature)
 
template<typename VectorType >
void evaluate_field_at_requested_location (const std::string &name, const VectorType &solution)
 
void start_new_dataset (const double key)
 
void push_back_independent (const std::vector< double > &independent_values)
 
void write_gnuplot (const std::string &base_name, const std::vector< Point< dim > > &postprocessor_locations=std::vector< Point< dim > >())
 
Vector< double > mark_support_locations ()
 
void get_support_locations (std::vector< std::vector< Point< dim > > > &locations)
 
void get_postprocessor_locations (const Quadrature< dim > &quadrature, std::vector< Point< dim > > &locations)
 
void close ()
 
void clear ()
 
void status (std::ostream &out)
 
bool deep_check (const bool strict)
 

Static Public Member Functions

static ::ExceptionBaseExcNoIndependent ()
 
static ::ExceptionBaseExcDataLostSync ()
 
static ::ExceptionBaseExcDoFHandlerRequired ()
 
static ::ExceptionBaseExcDoFHandlerChanged ()
 

Private Member Functions

void tria_change_listener ()
 

Private Attributes

std::vector< double > dataset_key
 
std::vector< std::vector< double > > independent_values
 
std::vector< std::string > indep_names
 
std::map< std::string, std::vector< std::vector< double > > > data_store
 
std::map< std::string, ComponentMaskcomponent_mask
 
std::map< std::string, std::vector< std::string > > component_names_map
 
std::vector< internal::PointValueHistoryImplementation::PointGeometryData< dim > > point_geometry_data
 
bool closed
 
bool cleared
 
SmartPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
 
bool triangulation_changed
 
bool have_dof_handler
 
boost::signals2::connection tria_listener
 
unsigned int n_indep
 

Detailed Description

template<int dim>
class PointValueHistory< dim >

PointValueHistory tackles the overhead of plotting time (or any other iterative process) graphs of solution values at specific points on the mesh. The user specifies the points which the solution should be monitored at ahead of time, as well as giving each solution vector that they want to record a mnemonic name. Then, for each step the user calls one of the three available "evaluate field" methods to store the data from each time step, and the class extracts data for the requested points to store it. Finally, once the computation is finished, the user can request output files to be generated; these files are in Gnuplot format but are basically just regular text and can easily be imported into other programs well, for example into spreadsheets.

The user can store extra variables which do not relate to mesh location specifying n_independent_variables. The class then expects a std::vector of size n_independent_variables to be added during each step using the method push_back_independent. This may be used for example for recording external input, logging solver performance data such as time taken to solve the step and solver steps before convergence, saving norms calculated, or simply saving the time, number of time step, or number of nonlinear iteration along with the data evaluated from the mesh.

The three "evaluate field" methods each have different strengths and weaknesses making each suitable for different contexts:

  1. Firstly, the evaluate_field version that does not take a DataPostprocessor object selects the nearest support point (see this entry in the glossary ) to a given point to extract data from. This makes the code that needs to be run at each time step very short, since looping over the mesh to extract the needed dof_index can be done just once at the start. However, this method is not suitable for FiniteElement objects that do not assign dofs to actual mesh locations (i.e. FEs without support points ) or if adaptive mesh refinement is used. The reason for the latter restriction is that the location of the closest support point to a given point may change upon mesh refinement. The class will throw an exception if any change to the triangulation is made (Although the nearest support point could be re- computed upon mesh refinement, the location of the support point will most likely change slightly, making the interpretation of the data difficult, hence this is not implemented currently.)

  2. Secondly, evaluate_field_at_requested_location calls VectorTools::point_value to compute values at the specific point requested. This method is valid for any FE that is supported by VectorTools::point_value. Specifically, this method can be called by codes using adaptive mesh refinement.

  3. Finally, the class offers a function evaluate_field that takes a DataPostprocessor object. This method allows the deal.II data postprocessor to be used to compute new quantities from the solution on the fly. The values are located at the nearest quadrature point to the requested point. If the mesh is refined between calls, this point will change, so care must be taken when using this method in code using adaptive refinement, but as the output will be meaningful (in the sense that the quadrature point selected is guaranteed to remain in the same vicinity, the class does not prevent the use of this method in adaptive codes. The class provides warnings in the output files if the mesh has changed. Note that one can reduce the error this procedure introduces by providing a quadrature formula that has more points, at the expense of performing more work since then the closest quadrature points is nearer to the point at which the evaluation is really supposed to happen. (As a sidenote: Why not do the evaluation at the requested point right away? The reason for this is that it would require setting up a new quadrature point object on each cell that has only a single point corresponding to the reference coordinates of the point you really want; then initializing a FEValues object with it; then evaluating the solution at this point; then handing the result to the DataPostprocessor object. This sequence of things is expensive – which is the reason why VectorTools::point_value() is expensive. Using the same quadrature formula on each cell on which we want to evaluate the solution and only having to initialize a FEValue object once is a much cheaper alternative, albeit of course at the expense of getting only an approximate result.)

When recording a new mnemonic name, the user must supply a component_mask (see this glossary entry ) to indicate the (vector) components to be extracted from the given input. If the user simply wants to extract all the components, the mask need not be explicitly supplied to the add_field_name method and the default value of the parameter is sufficient. If the evaluate_field with a DataPostprocessor object is used, the component_mask is interpreted as the mask of the DataPostprocessor return vector. The size of this mask can be different to that of the FE space, but must be provided when the add_field_name method is called. One variant of the add_field_name method allows an unsigned int input to construct a suitable mask, if all values from the DataPostprocessor are desired.

The class automatically generates names for the data stored based on the mnemonics supplied. The methods add_component_names and add_independent_names allow the user to provide lists of names to use instead if desired.

Following is a little code snippet that shows a common usage of this class:

//....
//... code to set up Triangulation, perform any refinement necessary
// and set up DoFHandler, sizing solution Vectors etc
// just one independent value, which happens to be an input
unsigned int n_inputs = 1;
// call the constructor
PointValueHistory<dim> node_monitor(dof_handler, n_inputs);
// set up fields and points required
node_monitor.add_field_name("Solution");
std::vector <Point <dim> > point_vector(2);
point_vector[0] = Point <dim>(0, 0);
point_vector[1] = Point <dim>(0.25, 0);
node_monitor.add_points(point_vector); // multiple points at once
node_monitor.add_point(Point<dim>(1, 0.2)); // add a single point
node_monitor.close(); // close the class once the setup is complete
node_monitor.status(std::cout); // print out status to check if desired
// ... more code ...
// ... in an iterative loop ...
// double time, vector <double> with size 1 input_value,
// and Vector <double> solution calculated in the loop
node_monitor.start_new_dataset(time);
node_monitor.push_back_independent(input_value);
node_monitor.evaluate_field("Solution", solution);
// ... end of iterative loop ...
node_monitor.write_gnuplot("node"); // write out data files
SmartPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
Definition point.h:111

Definition at line 212 of file point_value_history.h.

Constructor & Destructor Documentation

◆ PointValueHistory() [1/3]

template<int dim>
PointValueHistory< dim >::PointValueHistory ( const unsigned int  n_independent_variables = 0)

Provide a stripped down instance of the class which does not support adding points or mesh data. This may be used for example for recording external input or logging solver performance data.

Definition at line 56 of file point_value_history.cc.

◆ PointValueHistory() [2/3]

template<int dim>
PointValueHistory< dim >::PointValueHistory ( const DoFHandler< dim > &  dof_handler,
const unsigned int  n_independent_variables = 0 
)

Constructor linking the class to a specific DoFHandler. This class reads specific data from the DoFHandler and stores it internally for quick access (in particular dof indices of closest neighbors to requested points) the class is fairly intolerant to changes to the DoFHandler if data at support points is required. Mesh refinement and DoFRenumbering methods should be performed before the add_points method is called and adaptive grid refinement is only supported by some methods.

The user can store extra variables which do not relate to mesh location by specifying the number required using n_independent_variables and making calls to push_back_independent as needed. This may be used for example for recording external input or logging solver performance data.

Definition at line 77 of file point_value_history.cc.

◆ PointValueHistory() [3/3]

template<int dim>
PointValueHistory< dim >::PointValueHistory ( const PointValueHistory< dim > &  point_value_history)

Copy constructor. This constructor can be safely called with a PointValueHistory object that contains data, but this could be expensive and should be avoided.

Definition at line 103 of file point_value_history.cc.

◆ ~PointValueHistory()

template<int dim>
PointValueHistory< dim >::~PointValueHistory ( )

Deconstructor.

Definition at line 171 of file point_value_history.cc.

Member Function Documentation

◆ operator=()

template<int dim>
PointValueHistory< dim > & PointValueHistory< dim >::operator= ( const PointValueHistory< dim > &  point_value_history)

Assignment operator. This assignment operator can be safely called once the class is closed and data added, but this is provided primarily to allow a PointValueHistory object declared in a class to be reinitialized later in the class. Using the assignment operator when the object contains data could be expensive.

Definition at line 137 of file point_value_history.cc.

◆ add_point()

template<int dim>
void PointValueHistory< dim >::add_point ( const Point< dim > &  location)

Add a single point to the class. The support points (one per component) in the mesh that are closest to that point are found and their details stored for use when evaluate_field is called. If more than one point is required rather use the add_points method since this minimizes iterations over the mesh.

Definition at line 183 of file point_value_history.cc.

◆ add_points()

template<int dim>
void PointValueHistory< dim >::add_points ( const std::vector< Point< dim > > &  locations)

Add multiple points to the class. The support points (one per component) in the mesh that are closest to that point is found and their details stored for use when evaluate_field is called. If more than one point is required, rather call this method as it is more efficient than the add_point method since it minimizes iterations over the mesh. The points are added to the internal database in the order they appear in the list and there is always a one to one correspondence between the requested point and the added point, even if a point is requested multiple times.

Definition at line 319 of file point_value_history.cc.

◆ add_field_name() [1/2]

template<int dim>
void PointValueHistory< dim >::add_field_name ( const std::string &  vector_name,
const ComponentMask component_mask = {} 
)

Put another mnemonic string (and hence VectorType) into the class. This method adds storage space for variables equal to the number of true values in component_mask. This also adds extra entries for points that are already in the class, so add_field_name and add_points can be called in any order.

Definition at line 456 of file point_value_history.cc.

◆ add_field_name() [2/2]

template<int dim>
void PointValueHistory< dim >::add_field_name ( const std::string &  vector_name,
const unsigned int  n_components 
)

Put another mnemonic string (and hence VectorType) into the class. This method adds storage space for n_components variables. This also adds extra entries for points that are already in the class, so add_field_name and add_points can be called in any order. This method generates a std::vector 0, ..., n_components-1 and calls the previous function.

Definition at line 502 of file point_value_history.cc.

◆ add_component_names()

template<int dim>
void PointValueHistory< dim >::add_component_names ( const std::string &  vector_name,
const std::vector< std::string > &  component_names 
)

Provide optional names for each component of a field. These names will be used instead of names generated from the field name, if supplied.

Definition at line 512 of file point_value_history.cc.

◆ add_independent_names()

template<int dim>
void PointValueHistory< dim >::add_independent_names ( const std::vector< std::string > &  independent_names)

Provide optional names for the independent values. These names will be used instead of "Indep_...", if supplied.

Definition at line 535 of file point_value_history.cc.

◆ evaluate_field() [1/3]

template<int dim>
template<typename VectorType >
void PointValueHistory< dim >::evaluate_field ( const std::string &  name,
const VectorType &  solution 
)

Extract values at the stored points from the VectorType supplied and add them to the new dataset in vector_name. The component mask supplied when the field was added is used to select components to extract. If a DoFHandler is used, one (and only one) evaluate_field method must be called for each dataset (time step, iteration, etc) for each vector_name, otherwise a ExcDataLostSync error can occur.

Definition at line 578 of file point_value_history.cc.

◆ evaluate_field() [2/3]

template<int dim>
template<typename VectorType >
void PointValueHistory< dim >::evaluate_field ( const std::vector< std::string > &  names,
const VectorType &  solution,
const DataPostprocessor< dim > &  data_postprocessor,
const Quadrature< dim > &  quadrature 
)

Compute values using a DataPostprocessor object with the VectorType supplied and add them to the new dataset in vector_name. The component_mask supplied when the field was added is used to select components to extract from the DataPostprocessor return vector. This method takes a vector of field names to process and is preferred if many fields use the same DataPostprocessor object as each cell is only located once. The quadrature object supplied is used for all components of a vector field. Although this method will not throw an exception if the mesh has changed. (No internal data structures are invalidated as the quadrature points are repicked each time the function is called.) Nevertheless the user must be aware that if the mesh changes the point selected will also vary slightly, making interpretation of the data more difficult. If a DoFHandler is used, one (and only one) evaluate_field method must be called for each dataset (time step, iteration, etc) for each vector_name, otherwise a ExcDataLostSync error can occur.

Definition at line 646 of file point_value_history.cc.

◆ evaluate_field() [3/3]

template<int dim>
template<typename VectorType >
void PointValueHistory< dim >::evaluate_field ( const std::string &  name,
const VectorType &  solution,
const DataPostprocessor< dim > &  data_postprocessor,
const Quadrature< dim > &  quadrature 
)

Construct a std::vector <std::string> containing only vector_name and call the above function. The above function is more efficient if multiple fields use the same DataPostprocessor object.

Definition at line 878 of file point_value_history.cc.

◆ evaluate_field_at_requested_location()

template<int dim>
template<typename VectorType >
void PointValueHistory< dim >::evaluate_field_at_requested_location ( const std::string &  name,
const VectorType &  solution 
)

Extract values at the points actually requested from the VectorType supplied and add them to the new dataset in vector_name. Unlike the other evaluate_field methods this method does not care if the dof_handler has been modified because it uses calls to VectorTools::point_value to extract there data. Therefore, if only this method is used, the class is fully compatible with adaptive refinement. The component_mask supplied when the field was added is used to select components to extract. If a DoFHandler is used, one (and only one) evaluate_field method must be called for each dataset (time step, iteration, etc) for each vector_name, otherwise a ExcDataLostSync error can occur.

Definition at line 894 of file point_value_history.cc.

◆ start_new_dataset()

template<int dim>
void PointValueHistory< dim >::start_new_dataset ( const double  key)

Add the key for the current dataset to the dataset. Although calling this method first is sensible, the order in which this method, evaluate_field and push_back_independent is not important. It is however important that all the data for a give dataset is added to each dataset and that it is added before a new data set is started. This prevents a ExcDataLostSync.

Definition at line 963 of file point_value_history.cc.

◆ push_back_independent()

template<int dim>
void PointValueHistory< dim >::push_back_independent ( const std::vector< double > &  independent_values)

If independent values have been set up, this method stores these values. This should only be called once per dataset, and if independent values are used it must be called for every dataset. A ExcDataLostSync exception can be thrown if this method is not called.

Definition at line 978 of file point_value_history.cc.

◆ write_gnuplot()

template<int dim>
void PointValueHistory< dim >::write_gnuplot ( const std::string &  base_name,
const std::vector< Point< dim > > &  postprocessor_locations = std::vector<Point<dim>>() 
)

Write out a series of .gpl files named base_name + "-00.gpl", base_name + "-01.gpl" etc. The data file gives information about where the support points selected and interpreting the data. If n_indep != 0 an additional file base_name + "_indep.gpl" containing key and independent data. The file name agrees with the order the points were added to the class. The names of the data columns can be supplied using the functions add_component_names and add_independent_names. The support point information is only meaningful if the dof_handler has not been changed. Therefore, if adaptive mesh refinement has been used the support point data should not be used. The optional parameter postprocessor_locations is used to add the postprocessor locations to the output files. If this is desired, the data should be obtained from a call to get_postprocessor_locations while the dof_handler is usable. The default parameter is an empty vector of strings, and will suppress postprocessor locations output.

Definition at line 1000 of file point_value_history.cc.

◆ mark_support_locations()

template<int dim>
Vector< double > PointValueHistory< dim >::mark_support_locations ( )

Return a Vector with the indices of selected points flagged with a 1. This method is mainly for testing and verifying that the class is working correctly. By passing this vector to a DataOut object, the user can verify that the positions returned by get_support_locations agree with the positions that DataOut interprets from the Vector returned. The code snippet below demonstrates how this could be done:

// Make a DataOut object and attach the dof_handler
DataOut<dim> data_out;
// Call the mark_locations method to get the vector with indices flagged
Vector<double> support_point_locations = node_monitor.mark_locations();
// Add the vector to the data_out object and
// write out a file in the usual way
data_out.add_data_vector(support_point_locations, "Monitor_Locations");
data_out.build_patches(2);
std::ofstream output("locations.gpl");
data_out.write_gnuplot(output);
void write_gnuplot(std::ostream &out) const
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition data_out.cc:1062

Definition at line 1204 of file point_value_history.cc.

◆ get_support_locations()

template<int dim>
void PointValueHistory< dim >::get_support_locations ( std::vector< std::vector< Point< dim > > > &  locations)

Stores the actual location of each support point selected by the add_point(s) method. This can be used to compare with the point requested, for example by using the Point<dim>::distance function. For convenience, location is resized to the correct number of points by the method.

Definition at line 1232 of file point_value_history.cc.

◆ get_postprocessor_locations()

template<int dim>
void PointValueHistory< dim >::get_postprocessor_locations ( const Quadrature< dim > &  quadrature,
std::vector< Point< dim > > &  locations 
)

Stores the actual location of the points used by the data_postprocessor. This can be used to compare with the points requested, for example by using the Point<dim>::distance function. Unlike the support_locations, these locations are computed every time the evaluate_field method is called with a postprocessor. This method uses the same algorithm so can will find the same points. For convenience, location is resized to the correct number of points by the method.

Definition at line 1255 of file point_value_history.cc.

◆ close()

template<int dim>
void PointValueHistory< dim >::close ( )

Once datasets have been added to the class, requests to add additional points will make the data interpretation unclear. The boolean closed defines a state of the class and ensures this does not happen. Additional points or vectors can only be added while the class is not closed, and the class must be closed before datasets can be added or written to file. PointValueHistory::get_support_locations and PointValueHistory::status do not require the class to be closed. If a method that requires a class to be open or close is called while in the wrong state a ExcInvalidState exception is thrown.

Definition at line 547 of file point_value_history.cc.

◆ clear()

template<int dim>
void PointValueHistory< dim >::clear ( )

Delete the lock this object has to the DoFHandler used the last time the class was created. This method should not normally need to be called, but can be useful to ensure that the DoFHandler is released before it goes out of scope if the PointValueHistory class might live longer than it. Once this method has been called, the majority of methods will throw a ExcInvalidState exception, so if used this method should be the last call to the class.

Definition at line 556 of file point_value_history.cc.

◆ status()

template<int dim>
void PointValueHistory< dim >::status ( std::ostream &  out)

Print useful debugging information about the class, include details about which support points were selected for each point and sizes of the data stored.

Definition at line 1314 of file point_value_history.cc.

◆ deep_check()

template<int dim>
bool PointValueHistory< dim >::deep_check ( const bool  strict)

Check the internal data sizes to test for a loss of data sync. This is often used in Assert statements with the ExcDataLostSync exception. If strict is false this method returns true if all sizes are within 1 of each other (needed to allow data to be added), with strict = true they must be exactly equal.

Definition at line 1425 of file point_value_history.cc.

◆ tria_change_listener()

template<int dim>
void PointValueHistory< dim >::tria_change_listener ( )
private

A function that will be triggered through signals whenever the triangulation is modified.

It is currently used to check if the triangulation has changed, invalidating precomputed values.

Definition at line 1487 of file point_value_history.cc.

Member Data Documentation

◆ dataset_key

template<int dim>
std::vector<double> PointValueHistory< dim >::dataset_key
private

Stores keys, values on the abscissa. This will often be time, but possibly time step, iteration etc.

Definition at line 580 of file point_value_history.h.

◆ independent_values

template<int dim>
std::vector<std::vector<double> > PointValueHistory< dim >::independent_values
private

Values that do not depend on grid location.

Definition at line 585 of file point_value_history.h.

◆ indep_names

template<int dim>
std::vector<std::string> PointValueHistory< dim >::indep_names
private

Save a vector listing component names associated with a independent_values. This will be an empty vector if the user does not supplies names.

Definition at line 592 of file point_value_history.h.

◆ data_store

template<int dim>
std::map<std::string, std::vector<std::vector<double> > > PointValueHistory< dim >::data_store
private

Save data for each mnemonic entry. data_store: mnemonic -> [point_0_components point_1_components ... point_n-1_components][key] This format facilitates scalar mnemonics in a vector space, because scalar mnemonics will only have one component per point. Vector components are strictly FE.n_components () long.

Definition at line 601 of file point_value_history.h.

◆ component_mask

template<int dim>
std::map<std::string, ComponentMask> PointValueHistory< dim >::component_mask
private

Save a component mask for each mnemonic.

Definition at line 606 of file point_value_history.h.

◆ component_names_map

template<int dim>
std::map<std::string, std::vector<std::string> > PointValueHistory< dim >::component_names_map
private

Save a vector listing component names associated with a mnemonic. This will be an empty vector if the user does not supplies names.

Definition at line 613 of file point_value_history.h.

◆ point_geometry_data

template<int dim>
std::vector<internal::PointValueHistoryImplementation::PointGeometryData<dim> > PointValueHistory< dim >::point_geometry_data
private

Save the location and other mesh information about support points.

Definition at line 619 of file point_value_history.h.

◆ closed

template<int dim>
bool PointValueHistory< dim >::closed
private

Used to enforce closed state for some methods.

Definition at line 625 of file point_value_history.h.

◆ cleared

template<int dim>
bool PointValueHistory< dim >::cleared
private

Used to enforce !cleared state for some methods.

Definition at line 630 of file point_value_history.h.

◆ dof_handler

template<int dim>
SmartPointer<const DoFHandler<dim>, PointValueHistory<dim> > PointValueHistory< dim >::dof_handler
private

A smart pointer to the dof_handler supplied to the constructor. This can be released by calling clear().

Definition at line 637 of file point_value_history.h.

◆ triangulation_changed

template<int dim>
bool PointValueHistory< dim >::triangulation_changed
private

Variable to check if the triangulation has changed. If it has changed, certain data is out of date (especially the PointGeometryData::solution_indices.

Definition at line 645 of file point_value_history.h.

◆ have_dof_handler

template<int dim>
bool PointValueHistory< dim >::have_dof_handler
private

A boolean to record whether the class was initialized with a DoFHandler or not.

Definition at line 651 of file point_value_history.h.

◆ tria_listener

template<int dim>
boost::signals2::connection PointValueHistory< dim >::tria_listener
private

Used to detect signals from the Triangulation.

Definition at line 656 of file point_value_history.h.

◆ n_indep

template<int dim>
unsigned int PointValueHistory< dim >::n_indep
private

Stores the number of independent variables requested.

Definition at line 661 of file point_value_history.h.


The documentation for this class was generated from the following files: