Reference documentation for deal.II version 9.5.0
\(\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 Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
ChunkSparsityPattern Class Reference

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

Inheritance diagram for ChunkSparsityPattern:
[legend]

Public Types

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

Public Member Functions

 ChunkSparsityPattern ()
 
 ChunkSparsityPattern (const ChunkSparsityPattern &)
 
 ChunkSparsityPattern (const size_type m, const size_type n, const size_type max_chunks_per_row, const size_type chunk_size)
 
 ChunkSparsityPattern (const size_type m, const size_type n, const std::vector< size_type > &row_lengths, const size_type chunk_size)
 
 ChunkSparsityPattern (const size_type n, const size_type max_per_row, const size_type chunk_size)
 
 ChunkSparsityPattern (const size_type m, const std::vector< size_type > &row_lengths, const size_type chunk_size)
 
 ~ChunkSparsityPattern () override=default
 
ChunkSparsityPatternoperator= (const ChunkSparsityPattern &)
 
void reinit (const size_type m, const size_type n, const size_type max_per_row, const size_type chunk_size)
 
void reinit (const size_type m, const size_type n, const std::vector< size_type > &row_lengths, const size_type chunk_size)
 
void reinit (const size_type m, const size_type n, const ArrayView< const size_type > &row_lengths, const size_type chunk_size)
 
void compress ()
 
template<typename ForwardIterator >
void copy_from (const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end, const size_type chunk_size)
 
template<typename SparsityPatternType >
void copy_from (const SparsityPatternType &dsp, const size_type chunk_size)
 
template<typename number >
void copy_from (const FullMatrix< number > &matrix, const size_type chunk_size)
 
template<typename Sparsity >
void create_from (const size_type m, const size_type n, const Sparsity &sparsity_pattern_for_chunks, const size_type chunk_size, const bool optimize_diagonal=true)
 
bool empty () const
 
size_type get_chunk_size () const
 
size_type max_entries_per_row () const
 
void add (const size_type i, const size_type j)
 
void symmetrize ()
 
size_type n_rows () const
 
size_type n_cols () const
 
bool exists (const size_type i, const size_type j) const
 
size_type row_length (const size_type row) const
 
size_type bandwidth () const
 
size_type n_nonzero_elements () const
 
bool is_compressed () const
 
bool stores_only_added_elements () const
 
iterator begin () const
 
iterator end () const
 
iterator begin (const size_type r) const
 
iterator end (const size_type r) const
 
void block_write (std::ostream &out) const
 
void block_read (std::istream &in)
 
void print (std::ostream &out) const
 
void print_gnuplot (std::ostream &out) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
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 ::ExceptionBaseExcInvalidNumber (size_type arg1)
 
static ::ExceptionBaseExcInvalidIndex (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcNotEnoughSpace (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcNotCompressed ()
 
static ::ExceptionBaseExcMatrixIsCompressed ()
 
static ::ExceptionBaseExcEmptyObject ()
 
static ::ExceptionBaseExcIteratorRange (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcMETISNotInstalled ()
 
static ::ExceptionBaseExcInvalidNumberOfPartitions (size_type arg1)
 
static ::ExceptionBaseExcInvalidArraySize (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Static Public Attributes

static const size_type invalid_entry = SparsityPattern::invalid_entry
 

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

size_type rows
 
size_type cols
 
size_type chunk_size
 
SparsityPattern sparsity_pattern
 
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

template<typename >
class ChunkSparseMatrix
 
class ChunkSparsityPatternIterators::Accessor
 

Detailed Description

Structure representing the sparsity pattern of a sparse matrix. This class is an example of the "static" type of Sparsity. It uses the compressed row storage (CSR) format to store data.

The use of this class is demonstrated in step-51.

Definition at line 245 of file chunk_sparsity_pattern.h.

Member Typedef Documentation

◆ size_type

Declare the type for container size.

Definition at line 251 of file chunk_sparsity_pattern.h.

◆ const_iterator

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

Definition at line 256 of file chunk_sparsity_pattern.h.

◆ iterator

Typedef an 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 265 of file chunk_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 230 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 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ ChunkSparsityPattern() [1/6]

ChunkSparsityPattern::ChunkSparsityPattern ( )

Initialize the matrix empty, that is with no memory allocated. 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 25 of file chunk_sparsity_pattern.cc.

◆ ChunkSparsityPattern() [2/6]

ChunkSparsityPattern::ChunkSparsityPattern ( const ChunkSparsityPattern s)

Copy constructor. This constructor is only allowed to be called if the matrix 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 one wants to place a ChunkSparsityPattern in a container, e.g., to write such statements like v.push_back (ChunkSparsityPattern());, with v a vector of ChunkSparsityPattern objects.

Usually, it is sufficient to use the explicit keyword to disallow unwanted temporaries, but this does not work for std::vector. Since copying a structure like this is not useful anyway because multiple matrices can use the same sparsity structure, copies are only allowed for empty objects, as described above.

Definition at line 32 of file chunk_sparsity_pattern.cc.

◆ ChunkSparsityPattern() [3/6]

ChunkSparsityPattern::ChunkSparsityPattern ( const size_type  m,
const size_type  n,
const size_type  max_chunks_per_row,
const size_type  chunk_size 
)

Initialize a rectangular matrix.

  • m number of rows
  • n number of columns
  • max_per_row maximum number of nonzero entries per row

Definition at line 48 of file chunk_sparsity_pattern.cc.

◆ ChunkSparsityPattern() [4/6]

ChunkSparsityPattern::ChunkSparsityPattern ( const size_type  m,
const size_type  n,
const std::vector< size_type > &  row_lengths,
const size_type  chunk_size 
)

Initialize a rectangular matrix.

  • m number of rows
  • n number of columns
  • row_lengths possible number of nonzero entries for each row. This vector must have one entry for each row.

Definition at line 60 of file chunk_sparsity_pattern.cc.

◆ ChunkSparsityPattern() [5/6]

ChunkSparsityPattern::ChunkSparsityPattern ( const size_type  n,
const size_type  max_per_row,
const size_type  chunk_size 
)

Initialize a quadratic matrix of dimension n with at most max_per_row nonzero entries per row.

This constructor automatically enables optimized storage of diagonal elements. To avoid this, use the constructor taking row and column numbers separately.

Definition at line 73 of file chunk_sparsity_pattern.cc.

◆ ChunkSparsityPattern() [6/6]

ChunkSparsityPattern::ChunkSparsityPattern ( const size_type  m,
const std::vector< size_type > &  row_lengths,
const size_type  chunk_size 
)

Initialize a quadratic matrix.

  • m number of rows and columns
  • row_lengths possible number of nonzero entries for each row. This vector must have one entry for each row.

Definition at line 82 of file chunk_sparsity_pattern.cc.

◆ ~ChunkSparsityPattern()

ChunkSparsityPattern::~ChunkSparsityPattern ( )
overridedefault

Destructor.

Member Function Documentation

◆ operator=()

ChunkSparsityPattern & ChunkSparsityPattern::operator= ( const ChunkSparsityPattern 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 95 of file chunk_sparsity_pattern.cc.

◆ reinit() [1/3]

void ChunkSparsityPattern::reinit ( const size_type  m,
const size_type  n,
const size_type  max_per_row,
const size_type  chunk_size 
)

Reallocate memory and set up data structures for a new matrix with m rows and n columns, with at most max_per_row nonzero entries per row.

This function simply maps its operations to the other reinit function.

Definition at line 116 of file chunk_sparsity_pattern.cc.

◆ reinit() [2/3]

void ChunkSparsityPattern::reinit ( const size_type  m,
const size_type  n,
const std::vector< size_type > &  row_lengths,
const size_type  chunk_size 
)

Reallocate memory for a matrix of size m x n. The number of entries for each row is taken from the array row_lengths which has to give this number of each row i=1...m.

If m*n==0 all memory is freed, resulting in a total reinitialization of the object. If it is nonzero, new memory is only allocated if the new size extends the old one. This is done to save time and to avoid fragmentation of the heap.

If the number of rows equals the number of columns then diagonal elements are stored first in each row to allow optimized access in relaxation methods of SparseMatrix.

Definition at line 260 of file chunk_sparsity_pattern.cc.

◆ reinit() [3/3]

void ChunkSparsityPattern::reinit ( const size_type  m,
const size_type  n,
const ArrayView< const size_type > &  row_lengths,
const size_type  chunk_size 
)

Same as above, but with an ArrayView argument instead.

Definition at line 131 of file chunk_sparsity_pattern.cc.

◆ compress()

void ChunkSparsityPattern::compress ( )

This function compresses the sparsity structure that this object represents. It does so by eliminating unused entries and sorting the remaining ones to allow faster access by usage of binary search algorithms. A special sorting scheme is used for the diagonal entry of quadratic matrices, which is always the first entry of each row.

The memory which is no more needed is released.

SparseMatrix objects require the ChunkSparsityPattern objects they are initialized with to be compressed, to reduce memory requirements.

Definition at line 174 of file chunk_sparsity_pattern.cc.

◆ copy_from() [1/3]

template<typename ForwardIterator >
void ChunkSparsityPattern::copy_from ( const size_type  n_rows,
const size_type  n_cols,
const ForwardIterator  begin,
const ForwardIterator  end,
const size_type  chunk_size 
)

This function can be used as a replacement for reinit(), subsequent calls to add() and a final call to close() if you know exactly in advance the entries that will form the matrix sparsity pattern.

The first two parameters determine the size of the matrix. For the two last ones, note that a sparse matrix can be described by a sequence of rows, each of which is represented by a sequence of pairs of column indices and values. In the present context, the begin() and end() parameters designate iterators (of forward iterator type) into a container, one representing one row. The distance between begin() and end() should therefore be equal to n_rows(). These iterators may be iterators of std::vector, std::list, pointers into a C-style array, or any other iterator satisfying the requirements of a forward iterator. The objects pointed to by these iterators (i.e. what we get after applying operator* or operator-> to one of these iterators) must be a container itself that provides functions begin and end designating a range of iterators that describe the contents of one line. Dereferencing these inner iterators must either yield a pair of an integer as column index and a value of arbitrary type (such a type would be used if we wanted to describe a sparse matrix with one such object), or simply an integer (if we only wanted to describe a sparsity pattern). The function is able to determine itself whether an integer or a pair is what we get after dereferencing the inner iterators, through some template magic.

While the order of the outer iterators denotes the different rows of the matrix, the order of the inner iterator denoting the columns does not matter, as they are sorted internal to this function anyway.

Since that all sounds very complicated, consider the following example code, which may be used to fill a sparsity pattern:

std::vector<std::vector<size_type> > column_indices (n_rows);
for (size_type row=0; row<n_rows; ++row)
// generate necessary columns in this row
fill_row (column_indices[row]);
sparsity.copy_from (n_rows, n_cols,
column_indices.begin(),
column_indices.end());
size_type n_cols() const
size_type n_rows() const

Note that this example works since the iterators dereferenced yield containers with functions begin and end (namely std::vectors), and the inner iterators dereferenced yield integers as column indices. Note that we could have replaced each of the two std::vector occurrences by std::list, and the inner one by std::set as well.

Another example would be as follows, where we initialize a whole matrix, not only a sparsity pattern:

std::vector<std::map<size_type,double> > entries (n_rows);
for (size_type row=0; row<n_rows; ++row)
// generate necessary pairs of columns
// and corresponding values in this row
fill_row (entries[row]);
sparsity.copy_from (n_rows, n_cols,
column_indices.begin(),
column_indices.end());
matrix.reinit (sparsity);
matrix.copy_from (column_indices.begin(),
column_indices.end());

This example works because dereferencing iterators of the inner type yields a pair of integers and a value, the first of which we take as column index. As previously, the outer std::vector could be replaced by std::list, and the inner std::map<size_type,double> could be replaced by std::vector<std::pair<size_type,double> >, or a list or set of such pairs, as they all return iterators that point to such pairs.

◆ copy_from() [2/3]

template<typename SparsityPatternType >
void ChunkSparsityPattern::copy_from ( const SparsityPatternType &  dsp,
const size_type  chunk_size 
)

Copy data from an object of type DynamicSparsityPattern. Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 183 of file chunk_sparsity_pattern.cc.

◆ copy_from() [3/3]

template<typename number >
void ChunkSparsityPattern::copy_from ( const FullMatrix< number > &  matrix,
const size_type  chunk_size 
)

Take a full matrix and use its nonzero entries to generate a sparse matrix entry pattern for this object.

Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 225 of file chunk_sparsity_pattern.cc.

◆ create_from()

template<typename Sparsity >
void ChunkSparsityPattern::create_from ( const size_type  m,
const size_type  n,
const Sparsity &  sparsity_pattern_for_chunks,
const size_type  chunk_size,
const bool  optimize_diagonal = true 
)

Set the sparsity pattern of the chunk sparsity pattern to be given by chunk_size*chunksize blocks of the sparsity pattern for chunks specified. Note that the final number of rows m of the sparsity pattern will be approximately sparsity_pattern_for_chunks.n_rows() * chunk_size (modulo padding elements in the last chunk) and similarly for the number of columns n.

This is a special initialization option in case you can tell the position of the chunk already from the beginning without generating the sparsity pattern using make_sparsity_pattern calls. This bypasses the search for chunks but of course needs to be handled with care in order to give a correct sparsity pattern.

Previous content of this object is lost, and the sparsity pattern is in compressed mode afterwards.

Definition at line 295 of file chunk_sparsity_pattern.cc.

◆ empty()

bool ChunkSparsityPattern::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 320 of file chunk_sparsity_pattern.cc.

◆ get_chunk_size()

size_type ChunkSparsityPattern::get_chunk_size ( ) const

Return the chunk size given as argument when constructing this object.

◆ max_entries_per_row()

ChunkSparsityPattern::size_type ChunkSparsityPattern::max_entries_per_row ( ) const

Return the maximum number of entries per row. Before compression, this equals the number given to the constructor, while after compression, it equals the maximum number of entries actually allocated by the user.

Definition at line 328 of file chunk_sparsity_pattern.cc.

◆ add()

void ChunkSparsityPattern::add ( const size_type  i,
const size_type  j 
)

Add a nonzero entry to the matrix. This function may only be called for non-compressed sparsity patterns.

If the entry already exists, nothing bad happens.

Definition at line 336 of file chunk_sparsity_pattern.cc.

◆ symmetrize()

void ChunkSparsityPattern::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 quadratic matrix.

Definition at line 357 of file chunk_sparsity_pattern.cc.

◆ n_rows()

size_type ChunkSparsityPattern::n_rows ( ) const
inline

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

◆ n_cols()

size_type ChunkSparsityPattern::n_cols ( ) const
inline

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

◆ exists()

bool ChunkSparsityPattern::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 346 of file chunk_sparsity_pattern.cc.

◆ row_length()

ChunkSparsityPattern::size_type ChunkSparsityPattern::row_length ( const size_type  row) const

Number of entries in a specific row.

Definition at line 370 of file chunk_sparsity_pattern.cc.

◆ bandwidth()

ChunkSparsityPattern::size_type ChunkSparsityPattern::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. Consequently, the maximum bandwidth a \(n\times m\) matrix can have is \(\max\{n-1,m-1\}\).

Definition at line 520 of file chunk_sparsity_pattern.cc.

◆ n_nonzero_elements()

ChunkSparsityPattern::size_type ChunkSparsityPattern::n_nonzero_elements ( ) const

Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.

This function may only be called if the matrix struct is compressed. It does not make too much sense otherwise anyway.

Definition at line 398 of file chunk_sparsity_pattern.cc.

◆ is_compressed()

bool ChunkSparsityPattern::is_compressed ( ) const

Return whether the structure is compressed or not.

◆ stores_only_added_elements()

bool ChunkSparsityPattern::stores_only_added_elements ( ) const

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 true if and only if it is square because it then unconditionally stores the diagonal entries whether they have been added explicitly or not.

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 535 of file chunk_sparsity_pattern.cc.

◆ begin() [1/2]

iterator ChunkSparsityPattern::begin ( ) const

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.

◆ end() [1/2]

iterator ChunkSparsityPattern::end ( ) const

Final iterator.

◆ begin() [2/2]

iterator ChunkSparsityPattern::begin ( const size_type  r) const

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.

◆ end() [2/2]

iterator ChunkSparsityPattern::end ( const size_type  r) const

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.

◆ block_write()

void ChunkSparsityPattern::block_write ( std::ostream &  out) const

Write the data of this object en bloc to a file. This is done in a binary mode, so the output is neither readable by humans nor (probably) by other computers using a different operating system of number format.

The purpose of this function is that you can swap out matrices and sparsity pattern if you are short of memory, want to communicate between different programs, or allow objects to be persistent across different runs of the program.

Definition at line 546 of file chunk_sparsity_pattern.cc.

◆ block_read()

void ChunkSparsityPattern::block_read ( std::istream &  in)

Read data that has previously been written by block_write() from a file. This is done using the inverse operations to the above function, so it is reasonably fast because the bitstream is not interpreted except for a few numbers up front.

The object is resized on this operation, and all previous contents are lost.

A primitive form of error checking is performed which will recognize the bluntest attempts to interpret some data as a vector stored bitwise to a file, but not more.

Definition at line 562 of file chunk_sparsity_pattern.cc.

◆ print()

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

Print the sparsity of the matrix. 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 455 of file chunk_sparsity_pattern.cc.

◆ print_gnuplot()

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

Print the sparsity of the matrix 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 of this matrix, 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 486 of file chunk_sparsity_pattern.cc.

◆ memory_consumption()

std::size_t ChunkSparsityPattern::memory_consumption ( ) const

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

Definition at line 588 of file chunk_sparsity_pattern.cc.

◆ 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 136 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 156 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 204 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 53 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ ChunkSparseMatrix

template<typename >
friend class ChunkSparseMatrix
friend

Definition at line 855 of file chunk_sparsity_pattern.h.

◆ ChunkSparsityPatternIterators::Accessor

Definition at line 858 of file chunk_sparsity_pattern.h.

Member Data Documentation

◆ invalid_entry

const size_type ChunkSparsityPattern::invalid_entry = SparsityPattern::invalid_entry
static

Define a value which is used to indicate that a certain value in the colnums array is unused, i.e. does not represent a certain column number index.

Indices with this invalid value are used to insert new entries to the sparsity pattern using the add() member function, and are removed when calling compress().

You should not assume that the variable declared here has a certain value. The initialization is given here only to enable the compiler to perform some optimizations, but the actual value of the variable may change over time.

Definition at line 281 of file chunk_sparsity_pattern.h.

◆ rows

size_type ChunkSparsityPattern::rows
private

Number of rows that this sparsity structure shall represent.

Definition at line 835 of file chunk_sparsity_pattern.h.

◆ cols

size_type ChunkSparsityPattern::cols
private

Number of columns that this sparsity structure shall represent.

Definition at line 840 of file chunk_sparsity_pattern.h.

◆ chunk_size

size_type ChunkSparsityPattern::chunk_size
private

The size of chunks.

Definition at line 845 of file chunk_sparsity_pattern.h.

◆ sparsity_pattern

SparsityPattern ChunkSparsityPattern::sparsity_pattern
private

The reduced sparsity pattern. We store only which chunks exist, with each chunk a block in the matrix of size chunk_size by chunk_size.

Definition at line 851 of file chunk_sparsity_pattern.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 219 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 225 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 241 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 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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