Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+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
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
DynamicSparsityPattern Class Reference

#include <deal.II/lac/dynamic_sparsity_pattern.h>

Inheritance diagram for DynamicSparsityPattern:
Inheritance graph
[legend]

Classes

struct  Line
 

Public Types

using size_type = types::global_dof_index
 
using iterator = DynamicSparsityPatternIterators::Iterator
 
using const_iterator = DynamicSparsityPatternIterators::Iterator
 

Public Member Functions

 DynamicSparsityPattern ()
 
 DynamicSparsityPattern (const DynamicSparsityPattern &)
 
 DynamicSparsityPattern (const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
 
 DynamicSparsityPattern (const IndexSet &indexset)
 
 DynamicSparsityPattern (const size_type n)
 
DynamicSparsityPatternoperator= (const DynamicSparsityPattern &)
 
void reinit (const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
 
void compress ()
 
bool empty () const
 
size_type max_entries_per_row () const
 
void add (const size_type i, const size_type j)
 
template<typename ForwardIterator >
void add_entries (const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
 
virtual void add_row_entries (const size_type &row, const ArrayView< const size_type > &columns, const bool indices_are_sorted=false) override
 
bool exists (const size_type i, const size_type j) const
 
DynamicSparsityPattern get_view (const IndexSet &rows) const
 
void symmetrize ()
 
template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight >
void compute_mmult_pattern (const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
 
template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight >
void compute_Tmmult_pattern (const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
 
void print (std::ostream &out) const
 
void print_gnuplot (std::ostream &out) const
 
size_type row_length (const size_type row) const
 
void clear_row (const size_type row)
 
size_type column_number (const size_type row, const size_type index) const
 
size_type column_index (const size_type row, const size_type col) const
 
size_type bandwidth () const
 
size_type n_nonzero_elements () const
 
const IndexSetrow_index_set () const
 
IndexSet nonempty_cols () const
 
IndexSet nonempty_rows () const
 
size_type memory_consumption () const
 
virtual void add_entries (const ArrayView< const std::pair< size_type, size_type > > &entries)
 
size_type n_rows () const
 
size_type n_cols () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Iterators
iterator begin () const
 
iterator end () const
 
iterator begin (const size_type r) const
 
iterator end (const size_type r) const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer 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 bool stores_only_added_elements ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

virtual void resize (const size_type rows, const size_type cols)
 

Protected Attributes

size_type rows
 
size_type cols
 

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

bool have_entries
 
IndexSet rowset
 
std::vector< Linelines
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

class DynamicSparsityPatternIterators::Accessor
 

Detailed Description

This class acts as an intermediate form of the SparsityPattern class. From the interface it mostly represents a SparsityPattern object that is kept compressed at all times. However, since the final sparsity pattern is not known while constructing it, keeping the pattern compressed at all times can only be achieved at the expense of either increased memory or run time consumption upon use. The main purpose of this class is to avoid some memory bottlenecks, so we chose to implement it memory conservative. The chosen data format is too unsuited to be used for actual matrices, though. It is therefore necessary to first copy the data of this object over to an object of type SparsityPattern before using it in actual matrices.

Another viewpoint is that this class does not need up front allocation of a certain amount of memory, but grows as necessary. An extensive description of sparsity patterns can be found in the documentation of the Sparsity patterns module.

This class is an example of the "dynamic" type of Sparsity patterns. It is used in most tutorial programs in one way or another.

Interface

Since this class is intended as an intermediate replacement of the SparsityPattern class, it has mostly the same interface, with small changes where necessary. In particular, the add() function, and the functions inquiring properties of the sparsity pattern are the same.

Usage

Usage of this class is explained in step-2 (without constraints) and step-6 (with AffineConstraints) and typically looks as follows:

DynamicSparsityPattern dynamic_pattern (dof_handler.n_dofs());
dynamic_pattern,
constraints);
sp.copy_from (dynamic_pattern);
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)

Definition at line 321 of file dynamic_sparsity_pattern.h.

Member Typedef Documentation

◆ size_type

Declare the type for container size.

Definition at line 327 of file dynamic_sparsity_pattern.h.

◆ iterator

Typedef an for iterator class that allows to walk over all nonzero elements of a sparsity pattern.

Since the iterator does not allow to modify the sparsity pattern, this type is the same as that for const_iterator.

Definition at line 336 of file dynamic_sparsity_pattern.h.

◆ const_iterator

Typedef for an iterator class that allows to walk over all nonzero elements of a sparsity pattern.

Definition at line 342 of file dynamic_sparsity_pattern.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ DynamicSparsityPattern() [1/5]

DynamicSparsityPattern::DynamicSparsityPattern ( )

Initialize as an empty object. This is useful if you want such objects as member variables in other classes. You can make the structure usable by calling the reinit() function.

Definition at line 221 of file dynamic_sparsity_pattern.cc.

◆ DynamicSparsityPattern() [2/5]

DynamicSparsityPattern::DynamicSparsityPattern ( const DynamicSparsityPattern s)

Copy constructor. This constructor is only allowed to be called if the sparsity structure to be copied is empty. This is so in order to prevent involuntary copies of objects for temporaries, which can use large amounts of computing time. However, copy constructors are needed if you want to place a DynamicSparsityPattern in a container, e.g. to write such statements like v.push_back (DynamicSparsityPattern());, with v a vector of DynamicSparsityPattern objects.

Definition at line 229 of file dynamic_sparsity_pattern.cc.

◆ DynamicSparsityPattern() [3/5]

DynamicSparsityPattern::DynamicSparsityPattern ( const size_type  m,
const size_type  n,
const IndexSet rowset = IndexSet() 
)

Initialize a rectangular sparsity pattern with m rows and n columns. The rowset restricts the storage to elements in rows of this set. Adding elements outside of this set has no effect. The default argument keeps all entries.

Definition at line 244 of file dynamic_sparsity_pattern.cc.

◆ DynamicSparsityPattern() [4/5]

DynamicSparsityPattern::DynamicSparsityPattern ( const IndexSet indexset)

Create a square SparsityPattern using the given index set. The total size is given by the size of indexset and only rows corresponding to indices in indexset are stored on the current processor.

Definition at line 255 of file dynamic_sparsity_pattern.cc.

◆ DynamicSparsityPattern() [5/5]

DynamicSparsityPattern::DynamicSparsityPattern ( const size_type  n)

Initialize a square pattern of dimension n.

Definition at line 260 of file dynamic_sparsity_pattern.cc.

Member Function Documentation

◆ operator=()

DynamicSparsityPattern & DynamicSparsityPattern::operator= ( const DynamicSparsityPattern s)

Copy operator. For this the same holds as for the copy constructor: it is declared, defined and fine to be called, but the latter only for empty objects.

Definition at line 271 of file dynamic_sparsity_pattern.cc.

◆ reinit()

void DynamicSparsityPattern::reinit ( const size_type  m,
const size_type  n,
const IndexSet rowset = IndexSet() 
)

Reallocate memory and set up data structures for a new sparsity pattern with m rows and n columns. The rowset restricts the storage to elements in rows of this set. Adding elements outside of this set has no effect. The default argument keeps all entries.

Definition at line 290 of file dynamic_sparsity_pattern.cc.

◆ compress()

void DynamicSparsityPattern::compress ( )

Since this object is kept compressed at all times anyway, this function does nothing, but is declared to make the interface of this class as much alike as that of the SparsityPattern class.

Definition at line 315 of file dynamic_sparsity_pattern.cc.

◆ empty()

bool DynamicSparsityPattern::empty ( ) const

Return whether the object is empty. It is empty if no memory is allocated, which is the same as that both dimensions are zero.

Definition at line 321 of file dynamic_sparsity_pattern.cc.

◆ max_entries_per_row()

DynamicSparsityPattern::size_type DynamicSparsityPattern::max_entries_per_row ( ) const

Return the maximum number of entries per row. Note that this number may change as entries are added.

Definition at line 329 of file dynamic_sparsity_pattern.cc.

◆ add()

void DynamicSparsityPattern::add ( const size_type  i,
const size_type  j 
)
inline

Add a nonzero entry. If the entry already exists, this call does nothing.

Definition at line 1002 of file dynamic_sparsity_pattern.h.

◆ add_entries() [1/2]

template<typename ForwardIterator >
void DynamicSparsityPattern::add_entries ( const size_type  row,
ForwardIterator  begin,
ForwardIterator  end,
const bool  indices_are_unique_and_sorted = false 
)
inline

Add several nonzero entries to the specified row. Already existing entries are ignored.

Definition at line 1021 of file dynamic_sparsity_pattern.h.

◆ add_row_entries()

void DynamicSparsityPattern::add_row_entries ( const size_type row,
const ArrayView< const size_type > &  columns,
const bool  indices_are_sorted = false 
)
overridevirtual

Optimized function for adding new entries to a given row.

Implements SparsityPatternBase.

Definition at line 346 of file dynamic_sparsity_pattern.cc.

◆ exists()

bool DynamicSparsityPattern::exists ( const size_type  i,
const size_type  j 
) const

Check if a value at a certain position may be non-zero.

Definition at line 357 of file dynamic_sparsity_pattern.cc.

◆ get_view()

DynamicSparsityPattern DynamicSparsityPattern::get_view ( const IndexSet rows) const

Return a view of this sparsity pattern. That is, for all rows in rows extract non-empty columns. The resulting sparsity pattern will have number of rows equal rows.n_elements().

Definition at line 432 of file dynamic_sparsity_pattern.cc.

◆ symmetrize()

void DynamicSparsityPattern::symmetrize ( )

Make the sparsity pattern symmetric by adding the sparsity pattern of the transpose object.

This function throws an exception if the sparsity pattern does not represent a square matrix.

Definition at line 386 of file dynamic_sparsity_pattern.cc.

◆ compute_mmult_pattern()

template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight >
template void DynamicSparsityPattern::compute_mmult_pattern ( const SparsityPatternTypeLeft &  left,
const SparsityPatternTypeRight &  right 
)

Construct and store in this object the sparsity pattern corresponding to the product of left and right sparsity pattern.

Definition at line 496 of file dynamic_sparsity_pattern.cc.

◆ compute_Tmmult_pattern()

template<typename SparsityPatternTypeLeft , typename SparsityPatternTypeRight >
template void DynamicSparsityPattern::compute_Tmmult_pattern ( const SparsityPatternTypeLeft &  left,
const SparsityPatternTypeRight &  right 
)

Construct and store in this object the sparsity pattern corresponding to the product of transposed left and non-transpose right sparsity pattern.

Definition at line 455 of file dynamic_sparsity_pattern.cc.

◆ print()

void DynamicSparsityPattern::print ( std::ostream &  out) const

Print the sparsity pattern. The output consists of one line per row of the format [i,j1,j2,j3,...]. i is the row number and jn are the allocated columns in this row.

Definition at line 525 of file dynamic_sparsity_pattern.cc.

◆ print_gnuplot()

void DynamicSparsityPattern::print_gnuplot ( std::ostream &  out) const

Print the sparsity pattern in a format that gnuplot understands and which can be used to plot the sparsity pattern in a graphical way. The format consists of pairs i j of nonzero elements, each representing one entry, one per line of the output file. Indices are counted from zero on, as usual. Since sparsity patterns are printed in the same way as matrices are displayed, we print the negative of the column index, which means that the (0,0) element is in the top left rather than in the bottom left corner.

Print the sparsity pattern in gnuplot by setting the data style to dots or points and use the plot command.

Definition at line 543 of file dynamic_sparsity_pattern.cc.

◆ row_length()

types::global_dof_index DynamicSparsityPattern::row_length ( const size_type  row) const
inline

Number of entries in a specific row. This function can only be called if the given row is a member of the index set of rows that we want to store.

Definition at line 1042 of file dynamic_sparsity_pattern.h.

◆ clear_row()

void DynamicSparsityPattern::clear_row ( const size_type  row)

Clear all entries stored in a specific row.

Definition at line 413 of file dynamic_sparsity_pattern.cc.

◆ column_number()

types::global_dof_index DynamicSparsityPattern::column_number ( const size_type  row,
const size_type  index 
) const
inline

Access to column number field. Return the column number of the indexth entry in row.

Definition at line 1060 of file dynamic_sparsity_pattern.h.

◆ column_index()

types::global_dof_index DynamicSparsityPattern::column_index ( const size_type  row,
const size_type  col 
) const

Return index of column col in row row. If the column does not exist in this sparsity pattern, the returned value will be 'numbers::invalid_size_type'.

Definition at line 656 of file dynamic_sparsity_pattern.cc.

◆ begin() [1/2]

DynamicSparsityPattern::iterator DynamicSparsityPattern::begin ( ) const
inline

Iterator starting at the first entry of the matrix. The resulting iterator can be used to walk over all nonzero entries of the sparsity pattern.

Note the discussion in the general documentation of this class about the order in which elements are accessed.

Note
If the sparsity pattern has been initialized with an IndexSet that denotes which rows to store, then iterators will simply skip over rows that are not stored. In other words, they will look like empty rows, but no exception will be generated when iterating over such rows.

Definition at line 1075 of file dynamic_sparsity_pattern.h.

◆ end() [1/2]

DynamicSparsityPattern::iterator DynamicSparsityPattern::end ( ) const
inline

Final iterator.

Definition at line 1085 of file dynamic_sparsity_pattern.h.

◆ begin() [2/2]

DynamicSparsityPattern::iterator DynamicSparsityPattern::begin ( const size_type  r) const
inline

Iterator starting at the first entry of row r.

Note that if the given row is empty, i.e. does not contain any nonzero entries, then the iterator returned by this function equals end(r). Note also that the iterator may not be dereferenceable in that case.

Note also the discussion in the general documentation of this class about the order in which elements are accessed.

Note
If the sparsity pattern has been initialized with an IndexSet that denotes which rows to store, then iterators will simply skip over rows that are not stored. In other words, they will look like empty rows, but no exception will be generated when iterating over such rows.

Definition at line 1093 of file dynamic_sparsity_pattern.h.

◆ end() [2/2]

DynamicSparsityPattern::iterator DynamicSparsityPattern::end ( const size_type  r) const
inline

Final iterator of row r. It points to the first element past the end of line r, or past the end of the entire sparsity pattern.

Note that the end iterator is not necessarily dereferenceable. This is in particular the case if it is the end iterator for the last row of a matrix.

Definition at line 1151 of file dynamic_sparsity_pattern.h.

◆ bandwidth()

DynamicSparsityPattern::size_type DynamicSparsityPattern::bandwidth ( ) const

Compute the bandwidth of the matrix represented by this structure. The bandwidth is the maximum of \(|i-j|\) for which the index pair \((i,j)\) represents a nonzero entry of the matrix.

Definition at line 566 of file dynamic_sparsity_pattern.cc.

◆ n_nonzero_elements()

DynamicSparsityPattern::size_type DynamicSparsityPattern::n_nonzero_elements ( ) const

Return the number of nonzero elements allocated through this sparsity pattern.

Definition at line 586 of file dynamic_sparsity_pattern.cc.

◆ row_index_set()

const IndexSet & DynamicSparsityPattern::row_index_set ( ) const
inline

Return the IndexSet that sets which rows are active on the current processor. It corresponds to the IndexSet given to this class in the constructor or in the reinit function.

Definition at line 1165 of file dynamic_sparsity_pattern.h.

◆ nonempty_cols()

IndexSet DynamicSparsityPattern::nonempty_cols ( ) const

Return the IndexSet that contains entries for all columns in which at least one element exists in this sparsity pattern.

Note
In a parallel context, this only considers the locally stored rows.

Definition at line 603 of file dynamic_sparsity_pattern.cc.

◆ nonempty_rows()

IndexSet DynamicSparsityPattern::nonempty_rows ( ) const

Return the IndexSet that contains entries for all rows in which at least one element exists in this sparsity pattern.

Note
In a parallel context, this only considers the locally stored rows.

Definition at line 617 of file dynamic_sparsity_pattern.cc.

◆ stores_only_added_elements()

bool DynamicSparsityPattern::stores_only_added_elements ( )
inlinestatic

return whether this object stores only those entries that have been added explicitly, or if the sparsity pattern contains elements that have been added through other means (implicitly) while building it. For the current class, the result is always true.

This function mainly serves the purpose of describing the current class in cases where several kinds of sparsity patterns can be passed as template arguments.

Definition at line 1173 of file dynamic_sparsity_pattern.h.

◆ memory_consumption()

DynamicSparsityPattern::size_type DynamicSparsityPattern::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 641 of file dynamic_sparsity_pattern.cc.

◆ add_entries() [2/2]

void SparsityPatternBase::add_entries ( const ArrayView< const std::pair< size_type, size_type > > &  entries)
virtual

General function for adding new entries from an unsorted list of pairs.

This function is useful when multiple entries need to be added which do not correspond to a particular row, e.g., when assembling a flux sparsity pattern.

Reimplemented from SparsityPatternBase.

Definition at line 109 of file sparsity_pattern_base.cc.

◆ n_rows()

size_type SparsityPatternBase::n_rows ( ) const
inherited

Return number of rows of this matrix, which equals the dimension of the image space.

◆ n_cols()

size_type SparsityPatternBase::n_cols ( ) const
inherited

Return number of columns of this matrix, which equals the dimension of the range space.

◆ resize()

virtual void SparsityPatternBase::resize ( const size_type  rows,
const size_type  cols 
)
protectedvirtualinherited

Internal function for updating the stored size of the sparsity pattern.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
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 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
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 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
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.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ DynamicSparsityPatternIterators::Accessor

Definition at line 723 of file dynamic_sparsity_pattern.h.

Member Data Documentation

◆ have_entries

bool DynamicSparsityPattern::have_entries
private

A flag that stores whether any entries have been added so far.

Definition at line 670 of file dynamic_sparsity_pattern.h.

◆ rowset

IndexSet DynamicSparsityPattern::rowset
private

A set that contains the valid rows.

Definition at line 676 of file dynamic_sparsity_pattern.h.

◆ lines

std::vector<Line> DynamicSparsityPattern::lines
private

Actual data: store for each row the set of nonzero entries.

Definition at line 720 of file dynamic_sparsity_pattern.h.

◆ rows

size_type SparsityPatternBase::rows
protectedinherited

Number of rows that this sparsity pattern shall represent.

Definition at line 121 of file sparsity_pattern_base.h.

◆ cols

size_type SparsityPatternBase::cols
protectedinherited

Number of columns that this sparsity pattern shall represent.

Definition at line 126 of file sparsity_pattern_base.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
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 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
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 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
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 271 of file subscriptor.h.


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