deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40:00+00:00
|
#include <deal.II/base/table.h>
Public Types | |
using | value_type = T |
using | size_type = typename AlignedVector< T >::size_type |
Public Member Functions | |
TableBase ()=default | |
TableBase (const TableIndices< N > &sizes) | |
template<typename InputIterator > | |
TableBase (const TableIndices< N > &sizes, InputIterator entries, const bool C_style_indexing=true) | |
TableBase (const TableBase< N, T > &src) | |
template<typename T2 > | |
TableBase (const TableBase< N, T2 > &src) | |
TableBase (TableBase< N, T > &&src) noexcept | |
~TableBase () override=default | |
TableBase< N, T > & | operator= (const TableBase< N, T > &src) |
template<typename T2 > | |
TableBase< N, T > & | operator= (const TableBase< N, T2 > &src) |
TableBase< N, T > & | operator= (TableBase< N, T > &&src) noexcept |
bool | operator== (const TableBase< N, T > &T2) const |
void | reset_values () |
void | reinit (const TableIndices< N > &new_size, const bool omit_default_initialization=false) |
void | clear () |
size_type | size (const unsigned int i) const |
const TableIndices< N > & | size () const |
size_type | n_elements () const |
bool | empty () const |
template<typename InputIterator > | |
void | fill (InputIterator entries, const bool C_style_indexing=true) |
void | fill (const T &value) |
AlignedVector< T >::reference | operator() (const TableIndices< N > &indices) |
AlignedVector< T >::const_reference | operator() (const TableIndices< N > &indices) const |
void | replicate_across_communicator (const MPI_Comm communicator, const unsigned int root_process) |
void | swap (TableBase< N, T > &v) noexcept |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
EnableObserverPointer functionality | |
Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer class. | |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
size_type | position (const TableIndices< N > &indices) const |
AlignedVector< T >::reference | el (const TableIndices< N > &indices) |
AlignedVector< T >::const_reference | el (const TableIndices< N > &indices) const |
Protected Attributes | |
AlignedVector< T > | values |
TableIndices< N > | table_size |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | check_no_subscribers () const noexcept |
Private Attributes | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
Friends | |
template<int , typename > | |
class | TableBase |
A class holding a multi-dimensional array of objects of templated type. If the template parameter indicating the number of dimensions is one, then this class more or less represents a vector; if it is two then it is a matrix; and so on.
This class specifically replaces attempts at higher-dimensional arrays like std::vector<std::vector<T>>
, or even higher nested constructs. These constructs have the disadvantage that they are hard to initialize, and most importantly that they are very inefficient if all rows of a matrix or higher-dimensional table have the same size (which is the usual case), since then the memory for each row is allocated independently, both wasting time and memory. This can be made more efficient by allocating only one chunk of memory for the entire object, which is what the current class does.
In some way, this class is similar to the Tensor class, in that it templatizes on the number of dimensions. However, there are two major differences. The first is that the Tensor class stores only numeric values (as double
s), while the Table class stores arbitrary objects. The second is that the Tensor class has fixed sizes in each dimension, also given as a template argument, while this class can handle arbitrary and different sizes in each dimension.
This has two consequences. First, since the size is not known at compile time, it has to do explicit memory allocation. Second, the layout of individual elements is not known at compile time, so access is slower than for the Tensor class where the number of elements are their location is known at compile time and the compiler can optimize with this knowledge (for example when unrolling loops). On the other hand, this class is of course more flexible, for example when you want a two-dimensional table with the number of rows equal to the number of degrees of freedom on a cell, and the number of columns equal to the number of quadrature points. Both numbers may only be known at run-time, so a flexible table is needed here. Furthermore, you may want to store, say, the gradients of shape functions, so the data type is not a single scalar value, but a tensor itself.
The Table classes (derived from this class) are frequently used to store large data tables. A modest example is given in step-53 where we store a \(380 \times 220\) table of geographic elevation data for a region of Africa, and this data requires about 670 kB if memory; however, tables that store three- or more-dimensional data (say, information about the density, pressure, and temperature in the earth interior on a regular grid of (latitude, longitude, depth)
points) can easily run into hundreds of megabytes or more. These tables are then often provided to classes such as InterpolatedTensorProductGridData or InterpolatedUniformGridData.
If you need to load such tables on single-processor (or multi-threaded) jobs, then there is nothing you can do about the size of these tables: The table just has to fit into memory. But, if your program is parallelized via MPI, then a typical first implementation would create a table object on every process and fill it on every MPI process by reading the data from a file. This is inefficient from two perspectives:
Both of these use cases are enabled by the TableBase::replicate_across_communicator() function that is internally based on AlignedVector::replicate_across_communicator(). This function allows for workflows like the following where we put that MPI process with rank zero in charge of reading the data (but it could have been any other "root rank" as well):
The last call in this code snippet makes sure that the data is made available on all non-root processes, either by re-creating a copy of the table in the other processes' memory space or, if possible, by creating copies in shared memory once for all processes located on each of the machines used by the MPI job.
using TableBase< N, T >::size_type = typename AlignedVector<T>::size_type |
|
privateinherited |
The data type used in counter_map.
Definition at line 238 of file enable_observer_pointer.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 243 of file enable_observer_pointer.h.
|
default |
Default constructor. Set all dimensions to zero.
|
explicit |
Constructor. Initialize the array with the given dimensions in each index component.
TableBase< N, T >::TableBase | ( | const TableIndices< N > & | sizes, |
InputIterator | entries, | ||
const bool | C_style_indexing = true |
||
) |
Constructor. Initialize the array with the given dimensions in each index component, and then initialize the elements of the table using the second and third argument by calling fill(entries,C_style_indexing).
TableBase< N, T >::TableBase | ( | const TableBase< N, T > & | src | ) |
Copy constructor. Performs a deep copy.
TableBase< N, T >::TableBase | ( | const TableBase< N, T2 > & | src | ) |
Copy constructor. Performs a deep copy from a table object storing some other data type.
|
noexcept |
Move constructor. Transfers the contents of another Table.
Destructor. Free allocated memory.
TableBase< N, T > & TableBase< N, T >::operator= | ( | const TableBase< N, T > & | src | ) |
Assignment operator. Copy all elements of src
into the matrix. The size is adjusted if needed.
We can't use the other, templatized version since if we don't declare this one, the compiler will happily generate a predefined copy operator which is not what we want.
TableBase< N, T > & TableBase< N, T >::operator= | ( | const TableBase< N, T2 > & | src | ) |
Copy operator. Copy all elements of src
into the array. The size is adjusted if needed.
This function requires that the type T2
is convertible to T
.
|
noexcept |
Move assignment operator. Transfer all elements of src
into the table.
bool TableBase< N, T >::operator== | ( | const TableBase< N, T > & | T2 | ) | const |
Test for equality of two tables.
Set all entries to their default value (i.e. copy them over with default constructed objects). Do not change the size of the table, though.
void TableBase< N, T >::reinit | ( | const TableIndices< N > & | new_size, |
const bool | omit_default_initialization = false |
||
) |
Set the dimensions of this object to the sizes given in the first argument, and allocate the required memory for table entries to accommodate these sizes. If omit_default_initialization
is set to false
, all elements of the table are set to a default constructed object for the element type. Otherwise the memory is left in an uninitialized or otherwise undefined state.
Size of the table in direction i
.
const TableIndices< N > & TableBase< N, T >::size | ( | ) | const |
Return the sizes of this object in each direction.
Return the number of elements stored in this object, which is the product of the extensions in each dimension.
Return whether the object is empty, i.e. one of the directions is zero. This is equivalent to n_elements()==0
.
void TableBase< N, T >::fill | ( | InputIterator | entries, |
const bool | C_style_indexing = true |
||
) |
Fill this table (which is assumed to already have the correct size) from a source given by dereferencing the given forward iterator (which could, for example, be a pointer to the first element of an array, or an inserting std::istream_iterator). The second argument denotes whether the elements pointed to are arranged in a way that corresponds to the last index running fastest or slowest. The default is to use C-style indexing where the last index runs fastest (as opposed to Fortran-style where the first index runs fastest when traversing multidimensional arrays. For example, if you try to fill an object of type Table<2,T>, then calling this function with the default value for the second argument will result in the equivalent of doing
On the other hand, if the second argument to this function is false, then this would result in code of the following form:
Note the switched order in which we fill the table elements by traversing the given set of iterators.
entries | An iterator to a set of elements from which to initialize this table. It is assumed that iterator can be incremented and dereferenced a sufficient number of times to fill this table. |
C_style_indexing | If true, run over elements of the table with the last index changing fastest as we dereference subsequent elements of the input range. If false, change the first index fastest. |
Fill all table entries with the same value.
AlignedVector< T >::reference TableBase< N, T >::operator() | ( | const TableIndices< N > & | indices | ) |
Return a read-write reference to the indicated element.
AlignedVector< T >::const_reference TableBase< N, T >::operator() | ( | const TableIndices< N > & | indices | ) | const |
Return the value of the indicated element as a read-only reference.
We return the requested value as a constant reference rather than by value since this object may hold data types that may be large, and we don't know here whether copying is expensive or not.
void TableBase< N, T >::replicate_across_communicator | ( | const MPI_Comm | communicator, |
const unsigned int | root_process | ||
) |
This function replicates the state found on the process indicated by root_process
across all processes of the MPI communicator. The current state found on any of the processes other than root_process
is lost in this process. One can imagine this operation to act like a call to Utilities::MPI::broadcast() from the root process to all other processes, though in practice the function may try to move the data into shared memory regions on each of the machines that host MPI processes and let all MPI processes on this machine then access this shared memory region instead of keeping their own copy. See the general documentation of this class for a code example.
The intent of this function is to quickly exchange large arrays from one process to others, rather than having to compute or create it on all processes. This is specifically the case for data loaded from disk – say, large data tables – that are more easily dealt with by reading once and then distributing across all processes in an MPI universe, than letting each process read the data from disk itself. Specifically, the use of shared memory regions allows for replicating the data only once per multicore machine in the MPI universe, rather than replicating data once for each MPI process. This results in large memory savings if the data is large on today's machines that can easily house several dozen MPI processes per shared memory space.
This function does not imply a model of keeping data on different processes in sync, as LinearAlgebra::distributed::Vector and other vector classes do where there exists a notion of certain elements of the vector owned by each process and possibly ghost elements that are mirrored from its owning process to other processes. Rather, the elements of the current object are simply copied to the other processes, and it is useful to think of this operation as creating a set of const
TableBase objects on all processes that should not be changed any more after the replication operation, as this is the only way to ensure that the vectors remain the same on all processes. This is particularly true because of the use of shared memory regions where any modification of a vector element on one MPI process may also result in a modification of elements visible on other processes, assuming they are located within one shared memory node.
communicator
object, which is generally an expensive operation. Likewise, the generation of shared memory spaces is not a cheap operation. As a consequence, this function primarily makes sense when the goal is to share large read-only data tables among processes; examples are data tables that are loaded at start-up time and then used over the course of the run time of the program. In such cases, the start-up cost of running this function can be amortized over time, and the potential memory savings from not having to store the table on each process may be substantial on machines with large core counts on which many MPI processes run on the same machine.T
is "self-contained", i.e., all of its information is stored in its member variables, and if none of the member variables are pointers to other parts of the memory. This is because if a type T
does have pointers to other parts of memory, then moving T
into a shared memory space does not result in the other processes having access to data that the object points to with its member variable pointers: These continue to live only on one process, and are typically in memory areas not accessible to the other processes. As a consequence, the usual use case for this function is to share arrays of simple objects such as double
s or int
s.root_process
throw away their data, which is not a collective operation. Generally, these restrictions on what can and can not be done hint at the correctness of the comments above: You should treat an AlignedVector on which the current function has been called as const
, on which no further operations can be performed until the destructor is called. Swap the contents of this table and the other table v
. One could do this operation with a temporary variable and copying over the data elements, but this function is significantly more efficient since it only swaps the pointers to the data of the two vectors and therefore does not need to allocate temporary storage and move data around.
This function is analogous to the swap
function of all C++ standard containers. Also, there is a global function swap(u,v)
that simply calls u.swap(v)
, again in analogy to standard functions.
Determine an estimate for the memory consumption (in bytes) of this object.
void TableBase< N, T >::serialize | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
|
protected |
Return the position of the indicated element within the array of elements stored one after the other. This function does no index checking.
|
protected |
Return a read-write reference to the indicated element.
This function does no bounds checking and is only to be used internally and in functions already checked.
|
protected |
Return the value of the indicated element as a read-only reference.
This function does no bounds checking and is only to be used internally and in functions already checked.
We return the requested value as a constant reference rather than by value since this object may hold data types that may be large, and we don't know here whether copying is expensive or not.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 131 of file enable_observer_pointer.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 151 of file enable_observer_pointer.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 322 of file enable_observer_pointer.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 339 of file enable_observer_pointer.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 199 of file enable_observer_pointer.cc.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 53 of file enable_observer_pointer.cc.
|
protected |
|
protected |
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 227 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 233 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.
Definition at line 249 of file enable_observer_pointer.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 257 of file enable_observer_pointer.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 280 of file enable_observer_pointer.h.