Reference documentation for deal.II version 9.1.1
|
#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) | |
PointValueHistory & | operator= (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=ComponentMask()) |
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<class VectorType > | |
void | evaluate_field (const std::string &name, const VectorType &solution) |
template<class VectorType > | |
void | evaluate_field (const std::vector< std::string > &names, const VectorType &solution, const DataPostprocessor< dim > &data_postprocessor, const Quadrature< dim > &quadrature) |
template<class VectorType > | |
void | evaluate_field (const std::string &name, const VectorType &solution, const DataPostprocessor< dim > &data_postprocessor, const Quadrature< dim > &quadrature) |
template<class 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_points (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 ::ExceptionBase & | ExcNoIndependent () |
static ::ExceptionBase & | ExcDataLostSync () |
static ::ExceptionBase & | ExcDoFHandlerRequired () |
static ::ExceptionBase & | ExcDoFHandlerChanged () |
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, ComponentMask > | component_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 |
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:
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.)
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.
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:
Definition at line 212 of file point_value_history.h.
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 58 of file point_value_history.cc.
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 79 of file point_value_history.cc.
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 105 of file point_value_history.cc.
PointValueHistory< dim >::~PointValueHistory | ( | ) |
Deconstructor.
Definition at line 175 of file point_value_history.cc.
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 140 of file point_value_history.cc.
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 187 of file point_value_history.cc.
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 327 of file point_value_history.cc.
void PointValueHistory< dim >::add_field_name | ( | const std::string & | vector_name, |
const ComponentMask & | component_mask = ComponentMask() |
||
) |
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 468 of file point_value_history.cc.
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 514 of file point_value_history.cc.
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 524 of file point_value_history.cc.
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 547 of file point_value_history.cc.
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 590 of file point_value_history.cc.
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 658 of file point_value_history.cc.
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 885 of file point_value_history.cc.
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 901 of file point_value_history.cc.
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 970 of file point_value_history.cc.
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 985 of file point_value_history.cc.
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 1007 of file point_value_history.cc.
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_points
agree with the positions that DataOut
interprets from the Vector
returned. The code snippet below demonstrates how this could be done:
Definition at line 1219 of file point_value_history.cc.
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 1247 of file point_value_history.cc.
void PointValueHistory< dim >::get_points | ( | std::vector< std::vector< Point< dim >>> & | locations | ) |
This function only exists for backward compatibility as this is the interface provided by previous versions of the library. The function get_support_locations replaces it and reflects the fact that the points returned are actually the support points.
Definition at line 1269 of file point_value_history.cc.
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 1278 of file point_value_history.cc.
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_points
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 559 of file point_value_history.cc.
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 568 of file point_value_history.cc.
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 1332 of file point_value_history.cc.
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 1449 of file point_value_history.cc.
|
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 1515 of file point_value_history.cc.
|
private |
Stores keys, values on the abscissa. This will often be time, but possibly time step, iteration etc.
Definition at line 592 of file point_value_history.h.
|
private |
Values that do not depend on grid location.
Definition at line 597 of file point_value_history.h.
|
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 604 of file point_value_history.h.
|
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 613 of file point_value_history.h.
|
private |
Save a component mask for each mnemonic.
Definition at line 618 of file point_value_history.h.
|
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 625 of file point_value_history.h.
|
private |
Save the location and other mesh information about support points.
Definition at line 631 of file point_value_history.h.
|
private |
Used to enforce closed
state for some methods.
Definition at line 637 of file point_value_history.h.
|
private |
Used to enforce !cleared
state for some methods.
Definition at line 642 of file point_value_history.h.
|
private |
A smart pointer to the dof_handler supplied to the constructor. This can be released by calling clear()
.
Definition at line 649 of file point_value_history.h.
|
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 657 of file point_value_history.h.
|
private |
A boolean to record whether the class was initialized with a DoFHandler or not.
Definition at line 663 of file point_value_history.h.
|
private |
Used to detect signals from the Triangulation.
Definition at line 668 of file point_value_history.h.
|
private |
Stores the number of independent variables requested.
Definition at line 673 of file point_value_history.h.