GEOSX
DofManager.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 Total, S.A
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
19 #ifndef GEOSX_LINEARALGEBRA_DOFMANAGER_HPP_
20 #define GEOSX_LINEARALGEBRA_DOFMANAGER_HPP_
21 
23 #include "common/DataTypes.hpp"
26 #include "mesh/NodeManager.hpp"
27 
28 namespace geosx
29 {
30 
31 class DomainPartition;
32 class MeshLevel;
33 class ObjectManagerBase;
34 class FluxApproximationBase;
35 
43 {
44 public:
45 
51  enum class Location
52  {
53  Elem,
54  Face,
55  Edge,
56  Node
57  };
58 
64  enum class Connector
65  {
66  Elem,
67  Face,
68  Edge,
69  Node,
70  None,
71  Stencil
72  };
73 
77  static localIndex constexpr MAX_FIELDS = 10;
78 
83  {
84  string name;
86  std::vector< string > regions;
88  string key;
89  string docstring;
95  };
96 
101  {
103  std::vector< std::string > regions;
104  FluxApproximationBase const * stencils = nullptr;
105  };
106 
112  explicit DofManager( string name );
113 
117  DofManager( DofManager const & ) = delete;
118 
122  DofManager( DofManager && ) = default;
123 
128  DofManager & operator=( DofManager const & ) = delete;
129 
134  DofManager & operator=( DofManager && ) = delete;
135 
139  ~DofManager() = default;
140 
144  void clear();
145 
153  void setMesh( DomainPartition & domain,
154  localIndex const meshLevelIndex = 0,
155  localIndex const meshBodyIndex = 0 );
156 
163  void addField( string const & fieldName,
164  Location const location );
165 
173  void addField( string const & fieldName,
174  Location const location,
175  localIndex const components );
176 
184  void addField( string const & fieldName,
185  Location const location,
186  arrayView1d< string const > const & regions );
187 
213  void addField( string const & fieldName,
214  Location const location,
215  localIndex const components,
216  arrayView1d< string const > const & regions );
217 
225  void addCoupling( string const & rowFieldName,
226  string const & colFieldName,
227  Connector const connectivity );
228 
237  void addCoupling( string const & rowFieldName,
238  string const & colFieldName,
239  Connector const connectivity,
240  arrayView1d< string const > const & regions );
241 
250  void addCoupling( string const & rowFieldName,
251  string const & colFieldName,
252  Connector const connectivity,
253  bool const symmetric );
254 
273  void addCoupling( string const & rowFieldName,
274  string const & colFieldName,
275  Connector const connectivity,
276  arrayView1d< string const > const & regions,
277  bool const symmetric );
278 
287  void addCoupling( string const & fieldName,
288  FluxApproximationBase const & stencils );
289 
305  void reorderByRank();
306 
313  bool fieldExists( string const & name ) const;
314 
321  string const & getKey( string const & fieldName ) const;
322 
330  globalIndex numGlobalDofs( string const & fieldName = "" ) const;
331 
339  localIndex numLocalDofs( string const & fieldName = "" ) const;
340 
348 
349 
357 
365  globalIndex rankOffset( string const & fieldName = "" ) const;
366 
374  localIndex numComponents( string const & fieldName = "" ) const;
375 
383 
389  localIndex numLocalSupport( string const & fieldName ) const;
390 
396  globalIndex numGlobalSupport( string const & fieldName ) const;
397 
403  Location getLocation( string const & fieldName ) const;
404 
410  globalIndex globalOffset( string const & fieldName ) const;
411 
418  template< typename MATRIX >
419  void setSparsityPattern( MATRIX & matrix,
420  bool closePattern = true ) const;
421 
430  template< typename MATRIX >
431  void setSparsityPattern( MATRIX & matrix,
432  string const & rowFieldName,
433  string const & colFieldName,
434  bool closePattern = true ) const;
435 
440  void setSparsityPattern( SparsityPattern< globalIndex > & pattern ) const;
441 
449  string const & rowFieldName,
450  string const & colFieldName ) const;
451 
466  template< typename VECTOR >
467  void copyVectorToField( VECTOR const & vector,
468  string const & srcFieldName,
469  string const & dstFieldName,
470  real64 const scalingFactor,
471  localIndex const loCompIndex = 0,
472  localIndex const hiCompIndex = -1 ) const;
473 
487  void copyVectorToField( arrayView1d< real64 const > const & localVector,
488  string const & srcFieldName,
489  string const & dstFieldName,
490  real64 const scalingFactor,
491  localIndex const loCompIndex = 0,
492  localIndex const hiCompIndex = -1 ) const;
493 
508  template< typename VECTOR >
509  void addVectorToField( VECTOR const & vector,
510  string const & srcFieldName,
511  string const & dstFieldName,
512  real64 const scalingFactor,
513  localIndex const loCompIndex = 0,
514  localIndex const hiCompIndex = -1 ) const;
515 
529  void addVectorToField( arrayView1d< real64 const > const & localVector,
530  string const & srcFieldName,
531  string const & dstFieldName,
532  real64 const scalingFactor,
533  localIndex const loCompIndex = 0,
534  localIndex const hiCompIndex = -1 ) const;
535 
550  template< typename VECTOR >
551  void copyFieldToVector( VECTOR & vector,
552  string const & srcFieldName,
553  string const & dstFieldName,
554  real64 const scalingFactor,
555  localIndex const loCompIndex = 0,
556  localIndex const hiCompIndex = -1 ) const;
557 
572  void copyFieldToVector( arrayView1d< real64 > const & localVector,
573  string const & srcFieldName,
574  string const & dstFieldName,
575  real64 const scalingFactor,
576  localIndex const loCompIndex = 0,
577  localIndex const hiCompIndex = -1 ) const;
578 
593  template< typename VECTOR >
594  void addFieldToVector( VECTOR & vector,
595  string const & srcFieldName,
596  string const & dstFieldName,
597  real64 const scalingFactor,
598  localIndex const loCompIndex = 0,
599  localIndex const hiCompIndex = -1 ) const;
600 
615  void addFieldToVector( arrayView1d< real64 > const & localVector,
616  string const & srcFieldName,
617  string const & dstFieldName,
618  real64 const scalingFactor,
619  localIndex const loCompIndex = 0,
620  localIndex const hiCompIndex = -1 ) const;
621 
628  {
629  string fieldName;
632  };
633 
643  std::vector< SubComponent >
644  filterDofs( std::vector< SubComponent > const & excluded ) const;
645 
658  template< typename MATRIX >
659  void makeRestrictor( std::vector< SubComponent > const & selection,
660  MPI_Comm const & comm,
661  bool transpose,
662  MATRIX & restrictor ) const;
663 
669  void printFieldInfo( std::ostream & os = std::cout ) const;
670 
671 private:
672 
676  void initializeDataStructure();
677 
681  localIndex getFieldIndex( string const & name ) const;
682 
686  void createIndexArray( FieldDescription & field );
687 
691  void removeIndexArray( FieldDescription const & field );
692 
701  template< typename MATRIX >
702  void setSparsityPatternOneBlock( MATRIX & pattern,
703  localIndex const rowFieldIndex,
704  localIndex const colFieldIndex ) const;
705 
706  template< typename MATRIX >
707  void setSparsityPatternFromStencil( MATRIX & pattern,
708  localIndex const fieldIndex ) const;
709 
716  void countRowLengthsOneBlock( arrayView1d< localIndex > const & rowLengths,
717  localIndex const rowFieldIndex,
718  localIndex const colFieldIndex ) const;
719 
720  void countRowLengthsFromStencil( arrayView1d< localIndex > const & rowLengths,
721  localIndex const fieldIndex ) const;
722 
731  void setSparsityPatternOneBlock( SparsityPattern< globalIndex > & pattern,
732  localIndex const rowFieldIndex,
733  localIndex const colFieldIndex ) const;
734 
735  void setSparsityPatternFromStencil( SparsityPattern< globalIndex > & pattern,
736  localIndex const fieldIndex ) const;
737 
738  template< int DIMS_PER_DOF >
739  void setFiniteElementSparsityPattern( SparsityPattern< globalIndex > & pattern,
740  localIndex const fieldIndex ) const;
741 
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,
762  localIndex const loCompIndex,
763  localIndex const hiCompIndex ) const;
764 
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,
785  localIndex const loCompIndex,
786  localIndex const hiCompIndex ) const;
787 
789  string m_name;
790 
792  DomainPartition * m_domain = nullptr;
793 
795  MeshLevel * m_mesh = nullptr;
796 
798  std::vector< FieldDescription > m_fields;
799 
801  std::vector< std::vector< CouplingDescription > > m_coupling;
802 
804  bool m_reordered;
805 };
806 
807 } /* namespace geosx */
808 
809 #endif /*GEOSX_LINEARALGEBRA_DOFMANAGER_HPP_*/
Connector
Enumeration of geometric objects for connectivity type. Note that this enum is nearly identical to Lo...
Definition: DofManager.hpp:64
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).
Definition: DataTypes.hpp:128
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.
Definition: MeshLevel.hpp:38
localIndex numComponents
number of vector components
Definition: DofManager.hpp:87
location is element (like pressure in finite volumes)
Describes a selection of components from a DoF field.
Definition: DofManager.hpp:627
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&#39;s block in a block-wise ordered system
Definition: DofManager.hpp:92
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.
Definition: ArrayView.hpp:67
Location
Enumeration of geometric objects for support location. Note that this enum is nearly identical to Con...
Definition: DofManager.hpp:51
std::vector< std::string > regions
list of region names
Definition: DofManager.hpp:103
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns, and generally simplifying the interaction between PhysicsSolvers and linear algebra operations.
Definition: DofManager.hpp:42
localIndex hiComp
High component index (excluded from selection)
Definition: DofManager.hpp:631
localIndex numLocalSupport(string const &fieldName) const
Get the local number of support points on this processor.
string key
string key for index array
Definition: DofManager.hpp:88
double real64
64-bit floating point type.
Definition: DataTypes.hpp:136
string fieldName
Name of the DOF field in DofManager.
Definition: DofManager.hpp:629
std::vector< string > regions
list of support region names
Definition: DofManager.hpp:86
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&#39;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)
Definition: DofManager.hpp:630
DofManager(string name)
Constructor.
globalIndex rankOffset
field&#39;s first DoF on current processor (within its block, ignoring other fields)
Definition: DofManager.hpp:93
static localIndex constexpr MAX_FIELDS
Definition: DofManager.hpp:77
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
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
Definition: DofManager.hpp:89
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
Definition: DofManager.hpp:90
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&#39;s DOFs on current processor for multi-field problems
Definition: DofManager.hpp:94
Location location
support location
Definition: DofManager.hpp:85
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...
Definition: Array.hpp:55
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
Definition: DofManager.hpp:91
No Schur complement - just block-GS/block-Jacobi preconditioner.