MueLu
Version of the Day
|
Container class for aggregation information. More...
#include <MueLu_Aggregates_decl.hpp>
Public Types | |
using | local_ordinal_type = LocalOrdinal |
using | global_ordinal_type = GlobalOrdinal |
using | execution_space = typename Node::execution_space |
using | node_type = Node |
using | device_type = typename Node::device_type |
using | range_type = Kokkos::RangePolicy< local_ordinal_type, execution_space > |
using | LO_view = Kokkos::View< local_ordinal_type *, device_type > |
using | aggregates_sizes_type = Kokkos::View< LocalOrdinal *, device_type > |
using | local_graph_type = typename LWGraph_kokkos::local_graph_type |
using | colors_view_type = Kokkos::View< typename local_graph_type::entries_type::data_type, typename local_graph_type::device_type::memory_space > |
Public Member Functions | |
Aggregates (const LWGraph &graph) | |
Standard constructor for Aggregates structure. More... | |
Aggregates (LWGraph_kokkos graph) | |
Standard constructor for Aggregates structure. More... | |
Aggregates (const RCP< const Map > &map) | |
Constructor for Aggregates structure. More... | |
virtual | ~Aggregates () |
Destructor. More... | |
void | SetNumAggregates (LO nAggregates) |
Set number of local aggregates on current processor. More... | |
void | SetNumGlobalAggregates (GO nGlobalAggregates) |
Set number of global aggregates on current processor. More... | |
KOKKOS_INLINE_FUNCTION LO | GetNumAggregates () const |
KOKKOS_INLINE_FUNCTION void | AggregatesCrossProcessors (const bool &flag) |
Record whether aggregates include DOFs from other processes. More... | |
KOKKOS_INLINE_FUNCTION bool | AggregatesCrossProcessors () const |
Return false if and only if no aggregates include DOFs from other processes. More... | |
RCP< LOMultiVector > & | GetVertex2AggIdNonConst () |
Returns a nonconstant vector that maps local node IDs to local aggregates IDs. More... | |
RCP< LOVector > & | GetProcWinnerNonConst () |
Returns nonconstant vector that maps local node IDs to owning processor IDs. More... | |
const RCP< LOMultiVector > & | GetVertex2AggId () const |
Returns constant vector that maps local node IDs to local aggregates IDs. More... | |
const RCP< LOVector > & | GetProcWinner () const |
Returns constant vector that maps local node IDs to owning processor IDs. More... | |
bool | IsRoot (LO i) const |
Returns true if node with given local node id is marked to be a root node. More... | |
void | SetIsRoot (LO i, bool value=true) |
Set root node information. More... | |
const RCP< const Map > | GetMap () const |
returns (overlapping) map of aggregate/node distribution More... | |
aggregates_sizes_type::const_type | ComputeAggregateSizes (bool forceRecompute=false) const |
Compute sizes of aggregates. More... | |
Teuchos::ArrayRCP< LocalOrdinal > | ComputeAggregateSizesArrayRCP (bool forceRecompute=false) const |
Compute sizes of aggregates. More... | |
local_graph_type | GetGraph () const |
void | ComputeNodesInAggregate (LO_view &aggPtr, LO_view &aggNodes, LO_view &unaggregated) const |
Generates a compressed list of nodes in each aggregate, where the entries in aggNodes[aggPtr[i]] up to aggNodes[aggPtr[i+1]-1] contain the nodes in aggregate i. unaggregated contains the list of nodes which are, for whatever reason, not aggregated (e.g. Dirichlet) More... | |
GO | GetNumGlobalAggregatesComputeIfNeeded () |
Get global number of aggregates. More... | |
Public Member Functions inherited from MueLu::BaseClass | |
virtual | ~BaseClass () |
Destructor. More... | |
Public Member Functions inherited from MueLu::VerboseObject | |
VerbLevel | GetVerbLevel () const |
Get the verbosity level. More... | |
void | SetVerbLevel (const VerbLevel verbLevel) |
Set the verbosity level of this object. More... | |
int | GetProcRankVerbose () const |
Get proc rank used for printing. Do not use this information for any other purpose. More... | |
int | SetProcRankVerbose (int procRank) const |
Set proc rank used for printing. More... | |
bool | IsPrint (MsgType type, int thisProcRankOnly=-1) const |
Find out whether we need to print out information for a specific message type. More... | |
Teuchos::FancyOStream & | GetOStream (MsgType type, int thisProcRankOnly=0) const |
Get an output stream for outputting the input message type. More... | |
Teuchos::FancyOStream & | GetBlackHole () const |
VerboseObject () | |
virtual | ~VerboseObject () |
Destructor. More... | |
Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
virtual EVerbosityLevel | getVerbLevel () const |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT RCP< const ParameterList > | getValidVerboseObjectSublist () |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | setupVerboseObjectSublist (ParameterList *paramList) |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject) |
Public Member Functions inherited from Teuchos::VerboseObjectBase | |
virtual | ~VerboseObjectBase () |
VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const |
virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const |
virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) |
virtual RCP< FancyOStream > | getOStream () const |
virtual RCP< FancyOStream > | getOverridingOStream () const |
virtual std::string | getLinePrefix () const |
virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const |
Public Member Functions inherited from MueLu::Describable | |
virtual | ~Describable () |
Destructor. More... | |
virtual std::string | ShortClassName () const |
Return the class name of the object, without template parameters and without namespace. More... | |
virtual void | describe (Teuchos::FancyOStream &out_arg, const VerbLevel verbLevel=Default) const |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Print the object with some verbosity level to an FancyOStream object. More... | |
Public Member Functions inherited from Teuchos::Describable | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
LabeledObject () | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
virtual std::string | getObjectLabel () const |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Set/Get Methods for specific aggregation data | |
RCP< IndexManager_kokkos > & | GetIndexManagerKokkos () |
Get the index manager used by structured aggregation algorithms. This has to be done by the aggregation factory. More... | |
void | SetIndexManagerKokkos (RCP< IndexManager_kokkos > &geoDataKokkos) |
Set the index manager used by structured aggregation algorithms. This has to be done by the aggregation factory. More... | |
RCP< IndexManager > & | GetIndexManager () |
Get the index manager used by various aggregation algorithms. This has to be done by the aggregation factory. More... | |
void | SetIndexManager (RCP< IndexManager > &geoData) |
Set the index manager used by various aggregation algorithms. This has to be done by the aggregation factory. More... | |
colors_view_type & | GetGraphColors () |
Get a distance 2 coloring of the underlying graph. The coloring is computed and set during Phase1 of aggregation. More... | |
void | SetGraphColors (colors_view_type graphColors) |
Set a distance 2 coloring of the underlying graph. The coloring is computed and set during Phase1 of aggregation. More... | |
LO | GetGraphNumColors () |
Get the number of colors needed by the distance 2 coloring. More... | |
void | SetGraphNumColors (const LO graphNumColors) |
Set the number of colors needed by the distance 2 coloring. More... | |
Overridden from Teuchos::Describable | |
LO | numAggregates_ |
Number of aggregates on this processor. More... | |
GO | numGlobalAggregates_ |
Number of global aggregates. More... | |
RCP< LOMultiVector > | vertex2AggId_ |
RCP< LOVector > | procWinner_ |
RCP< IndexManager_kokkos > | geoDataKokkos_ |
RCP< IndexManager > | geoData_ |
colors_view_type | graphColors_ |
LO | graphNumColors_ |
Teuchos::ArrayRCP< bool > | isRoot_ |
An ArrayRCP of booleans specifying if a local entry is an aggregate root. More... | |
bool | aggregatesIncludeGhosts_ |
Set to false iff aggregates do not include any DOFs belong to other processes. More... | |
aggregates_sizes_type | aggregateSizes_ |
Array of sizes of each local aggregate. More... | |
aggregates_sizes_type::HostMirror | aggregateSizesHost_ |
local_graph_type | graph_ |
Aggregates represented as Kokkos graph type. More... | |
std::string | description () const |
Return a simple one-line description of this object. More... | |
void | print (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=verbLevel_default) const |
Print the object with some verbosity level to an FancyOStream object. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from MueLu::VerboseObject | |
static void | SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream) |
static void | SetMueLuOFileStream (const std::string &filename) |
static Teuchos::RCP < Teuchos::FancyOStream > | GetMueLuOStream () |
static void | SetDefaultVerbLevel (const VerbLevel defaultVerbLevel) |
Set the default (global) verbosity level. More... | |
static VerbLevel | GetDefaultVerbLevel () |
Get the default (global) verbosity level. More... | |
Static Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
static void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
Static Public Member Functions inherited from Teuchos::VerboseObjectBase | |
static void | setDefaultOStream (const RCP< FancyOStream > &defaultOStream) |
static RCP< FancyOStream > | getDefaultOStream () |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
Protected Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
Protected Member Functions inherited from Teuchos::VerboseObjectBase | |
void | initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual void | informUpdatedVerbosityState () const |
Container class for aggregation information.
Structure holding aggregate information. Right now, nAggregates, IsRoot, Vertex2AggId, procWinner are populated. This allows us to look at a node and determine the aggregate to which it has been assigned and the id of the processor that owns this aggregate. It is not so easy to determine vertices within the kth aggregate or the size of the kth aggregate. Thus, it might be useful to have a secondary structure which would be a rectangular CrsGraph where rows (or vertices) correspond to aggregates and colunmns (or edges) correspond to nodes. While not strictly necessary, it might be convenient.
Definition at line 69 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type = LocalOrdinal |
Definition at line 71 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type = GlobalOrdinal |
Definition at line 72 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::execution_space = typename Node::execution_space |
Definition at line 73 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::node_type = Node |
Definition at line 74 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::device_type = typename Node::device_type |
Definition at line 75 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::range_type = Kokkos::RangePolicy<local_ordinal_type, execution_space> |
Definition at line 76 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::LO_view = Kokkos::View<local_ordinal_type*, device_type> |
Definition at line 77 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::aggregates_sizes_type = Kokkos::View<LocalOrdinal*, device_type> |
Definition at line 79 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type = typename LWGraph_kokkos::local_graph_type |
Definition at line 87 of file MueLu_Aggregates_decl.hpp.
using MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::colors_view_type = Kokkos::View<typename local_graph_type::entries_type::data_type, typename local_graph_type::device_type::memory_space> |
Definition at line 89 of file MueLu_Aggregates_decl.hpp.
MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::Aggregates | ( | const LWGraph & | graph | ) |
Standard constructor for Aggregates structure.
Standard constructor of aggregates takes a Graph object as parameter. Uses the graph.GetImportMap() to initialize the internal vector for mapping nodes to (local) aggregate ids as well as the mapping of node to the owning processor id.
Definition at line 27 of file MueLu_Aggregates_def.hpp.
MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::Aggregates | ( | LWGraph_kokkos | graph | ) |
Standard constructor for Aggregates structure.
Standard constructor of aggregates takes a LWGraph object as parameter. Uses the graph.GetImportMap() to initialize the internal vector for mapping nodes to (local) aggregate ids as well as the mapping of node to the owning processor id.
Definition at line 45 of file MueLu_Aggregates_def.hpp.
MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::Aggregates | ( | const RCP< const Map > & | map | ) |
Constructor for Aggregates structure.
This constructor takes a RCP pointer to a map which is used for the internal mappings of nodes to the (local) aggregate ids and the owning processor.
Definition at line 63 of file MueLu_Aggregates_def.hpp.
|
inlinevirtual |
Destructor.
Definition at line 119 of file MueLu_Aggregates_decl.hpp.
|
inline |
Get the index manager used by structured aggregation algorithms. This has to be done by the aggregation factory.
Definition at line 127 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set the index manager used by structured aggregation algorithms. This has to be done by the aggregation factory.
Definition at line 132 of file MueLu_Aggregates_decl.hpp.
|
inline |
Get the index manager used by various aggregation algorithms. This has to be done by the aggregation factory.
Definition at line 137 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set the index manager used by various aggregation algorithms. This has to be done by the aggregation factory.
Definition at line 142 of file MueLu_Aggregates_decl.hpp.
|
inline |
Get a distance 2 coloring of the underlying graph. The coloring is computed and set during Phase1 of aggregation.
Definition at line 147 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set a distance 2 coloring of the underlying graph. The coloring is computed and set during Phase1 of aggregation.
Definition at line 152 of file MueLu_Aggregates_decl.hpp.
|
inline |
Get the number of colors needed by the distance 2 coloring.
Definition at line 156 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set the number of colors needed by the distance 2 coloring.
Definition at line 160 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set number of local aggregates on current processor.
This has to be done by the aggregation routines.
Definition at line 168 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set number of global aggregates on current processor.
This has to be done by the aggregation routines.returns the number of aggregates of the current processor. Note: could/should be renamed to GetNumLocalAggregates?
Definition at line 174 of file MueLu_Aggregates_decl.hpp.
|
inline |
Definition at line 177 of file MueLu_Aggregates_decl.hpp.
|
inline |
Record whether aggregates include DOFs from other processes.
Definition at line 182 of file MueLu_Aggregates_decl.hpp.
|
inline |
Return false if and only if no aggregates include DOFs from other processes.
Used in construction of tentative prolongator to skip a communication phase.
Definition at line 190 of file MueLu_Aggregates_decl.hpp.
|
inline |
Returns a nonconstant vector that maps local node IDs to local aggregates IDs.
For local node ID i, the corresponding vector entry v[i] is the local aggregate id to which i belongs on the current processor.
Definition at line 198 of file MueLu_Aggregates_decl.hpp.
|
inline |
Returns nonconstant vector that maps local node IDs to owning processor IDs.
For local node ID i, the corresponding vector entry v[i] is the owning processor ID.
Definition at line 204 of file MueLu_Aggregates_decl.hpp.
|
inline |
Returns constant vector that maps local node IDs to local aggregates IDs.
For local node ID i, the corresponding vector entry v[i] is the local aggregate id to which i belongs on the current processor.
Definition at line 209 of file MueLu_Aggregates_decl.hpp.
|
inline |
Returns constant vector that maps local node IDs to owning processor IDs.
For local node ID i, the corresponding vector entry v[i] is the owning processor ID.
Definition at line 215 of file MueLu_Aggregates_decl.hpp.
|
inline |
Returns true if node with given local node id is marked to be a root node.
Definition at line 218 of file MueLu_Aggregates_decl.hpp.
|
inline |
Set root node information.
Used by aggregation methods only.
Definition at line 224 of file MueLu_Aggregates_decl.hpp.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::GetMap | ( | ) | const |
returns (overlapping) map of aggregate/node distribution
Definition at line 278 of file MueLu_Aggregates_def.hpp.
Aggregates< LocalOrdinal, GlobalOrdinal, Node >::aggregates_sizes_type::const_type MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::ComputeAggregateSizes | ( | bool | forceRecompute = false | ) | const |
Compute sizes of aggregates.
Returns the number of nodes in each aggregate in an array. If the aggregate sizes are not stored internally (which is the default), they are computed and returned. If the aggregate sizes have been stored internally, then they are not recomputed, but instead the stored sizes are returned.
[in] | forceRecompute | if true, force recomputation of the aggregate sizes. |
Definition at line 81 of file MueLu_Aggregates_def.hpp.
Teuchos::ArrayRCP< LocalOrdinal > MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::ComputeAggregateSizesArrayRCP | ( | bool | forceRecompute = false | ) | const |
Compute sizes of aggregates.
Returns the number of nodes in each aggregate in an array. If the aggregate sizes are not stored internally (which is the default), they are computed and returned. If the aggregate sizes have been stored internally, then they are not recomputed, but instead the stored sizes are returned.
[in] | forceRecompute | if true, force recomputation of the aggregate sizes. |
Definition at line 111 of file MueLu_Aggregates_def.hpp.
Aggregates< LocalOrdinal, GlobalOrdinal, Node >::local_graph_type MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::GetGraph | ( | ) | const |
Definition at line 132 of file MueLu_Aggregates_def.hpp.
void MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::ComputeNodesInAggregate | ( | LO_view & | aggPtr, |
LO_view & | aggNodes, | ||
LO_view & | unaggregated | ||
) | const |
Generates a compressed list of nodes in each aggregate, where the entries in aggNodes[aggPtr[i]] up to aggNodes[aggPtr[i+1]-1] contain the nodes in aggregate i. unaggregated contains the list of nodes which are, for whatever reason, not aggregated (e.g. Dirichlet)
Definition at line 192 of file MueLu_Aggregates_def.hpp.
GlobalOrdinal MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::GetNumGlobalAggregatesComputeIfNeeded | ( | ) |
Get global number of aggregates.
Definition at line 266 of file MueLu_Aggregates_def.hpp.
|
virtual |
Return a simple one-line description of this object.
Reimplemented from MueLu::Describable.
Definition at line 246 of file MueLu_Aggregates_def.hpp.
void MueLu::Aggregates< LocalOrdinal, GlobalOrdinal, Node >::print | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = verbLevel_default |
||
) | const |
Print the object with some verbosity level to an FancyOStream object.
Definition at line 254 of file MueLu_Aggregates_def.hpp.
|
private |
Number of aggregates on this processor.
Definition at line 273 of file MueLu_Aggregates_decl.hpp.
|
private |
Number of global aggregates.
Definition at line 274 of file MueLu_Aggregates_decl.hpp.
|
private |
vertex2AggId[k] gives a local id corresponding to the aggregate to which local id k has been assigned. While k is the local id on my processor (MyPID), vertex2AggId[k] is the local id on the processor which actually owns the aggregate.
Definition at line 280 of file MueLu_Aggregates_decl.hpp.
|
private |
If k is the local id on my processor (MyPID), the owning processor has the id given by procWinner[k]
Definition at line 286 of file MueLu_Aggregates_decl.hpp.
|
private |
geoData stores an index manager object that is used to perform structured aggreation on a problem.
Definition at line 291 of file MueLu_Aggregates_decl.hpp.
|
private |
geoData stores an index manager object that is used to perform structured aggreation on a problem.
Definition at line 296 of file MueLu_Aggregates_decl.hpp.
|
private |
graphColors_ stores a view that assigns a color to each node in the graph These colors are used to parallelize the aggregation process in UncoupledAggregation
Definition at line 301 of file MueLu_Aggregates_decl.hpp.
|
private |
graphNumColors_ stores the number of colors that are needed to perform a distance 2 coloring of the underlying graph.
Definition at line 306 of file MueLu_Aggregates_decl.hpp.
|
private |
An ArrayRCP of booleans specifying if a local entry is an aggregate root.
Definition at line 309 of file MueLu_Aggregates_decl.hpp.
|
private |
Set to false iff aggregates do not include any DOFs belong to other processes.
Definition at line 312 of file MueLu_Aggregates_decl.hpp.
|
mutableprivate |
Array of sizes of each local aggregate.
Definition at line 315 of file MueLu_Aggregates_decl.hpp.
|
mutableprivate |
aggragateSizesHost_ is a host copy of aggregate sizes, which helps slightly reduce the cost of calling ComputeAggregateSizes from different parts of MueLu that require such data on the host device.
Definition at line 322 of file MueLu_Aggregates_decl.hpp.
|
mutableprivate |
Aggregates represented as Kokkos graph type.
Definition at line 325 of file MueLu_Aggregates_decl.hpp.