19 #ifndef GEOSX_LINEARALGEBRA_DOFMANAGER_HPP_ 20 #define GEOSX_LINEARALGEBRA_DOFMANAGER_HPP_ 31 class DomainPartition;
33 class ObjectManagerBase;
34 class FluxApproximationBase;
163 void addField(
string const & fieldName,
173 void addField(
string const & fieldName,
184 void addField(
string const & fieldName,
213 void addField(
string const & fieldName,
226 string const & colFieldName,
238 string const & colFieldName,
251 string const & colFieldName,
253 bool const symmetric );
274 string const & colFieldName,
277 bool const symmetric );
321 string const &
getKey(
string const & fieldName )
const;
418 template<
typename MATRIX >
420 bool closePattern =
true )
const;
430 template<
typename MATRIX >
432 string const & rowFieldName,
433 string const & colFieldName,
434 bool closePattern =
true )
const;
449 string const & rowFieldName,
450 string const & colFieldName )
const;
466 template<
typename VECTOR >
468 string const & srcFieldName,
469 string const & dstFieldName,
470 real64 const scalingFactor,
488 string const & srcFieldName,
489 string const & dstFieldName,
490 real64 const scalingFactor,
508 template<
typename VECTOR >
510 string const & srcFieldName,
511 string const & dstFieldName,
512 real64 const scalingFactor,
530 string const & srcFieldName,
531 string const & dstFieldName,
532 real64 const scalingFactor,
550 template<
typename VECTOR >
552 string const & srcFieldName,
553 string const & dstFieldName,
554 real64 const scalingFactor,
573 string const & srcFieldName,
574 string const & dstFieldName,
575 real64 const scalingFactor,
593 template<
typename VECTOR >
595 string const & srcFieldName,
596 string const & dstFieldName,
597 real64 const scalingFactor,
616 string const & srcFieldName,
617 string const & dstFieldName,
618 real64 const scalingFactor,
643 std::vector< SubComponent >
644 filterDofs( std::vector< SubComponent >
const & excluded )
const;
658 template<
typename MATRIX >
659 void makeRestrictor( std::vector< SubComponent >
const & selection,
660 MPI_Comm
const & comm,
662 MATRIX & restrictor )
const;
676 void initializeDataStructure();
681 localIndex getFieldIndex(
string const & name )
const;
701 template<
typename MATRIX >
702 void setSparsityPatternOneBlock( MATRIX & pattern,
706 template<
typename MATRIX >
707 void setSparsityPatternFromStencil( MATRIX & pattern,
738 template<
int DIMS_PER_DOF >
757 template<
typename FIELD_OP,
typename POLICY,
typename LOCAL_VECTOR >
758 void vectorToField( LOCAL_VECTOR
const localVector,
759 string const & srcFieldName,
760 string const & dstFieldName,
761 real64 const scalingFactor,
780 template<
typename FIELD_OP,
typename POLICY,
typename LOCAL_VECTOR >
781 void fieldToVector( LOCAL_VECTOR localVector,
782 string const & srcFieldName,
783 string const & dstFieldName,
784 real64 const scalingFactor,
798 std::vector< FieldDescription > m_fields;
801 std::vector< std::vector< CouplingDescription > > m_coupling;
Connector
Enumeration of geometric objects for connectivity type. Note that this enum is nearly identical to Lo...
std::vector< SubComponent > filterDofs(std::vector< SubComponent > const &excluded) const
Create a dof selection by filtering out excluded components.
DofManager & operator=(DofManager const &)=delete
Deleted copy assignment.
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
constexpr void transpose(DST_MATRIX &&LVARRAY_RESTRICT_REF dstMatrix, SRC_MATRIX const &LVARRAY_RESTRICT_REF srcMatrix)
Store the transpose of the NxM matrix srcMatrix in dstMatrix.
location is face (like flux in mixed finite elements)
void makeRestrictor(std::vector< SubComponent > const &selection, MPI_Comm const &comm, bool transpose, MATRIX &restrictor) const
Create a matrix that restricts vectors and matrices to a subset of DOFs.
Contains the implementation of LvArray:SparsityPattern.
localIndex numLocalDofs(string const &fieldName="") const
Return local number of dofs on this processor. If field argument is empty, return monolithic size...
array1d< localIndex > numComponentsPerField() const
Return an array of number of components per field, sorted by field registration order.
Class facilitating the representation of a multi-level discretization of a MeshBody.
localIndex numComponents
number of vector components
location is element (like pressure in finite volumes)
Describes a selection of components from a DoF field.
void setMesh(DomainPartition &domain, localIndex const meshLevelIndex=0, localIndex const meshBodyIndex=0)
Assign a mesh.
location is edge (like flux between fracture elements)
globalIndex blockOffset
offset of this field's block in a block-wise ordered system
globalIndex numGlobalDofs(string const &fieldName="") const
Return global number of dofs across all processors. If field argument is empty, return monolithic siz...
void addField(string const &fieldName, Location const location)
Just an interface to allow only three parameters.
void addVectorToField(VECTOR const &vector, string const &srcFieldName, string const &dstFieldName, real64 const scalingFactor, localIndex const loCompIndex=0, localIndex const hiCompIndex=-1) const
Add values from LA vectors to simulation data arrays.
globalIndex numGlobalSupport(string const &fieldName) const
Get the local number of support points across all processors.
void setSparsityPattern(MATRIX &matrix, bool closePattern=true) const
Populate sparsity pattern of the entire system matrix.
~DofManager()=default
Destructor.
This class serves to provide a "view" of a multidimensional array.
Location
Enumeration of geometric objects for support location. Note that this enum is nearly identical to Con...
std::vector< std::string > regions
list of region names
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns, and generally simplifying the interaction between PhysicsSolvers and linear algebra operations.
localIndex hiComp
High component index (excluded from selection)
localIndex numLocalSupport(string const &fieldName) const
Get the local number of support points on this processor.
string key
string key for index array
double real64
64-bit floating point type.
string fieldName
Name of the DOF field in DofManager.
std::vector< string > regions
list of support region names
void copyFieldToVector(VECTOR &vector, string const &srcFieldName, string const &dstFieldName, real64 const scalingFactor, localIndex const loCompIndex=0, localIndex const hiCompIndex=-1) const
Copy values from nodes to DOFs.
string const & getKey(string const &fieldName) const
Return the key used to record the field in the DofManager.
there is no connectivity (self connected field, like a lumped mass matrix)
void addFieldToVector(VECTOR &vector, string const &srcFieldName, string const &dstFieldName, real64 const scalingFactor, localIndex const loCompIndex=0, localIndex const hiCompIndex=-1) const
Add values from a simulation data array to a DOF vector.
Location getLocation(string const &fieldName) const
Get the support location type of the field.
globalIndex globalOffset(string const &fieldName) const
Get global offset of field's block on current processor in the system matrix.
globalIndex rankOffset(string const &fieldName="") const
Return the sum of local dofs across all previous processors w.r.t. to the calling one for the specifi...
void addCoupling(string const &rowFieldName, string const &colFieldName, Connector const connectivity)
Just an interface to allow only three parameters.
localIndex loComp
Low component index (included in selection)
DofManager(string name)
Constructor.
globalIndex rankOffset
field's first DoF on current processor (within its block, ignoring other fields)
static localIndex constexpr MAX_FIELDS
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
bool fieldExists(string const &name) const
Check if string key is already being used.
array1d< localIndex > numLocalDofsPerField() const
Return an array of local number of dofs on this processor sorted by field registration order...
string docstring
documentation string
array1d< localIndex > getLocalDofComponentLabels() const
Computes an array of size equal to sum of all field local number of dofs containing unique integer la...
localIndex numLocalDof
number of local rows
void copyVectorToField(VECTOR const &vector, string const &srcFieldName, string const &dstFieldName, real64 const scalingFactor, localIndex const loCompIndex=0, localIndex const hiCompIndex=-1) const
Copy values from LA vectors to simulation data arrays.
globalIndex globalOffset
global offset of field's DOFs on current processor for multi-field problems
Location location
support location
localIndex numComponents(string const &fieldName="") const
Get the number of components in a field. If field argument is empty, return the total number of compo...
location is node (like displacements in finite elements)
void printFieldInfo(std::ostream &os=std::cout) const
Print the summary of declared fields and coupling.
This class provides a fixed dimensional resizeable array interface in addition to an interface simila...
void reorderByRank()
Finish populating fields and apply appropriate dof renumbering.
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
globalIndex numGlobalDof
number of global rows
No Schur complement - just block-GS/block-Jacobi preconditioner.