GEOS
|
#include "common/DataTypes.hpp"
#include "linearAlgebra/multiscale/mesh/MeshObjectManager.hpp"
#include "mesh/ObjectManagerBase.hpp"
Go to the source code of this file.
Classes | |
class | geos::multiscale::meshUtils::ScopedDataRegistrar |
Utility class to manage registration of temporary scope-bound data in the data repository. More... | |
Namespaces | |
geos | |
Functions | |
template<typename T , typename U = T> | |
void | geos::multiscale::meshUtils::mapIndexArray (arrayView1d< T const > const &src, arrayView1d< U const > const &map, array1d< U > &dst) |
Transform an array of mesh indices using a map. More... | |
template<typename T , typename U = T> | |
void | geos::multiscale::meshUtils::mapIndexArrayUnique (arrayView1d< T const > const &src, arrayView1d< U const > const &map, array1d< U > &dst) |
Transform an array of mesh indices using a map, removing duplicates. More... | |
template<typename T , typename U = T> | |
void | geos::multiscale::meshUtils::mapIndexSet (SortedArrayView< T const > const &src, arrayView1d< U const > const &map, SortedArray< U > &dst) |
Transform a set of mesh indices using a map. More... | |
template<typename POLICY , typename T , int NDIM, int USD, typename INDEX , typename FUNC > | |
void | geos::multiscale::meshUtils::fillArrayByDstIndex (arrayView1d< INDEX const > const &map, ArrayView< T, NDIM, USD > const &dst, FUNC &&src) |
Populate an array using a function source and an index map. More... | |
template<typename POLICY , typename T , int NDIM, int USD, typename INDEX , typename FUNC > | |
void | geos::multiscale::meshUtils::fillArrayBySrcIndex (arrayView1d< INDEX const > const &map, ArrayView< T, NDIM, USD > const &dst, FUNC &&src) |
Populate an array using a function source and an index map. More... | |
template<typename FUNC > | |
void | geos::multiscale::meshUtils::copyNeighborData (ObjectManagerBase const &srcManager, string const &mapKey, std::vector< integer > const &ranks, ObjectManagerBase &dstManager, FUNC &©Func) |
Sets up NeighborData of a new object manager mapped from existing one. More... | |
void | geos::multiscale::meshUtils::copySets (ObjectManagerBase const &srcManager, string const &mapKey, ObjectManagerBase &dstManager) |
Copy and update sets from an existing to a new object manager. More... | |
template<typename ITER , typename FUNC > | |
void | geos::multiscale::meshUtils::forUniqueValues (ITER first, ITER const last, FUNC &&func) |
Call a function on unique values from a previously collected range. More... | |
template<integer MAX_NEIGHBORS, typename L2C_MAP , typename C2L_MAP , typename PRED , typename FUNC > | |
void | geos::multiscale::meshUtils::forUniqueNeighbors (localIndex const locIdx, L2C_MAP const &locToConn, C2L_MAP const &connToLoc, PRED &&connPred, FUNC &&func) |
Call a function on unique indices of topological neighbors visited through location-connector adjacency maps. More... | |
template<integer MAX_NEIGHBORS, typename L2C_MAP , typename C2L_MAP , typename FUNC > | |
void | geos::multiscale::meshUtils::forUniqueNeighbors (localIndex const locIdx, L2C_MAP const &locToConn, C2L_MAP const &connToLoc, FUNC &&func) |
Call a function on unique indices of topological neighbors visited through location-connector adjacency maps. More... | |
template<integer MAX_NEIGHBORS, typename NBR_MAP , typename VAL_FUNC , typename VAL_PRED , typename FUNC > | |
void | geos::multiscale::meshUtils::forUniqueNeighborValues (localIndex const locIdx, NBR_MAP const &neighbors, VAL_FUNC &&valueFunc, VAL_PRED &&pred, FUNC &&func) |
Call a function on unique values of topological neighbors visited through location-connector adjacency maps. More... | |
template<integer MAX_NEIGHBORS, typename NBR_MAP , typename VAL_FUNC , typename FUNC > | |
void | geos::multiscale::meshUtils::forUniqueNeighborValues (localIndex const locIdx, NBR_MAP const &neighbors, VAL_FUNC &&valueFunc, FUNC &&func) |
Call a function on unique values of topological neighbors visited through location-connector adjacency maps. More... | |
template<typename INDEX_TYPE > | |
ArrayOfSets< INDEX_TYPE > | geos::multiscale::meshUtils::buildFineObjectToSubdomainMap (MeshObjectManager const &fineObjectManager, arrayView1d< INDEX_TYPE const > const &subdomains, string_array const &boundaryObjectSets) |
Build a map from mesh objects (nodes/cells) to subdomains defined by a partitioning of dual objects (cells/nodes). More... | |
template<typename INDEX_TYPE > | |
ArrayOfSets< INDEX_TYPE > | geos::multiscale::meshUtils::addBoundarySubdomains (MeshObjectManager const &manager, ArrayOfSetsView< INDEX_TYPE const > const &inputMap, string_array const &boundaryObjectSets) |
Insert virtual boundary subdomains into an adjacency map. More... | |
array1d< localIndex > | geos::multiscale::meshUtils::findCoarseNodesByDualPartition (MeshObjectManager::MapViewConst const &nodeToDual, MeshObjectManager::MapViewConst const &dualToNode, ArrayOfSetsView< globalIndex const > const &nodeToSubdomain, integer const minSubdomains, bool allowMultiNodes=true) |
Find "coarse" nodes in a mesh in which dual objects have been partitioned into subdomains. More... | |
ArrayOfSets< INDEX_TYPE > geos::multiscale::meshUtils::addBoundarySubdomains | ( | MeshObjectManager const & | manager, |
ArrayOfSetsView< INDEX_TYPE const > const & | inputMap, | ||
string_array const & | boundaryObjectSets | ||
) |
Insert virtual boundary subdomains into an adjacency map.
INDEX_TYPE | type of index in the adjacency map |
manager | multiscale mesh object manager |
inputMap | source map |
boundaryObjectSets | names of boundary object sets |
Boundary (virtual) subdomains are assigned negative indices [-1 ... boundaryObjectSets.size()] in the map. For example, if the input map contains cells [1,4,6] for a given node, and the node is found in the first and third boundary sets (in order listed), the output map will contain entries [-3,-1,1,4,6].
Definition at line 541 of file MeshUtils.hpp.
ArrayOfSets< INDEX_TYPE > geos::multiscale::meshUtils::buildFineObjectToSubdomainMap | ( | MeshObjectManager const & | fineObjectManager, |
arrayView1d< INDEX_TYPE const > const & | subdomains, | ||
string_array const & | boundaryObjectSets | ||
) |
Build a map from mesh objects (nodes/cells) to subdomains defined by a partitioning of dual objects (cells/nodes).
INDEX_TYPE | type of subdomain index |
fineObjectManager | |
subdomains | array of subdomain indices of dual objects |
boundaryObjectSets | (optional) list of boundary set names |
Boundaries (if present) are treated as additional "virtual" subdomains. They are assigned unique negative indices to distinguish them from actual subdomains. Pass an empty list of set names to ignore this feature.
Definition at line 473 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::copyNeighborData | ( | ObjectManagerBase const & | srcManager, |
string const & | mapKey, | ||
std::vector< integer > const & | ranks, | ||
ObjectManagerBase & | dstManager, | ||
FUNC && | copyFunc | ||
) |
Sets up NeighborData of a new object manager mapped from existing one.
FUNC | type of copy utility function |
srcManager | source object manager |
mapKey | key to extract the index map from source manager |
ranks | ranks of neighbors |
dstManager | target object manager |
copyFunc | copy utility function to use on each array |
Definition at line 248 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::copySets | ( | ObjectManagerBase const & | srcManager, |
string const & | mapKey, | ||
ObjectManagerBase & | dstManager | ||
) |
Copy and update sets from an existing to a new object manager.
srcManager | source object manager |
mapKey | key to extract the index map from source manager |
dstManager | target object manager |
void geos::multiscale::meshUtils::fillArrayByDstIndex | ( | arrayView1d< INDEX const > const & | map, |
ArrayView< T, NDIM, USD > const & | dst, | ||
FUNC && | src | ||
) |
Populate an array using a function source and an index map.
POLICY | loop execution policy |
T | type of value |
NDIM | number of array dimensions |
USD | array unit stride dimension |
INDEX | array index type |
FUNC | type of source function |
map | index map from source to destination |
dst | destination array |
src | source function called with source indices |
Elements for which index in the map is negative are ignored.
Definition at line 188 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::fillArrayBySrcIndex | ( | arrayView1d< INDEX const > const & | map, |
ArrayView< T, NDIM, USD > const & | dst, | ||
FUNC && | src | ||
) |
Populate an array using a function source and an index map.
POLICY | loop execution policy |
T | type of value |
NDIM | number of array dimensions |
USD | array unit stride dimension |
INDEX | array index type |
FUNC | type of source function |
map | index map from destination to source |
dst | destination array |
src | source function called with source indices |
Elements for which index in the map is negative are ignored.
Definition at line 220 of file MeshUtils.hpp.
array1d< localIndex > geos::multiscale::meshUtils::findCoarseNodesByDualPartition | ( | MeshObjectManager::MapViewConst const & | nodeToDual, |
MeshObjectManager::MapViewConst const & | dualToNode, | ||
ArrayOfSetsView< globalIndex const > const & | nodeToSubdomain, | ||
integer const | minSubdomains, | ||
bool | allowMultiNodes = true |
||
) |
Find "coarse" nodes in a mesh in which dual objects have been partitioned into subdomains.
nodeToDual | adjacency map of nodes to its dual object |
dualToNode | adjacency map of dual objects to nodes |
nodeToSubdomain | adjacency of nodes to subdomains (must be built by the caller) |
minSubdomains | minimum number of adjacent subdomains required for a coarse node (reduces search space) |
allowMultiNodes | whether clusters of similar nodes should produce multiple coarse nodes (if false, just one will be chosen) |
"Coarse" nodes are defined as those that, given a partitioning of dual objects (e.g. mesh volumes or cells) into subdomains, are adjacent to the locally maximal (among its neighbors) number of such subdomains.
The implementation groups nodes into "features" (that are groups of nodes sharing the same subdomain list) and finds features with largest list (or "key") among adjacent features. Such features typically consist of just a single node. Occasionally, in complex topologies, coarse features with multiple nodes may be found. In that case, any of the nodes in the feature can be chosen as coarse (since they all have the same adjacent subdomains), or each of them can be made a separate coarse node (this is the default behavior).
The analysis is fully local. In order for coarse nodes to be chosen consistently across processes, adjacency maps must include ghosted node/dual objects and correct global subdomain indices.
void geos::multiscale::meshUtils::forUniqueNeighbors | ( | localIndex const | locIdx, |
L2C_MAP const & | locToConn, | ||
C2L_MAP const & | connToLoc, | ||
FUNC && | func | ||
) |
Call a function on unique indices of topological neighbors visited through location-connector adjacency maps.
MAX_NEIGHBORS | maximum number of total non-unique neighbor indices |
L2C_MAP | type of location to connector map |
C2L_MAP | type of connector to location map |
FUNC | type of function to call |
locIdx | location index |
locToConn | location to connector map |
connToLoc | connector to location map |
func | function to call |
Version that does not use a connector predicate.
Definition at line 386 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::forUniqueNeighbors | ( | localIndex const | locIdx, |
L2C_MAP const & | locToConn, | ||
C2L_MAP const & | connToLoc, | ||
PRED && | connPred, | ||
FUNC && | func | ||
) |
Call a function on unique indices of topological neighbors visited through location-connector adjacency maps.
MAX_NEIGHBORS | maximum number of total non-unique neighbor indices |
L2C_MAP | type of location to connector map |
C2L_MAP | type of connector to location map |
PRED | type of predicate function |
FUNC | type of function to call |
locIdx | location index |
locToConn | location to connector map |
connToLoc | connector to location map |
connPred | predicate used to filter connector indices |
func | function to call |
Definition at line 350 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::forUniqueNeighborValues | ( | localIndex const | locIdx, |
NBR_MAP const & | neighbors, | ||
VAL_FUNC && | valueFunc, | ||
FUNC && | func | ||
) |
Call a function on unique values of topological neighbors visited through location-connector adjacency maps.
MAX_NEIGHBORS | maximum number of total non-unique neighbor indices |
NBR_MAP | type of connectivity map |
VAL_FUNC | type of value function |
FUNC | type of target function |
locIdx | location index |
neighbors | neighbor map |
valueFunc | function called with neighbor indices that returns target values |
func | function to call |
Version that does not use a value predicate.
Definition at line 447 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::forUniqueNeighborValues | ( | localIndex const | locIdx, |
NBR_MAP const & | neighbors, | ||
VAL_FUNC && | valueFunc, | ||
VAL_PRED && | pred, | ||
FUNC && | func | ||
) |
Call a function on unique values of topological neighbors visited through location-connector adjacency maps.
MAX_NEIGHBORS | maximum number of total non-unique neighbor indices |
NBR_MAP | type of connectivity map |
VAL_FUNC | type of value function |
VAL_PRED | type of value predicate |
FUNC | type of target function |
locIdx | location index |
neighbors | neighbor map |
valueFunc | function called with neighbor indices that returns target values |
pred | predicate used to filter values |
func | function to call |
Definition at line 412 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::forUniqueValues | ( | ITER | first, |
ITER const | last, | ||
FUNC && | func | ||
) |
Call a function on unique values from a previously collected range.
ITER | type of range iterator |
FUNC | type of function to call |
first | start of the range |
last | end of the range |
func | the function to call |
ITER
must not be a const iterator. Definition at line 321 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::mapIndexArray | ( | arrayView1d< T const > const & | src, |
arrayView1d< U const > const & | map, | ||
array1d< U > & | dst | ||
) |
Transform an array of mesh indices using a map.
T | source index type |
U | destination index type |
src | source index array |
map | index map |
dst | destination array |
This function assigns dst[i] = map[src[i]] for each i s.t. map[src[i]] >= 0. Thus entries with negative indices in the map are filtered out.
Definition at line 98 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::mapIndexArrayUnique | ( | arrayView1d< T const > const & | src, |
arrayView1d< U const > const & | map, | ||
array1d< U > & | dst | ||
) |
Transform an array of mesh indices using a map, removing duplicates.
T | source index type |
U | destination index type |
src | source index array |
map | index map |
dst | destination array |
This function assigns dst[i] = map[src[i]] for each i s.t. map[src[i]] >= 0. Thus entries with negative indices in the map are filtered out. In addition, destination array only contains unique values.
Definition at line 127 of file MeshUtils.hpp.
void geos::multiscale::meshUtils::mapIndexSet | ( | SortedArrayView< T const > const & | src, |
arrayView1d< U const > const & | map, | ||
SortedArray< U > & | dst | ||
) |
Transform a set of mesh indices using a map.
T | source index type |
U | destination index type |
src | source index set |
map | index map |
dst | destination set |
Definition at line 158 of file MeshUtils.hpp.