Public Member Functions | |
SubDomain (Domain &domain, PetscInt dsNumber, const std::vector< std::shared_ptr< FieldDescription >> &allAuxFields) | |
SubDomain (const SubDomain &temp_obj)=delete | |
SubDomain & | operator= (const SubDomain &temp_obj)=delete |
void | CreateSubDomainStructures () |
const Field & | GetField (const std::string &fieldName) const |
const Field & | GetField (PetscInt id, FieldLocation location=FieldLocation::SOL) const |
const std::vector< Field > & | GetFields (FieldLocation type=FieldLocation::SOL) const |
std::vector< Field > | GetFields (FieldLocation type, std::string_view tag) const |
DM | GetFieldDM (const Field &field) const noexcept |
Vec | GetVec (const Field &field) const noexcept |
Vec | GetGlobalVec (const Field &field) noexcept |
auto | GetSolutionAccessor (const Field &field) const |
auto | GetConstSolutionAccessor (const Field &field) const |
auto | GetAuxAccessor (const Field &field) const |
auto | GetConstAuxAccessor (const Field &field) const |
auto | GetSolutionAccessor (const std::string &fieldName) const |
auto | GetConstSolutionAccessor (const std::string &fieldName) const |
auto | GetAuxAccessor (const std::string &fieldName) |
auto | GetConstAuxAccessor (const std::string &fieldName) |
bool | ContainsField (const std::string &fieldName) |
bool | InRegion (const domain::Region &) const |
bool | InRegion (PetscInt point) const |
MPI_Comm | GetComm () const |
DMLabel | GetLabel () |
const std::string & | GetId () const override |
PetscInt | GetDimensions () const noexcept |
PetscObject | GetPetscFieldObject (const Field &field) |
DM & | GetDM () const noexcept |
DM | GetAuxDM () const noexcept |
Vec | GetSolutionVector () const noexcept |
Vec | GetAuxVector () const noexcept |
Vec | GetAuxGlobalVector () |
DM | GetSubDM () |
Vec | GetSubSolutionVector () |
DM | GetSubAuxDM () |
Vec | GetSubAuxVector () |
PetscDS | GetDiscreteSystem () |
PetscDS | GetAuxDiscreteSystem () |
PetscInt | GetNumberFields () const |
void | ProjectFieldFunctionsToLocalVector (const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> &fieldFunctions, Vec locVec, PetscReal time=0.0) const |
PetscErrorCode | ProjectFieldFunctionsToSubDM (const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> &initialization, Vec globVec, PetscReal time=0.0) |
void | SetsExactSolutions (const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> &exactSolutions) |
PetscErrorCode | GetFieldGlobalVector (const Field &, IS *vecIs, Vec *vec, DM *subdm) |
PetscErrorCode | RestoreFieldGlobalVector (const Field &, IS *vecIs, Vec *vec, DM *subdm) |
PetscErrorCode | GetFieldLocalVector (const Field &, PetscReal time, IS *vecIs, Vec *vec, DM *subdm) |
PetscErrorCode | RestoreFieldLocalVector (const Field &, IS *vecIs, Vec *vec, DM *subdm) |
PetscErrorCode | Save (PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override |
PetscErrorCode | Restore (PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override |
void | CreateEmptySubDM (DM *inDM, std::shared_ptr< domain::Region > region={}) |
void | GetRange (const std::shared_ptr< ablate::domain::Region > ®ion, PetscInt depth, ablate::domain::Range &range) const |
void | GetCellRange (const std::shared_ptr< ablate::domain::Region > ®ion, ablate::domain::Range &cellRange) const |
void | GetFaceRange (const std::shared_ptr< ablate::domain::Region > ®ion, ablate::domain::Range &faceRange) const |
void | RestoreRange (ablate::domain::Range &range) const |
Public Member Functions inherited from ablate::io::Serializable | |
virtual SerializerType | Serialize () const |
Additional Inherited Members | |
Public Types inherited from ablate::io::Serializable | |
enum class | SerializerType { none , collective , serial } |
Static Protected Member Functions inherited from ablate::io::Serializable | |
static PetscErrorCode | SaveKeyValue (PetscViewer viewer, const char *name, PetscScalar value) |
static PetscErrorCode | RestoreKeyValue (PetscViewer viewer, const char *name, PetscScalar &value) |
template<class T > | |
static PetscErrorCode | SaveKeyValue (PetscViewer viewer, const char *name, T value) |
template<class T > | |
static PetscErrorCode | RestoreKeyValue (PetscViewer viewer, const char *name, T &value) |
template<class T > | |
static SerializerType | DetermineSerializerType (const T &types) |
ablate::domain::SubDomain::SubDomain | ( | Domain & | domain, |
PetscInt | dsNumber, | ||
const std::vector< std::shared_ptr< FieldDescription >> & | allAuxFields | ||
) |
Create a subdomain based upon a domain
domain | |
dsNumber | ds number in the dm, |
allAuxFields | and a list of aux fields to create in this subdomain |
|
inline |
Determine if the field is available as either a sol or aux field
fieldName |
void ablate::domain::SubDomain::CreateEmptySubDM | ( | DM * | inDM, |
std::shared_ptr< domain::Region > | region = {} |
||
) |
This checks for whether the label describing the subdomain exists. If it does, use DMPlexFilter. If not, use DMClone to return new DM.
inDM |
void ablate::domain::SubDomain::CreateSubDomainStructures | ( | ) |
Create the auxDM, auxVec, and other structures on the subDomain
|
inline |
Returns the accessor for the specified aux field
field |
|
inline |
Returns the accessor for the specified aux field
field |
|
inline |
Get an aux DS if it is available
|
inlinenoexcept |
Returns the dm describing the aux fields living in this subdomain. The dm is defined across the entire mesh, but the fields are only define under this subdomain
Vec ablate::domain::SubDomain::GetAuxGlobalVector | ( | ) |
Return the global aux vector with information updated from the auxLocVec
|
inlinenoexcept |
Returns the local aux vector
|
inline |
Return the range of cells in the subDomain and region.
region | |
cellRange |
|
inline |
The comm used to define this subDomain and resulting solvers
|
inline |
Returns the accessor for the specified aux field
field |
|
inline |
Returns the accessor for the specified aux field
field |
|
inline |
Returns the accessor for the specified solution field by name
field |
|
inline |
Returns the accessor for the specified solution field
field |
|
inlinenoexcept |
Returns the number of physical dimensions defining the dm
|
inline |
Get the discreteSystem describe system for this subDomain
|
inlinenoexcept |
Returns raw access to the global dm
|
inline |
Return the range of faces/edges in the subDomain and region.
region | |
cellRange |
|
inline |
returns a references to the field (sol/aux) for a given field name
fieldName | the string name of the field |
|
inline |
returns a references to the field for a given id and location
id | |
location |
|
inlinenoexcept |
gets the dm corresponding to a field location (aux/sol)
field |
PetscErrorCode ablate::domain::SubDomain::GetFieldGlobalVector | ( | const Field & | field, |
IS * | vecIs, | ||
Vec * | vec, | ||
DM * | subdm | ||
) |
Get a global vector with only a single field
vecIs | |
vec | |
subdm |
PetscErrorCode ablate::domain::SubDomain::GetFieldLocalVector | ( | const Field & | field, |
PetscReal | time, | ||
IS * | vecIs, | ||
Vec * | vec, | ||
DM * | subdm | ||
) |
Get a local vector (with boundary values) with only a single field
vecIs | |
time | time is ued to insert boundary conditions for the global solution vector |
vec | |
subdm |
std::vector< ablate::domain::Field > ablate::domain::SubDomain::GetFields | ( | FieldLocation | type, |
std::string_view | tag | ||
) | const |
returns all fields of a certain type with a specific tag. This can be costly and should be used only for setup
type | |
tag |
|
inline |
returns all fields of a certain type
type | (defaults to SOL) |
|
inlinenoexcept |
Returns always returns a global vector of sol or aux
field |
|
inlineoverridevirtual |
|
inline |
The label (if any) used to define this subDomain
|
inline |
The number of aux and solution fields in this subdomain
PetscObject ablate::domain::SubDomain::GetPetscFieldObject | ( | const Field & | field | ) |
Get the petscField object from the dm or auxDm for this region
fieldName |
|
inline |
Return the range of DMPlex objects at a given depth in the subDomain and region.
region | |
depth | |
range |
|
inline |
Returns the accessor for the specified solution field by name
field |
|
inline |
Returns the accessor for the specified solution field
field |
|
inlinenoexcept |
Returns the global solution vector
DM ablate::domain::SubDomain::GetSubAuxDM | ( | ) |
The SubDM is defined as a dm that lives only over this subdomain region. It may be useful for outputting. This function will create (if needed) and return the subdomain for the aux variables. If a subdomain is not needed (same ds over entire dm) the global aux dm is returned.
Vec ablate::domain::SubDomain::GetSubAuxVector | ( | ) |
Returns a "global" aux vector defined over the subDomain
DM ablate::domain::SubDomain::GetSubDM | ( | ) |
The SubDM is defined as a dm that lives only over this subdomain region. It may be useful for outputting. This function will create (if needed) and return the subdomain. If a subdomain is not needed (same ds over entire dm) the global dm is returned.
Vec ablate::domain::SubDomain::GetSubSolutionVector | ( | ) |
Returns a "global" solution vector defined over the subDomain
|
inlinenoexcept |
Returns the global solution field or the local aux vector depending upon the field
field |
bool ablate::domain::SubDomain::InRegion | ( | const domain::Region & | region | ) | const |
Determine if the provided region lives inside of subdomain region
|
inline |
determines if this point is in this region as defined by the label and labelID
point |
void ablate::domain::SubDomain::ProjectFieldFunctionsToLocalVector | ( | const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> & | fieldFunctions, |
Vec | locVec, | ||
PetscReal | time = 0.0 |
||
) | const |
project the list of field function into the provided local vector. Allows solution and aux vectors
fieldFunctions | |
globVec |
PetscErrorCode ablate::domain::SubDomain::ProjectFieldFunctionsToSubDM | ( | const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> & | initialization, |
Vec | globVec, | ||
PetscReal | time = 0.0 |
||
) |
Support function to project the fields on to vector that lives only on the subDM
initialization | |
globVec | |
time |
|
overridevirtual |
PetscErrorCode ablate::domain::SubDomain::RestoreFieldGlobalVector | ( | const Field & | field, |
IS * | vecIs, | ||
Vec * | vec, | ||
DM * | subdm | ||
) |
Restore a global vector with only a single field
vecIs | |
vec | |
subdm |
PetscErrorCode ablate::domain::SubDomain::RestoreFieldLocalVector | ( | const Field & | field, |
IS * | vecIs, | ||
Vec * | vec, | ||
DM * | subdm | ||
) |
Restore a local vector with only a single field
vecIs | |
vec | |
subdm |
|
inline |
Restore the range of cells in the subDomain and region.
range |
|
overridevirtual |
void ablate::domain::SubDomain::SetsExactSolutions | ( | const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> & | exactSolutions | ) |
set exactSolutions if the fields live in the ds
exactSolutions |