KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
|
#include <partitioned_fsi_utilities.hpp>
Public Types | |
typedef TSpace::VectorType | VectorType |
typedef TSpace::MatrixType | MatrixType |
typedef TSpace::VectorPointerType | VectorPointerType |
typedef TSpace::MatrixPointerType | MatrixPointerType |
Public Member Functions | |
KRATOS_CLASS_POINTER_DEFINITION (PartitionedFSIUtilities) | |
PartitionedFSIUtilities () | |
PartitionedFSIUtilities (const PartitionedFSIUtilities &Other)=delete | |
virtual | ~PartitionedFSIUtilities ()=default |
Public Operators | |
void | CreateCouplingSkin (const ModelPart &rOriginInterfaceModelPart, ModelPart &rDestinationInterfaceModelPart) |
Create a coupling condition-based skin object This method creates an condition-based skin model part by copying the conditions of a given skin model part. More... | |
int | GetInterfaceResidualSize (ModelPart &rInterfaceModelPart) |
double | GetInterfaceArea (ModelPart &rInterfaceModelPart) |
virtual VectorPointerType | SetUpInterfaceVector (ModelPart &rInterfaceModelPart) |
void | InitializeInterfaceVector (const ModelPart &rInterfaceModelPart, const Variable< TValueType > &rOriginVariable, VectorType &rInterfaceVector) |
virtual void | ComputeInterfaceResidualVector (ModelPart &rInterfaceModelPart, const Variable< TValueType > &rOriginalVariable, const Variable< TValueType > &rModifiedVariable, const Variable< TValueType > &rResidualVariable, VectorType &rInterfaceResidual, const std::string ResidualType="nodal", const Variable< double > &rResidualNormVariable=FSI_INTERFACE_RESIDUAL_NORM) |
Compute array variable interface residual vector This function computes (and stores in a vector) the residual of a vector variable over the fluid interface. The residual is defined as the OriginalVariable value minus the ModifiedVariable value. The nodal values of the residual are stored in the FSI_INTERFACE_RESIDUAL variable. Besides, the norm of the residual vector is stored in the ProcessInfo. More... | |
double | ComputeInterfaceResidualNorm (ModelPart &rInterfaceModelPart, const Variable< TValueType > &rOriginalVariable, const Variable< TValueType > &rModifiedVariable, const Variable< TValueType > &rResidualVariable, const std::string ResidualType="nodal") |
Computes the interface residual norm This method computes the interface residual norm. To do that it firstly computes the residual vector as is done in ComputeInterfaceResidualVector. More... | |
virtual void | AuxSetLocalValue (VectorType &rInterfaceResidual, const double &rResidualValue, const int AuxPosition) |
Auxiliar call to SetLocalValue for double type This method serves as an auxiliar call to the SetLocalValue function for double variables. More... | |
virtual void | AuxSetLocalValue (VectorType &rInterfaceResidual, const array_1d< double, 3 > &rResidualValue, const int AuxPosition) |
Auxiliar call to SetLocalValue for array type This method serves as an auxiliar call to the SetLocalValue function for array variables. More... | |
virtual void | UpdateInterfaceValues (ModelPart &rInterfaceModelPart, const Variable< TValueType > &rSolutionVariable, const VectorType &rCorrectedGuess) |
void | UpdateInterfaceLocalValue (const VectorType &rCorrectedGuess, double &rValueToUpdate, const int AuxPosition) |
Corrected guess local double value update This auxiliar method helps to update the nodal database using the values from a corrected guess vector values. More... | |
void | UpdateInterfaceLocalValue (const VectorType &rCorrectedGuess, array_1d< double, 3 > &rValueToUpdate, const int AuxPosition) |
Corrected guess local array value update This auxiliar method helps to update the nodal database using the values from a corrected guess vector values. Note that it performs the array components loop. More... | |
virtual void | ComputeAndPrintFluidInterfaceNorms (ModelPart &rInterfaceModelPart) |
virtual void | ComputeAndPrintStructureInterfaceNorms (ModelPart &rInterfaceModelPart) |
virtual void | CheckCurrentCoordinatesFluid (ModelPart &rModelPart, const double tolerance) |
virtual void | CheckCurrentCoordinatesStructure (ModelPart &rModelPart, const double tolerance) |
void | EmbeddedPressureToPositiveFacePressureInterpolator (ModelPart &rFluidModelPart, ModelPart &rStructureSkinModelPart) |
Pressure to positive face pressure interpolator This method interpolates the pressure in the embedded fluid background mesh to the positive face pressure of all nodes in the provided structure interface model part. More... | |
void | CalculateTractionFromPressureValues (ModelPart &rModelPart, const Variable< double > &rPressureVariable, const Variable< array_1d< double, 3 >> &rTractionVariable, const bool SwapTractionSign) |
void | CalculateTractionFromPressureValues (ModelPart &rModelPart, const Variable< double > &rPositivePressureVariable, const Variable< double > &rNegativePressureVariable, const Variable< array_1d< double, 3 >> &rTractionVariable, const bool SwapTractionSign) |
Protected Member Functions | |
Protected Operations | |
std::string | GetSkinElementName () |
Get the skin element name Auxiliary method that returns the auxiliary embedded skin element type name. More... | |
std::string | GetSkinConditionName () |
Get the skin condition name Auxiliary method that returns the auxiliary embedded skin condition type name. More... | |
void | ComputeConsistentResidual (ModelPart &rInterfaceModelPart, const Variable< TValueType > &rOriginalVariable, const Variable< TValueType > &rModifiedVariable, const Variable< TValueType > &rErrorStorageVariable) |
void | ComputeNodeByNodeResidual (ModelPart &rInterfaceModelPart, const Variable< TValueType > &rOriginalVariable, const Variable< TValueType > &rModifiedVariable, const Variable< TValueType > &rErrorStorageVariable) |
virtual void | SetLocalValue (VectorType &rVector, int LocalRow, double Value) const |
virtual double | GetLocalValue (const VectorType &rVector, int LocalRow) const |
Short class definition. Detail class definition.
typedef TSpace::MatrixPointerType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::MatrixPointerType |
typedef TSpace::MatrixType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::MatrixType |
typedef TSpace::VectorPointerType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::VectorPointerType |
typedef TSpace::VectorType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::VectorType |
Type Definitions
|
inline |
Constructor. Empty constructor
|
delete |
Copy constructor.
|
virtualdefault |
Destructor.
|
inlinevirtual |
Auxiliar call to SetLocalValue for array type This method serves as an auxiliar call to the SetLocalValue function for array variables.
rInterfaceResidual | Interface residual vector in where the value is set |
rResidualValue | Double residual value |
AuxPosition | Position in rInterfaceResidual where the value is to be set. Note that the final position is computed using the TDim template value, since each entry in the residual vector contains a component of the residual variable value. |
|
inlinevirtual |
Auxiliar call to SetLocalValue for double type This method serves as an auxiliar call to the SetLocalValue function for double variables.
rInterfaceResidual | Interface residual vector in where the value is set |
rResidualValue | Double residual value |
AuxPosition | Position in rInterfaceResidual where the value is to be set |
|
inline |
|
inline |
|
inlinevirtual |
Checks if X = X0 + deltaX
rModelPart | reference to the model part in where the check has to be performed |
|
inlinevirtual |
Checks if X = X0 + deltaX
rModelPart | reference to the model part in where the check has to be performed |
|
inlinevirtual |
Computes and prints the fluid interface residual norms for debugging purposes
rInterfaceModelPart | interface modelpart in where the residual is computed |
|
inlinevirtual |
Computes and prints the structure interface residual norms for debugging purposes
rInterfaceModelPart | interface modelpart in where the residual is computed |
|
inlineprotected |
This function computes the nodal error of a vector magnitude in a consistent manner. The error is defined as the integral over the interface of a tests function times the difference between rOriginalVariable and rModifiedVariable.
rInterfaceModelPart | interface modelpart in where the residual is computed |
rOriginalVariable | variable with the reference value |
rModifiedVariable | variable with the computed vvalue |
rErrorStorageVariable | variable to store the error nodal value |
|
inline |
Computes the interface residual norm This method computes the interface residual norm. To do that it firstly computes the residual vector as is done in ComputeInterfaceResidualVector.
rInterfaceModelPart | interface modelpart in where the residual is computed |
rOriginalVariable | origin variable to compute the residual |
rModifiedVariable | end variable to compute the residual |
rInterfaceResidual | reference to the residual vector |
ResidualType | residual computation type (nodal or consistent) |
|
inlinevirtual |
Compute array variable interface residual vector This function computes (and stores in a vector) the residual of a vector variable over the fluid interface. The residual is defined as the OriginalVariable value minus the ModifiedVariable value. The nodal values of the residual are stored in the FSI_INTERFACE_RESIDUAL variable. Besides, the norm of the residual vector is stored in the ProcessInfo.
rInterfaceModelPart | interface modelpart in where the residual is computed |
rOriginalVariable | origin variable to compute the residual |
rModifiedVariable | end variable to compute the residual |
rInterfaceResidual | reference to the residual vector |
ResidualType | residual computation type (nodal or consistent) |
rArrayResidualVariable | variable to save the residual values |
rResidualNormVariable | variable to save the residual norm |
|
inlineprotected |
This function computes the nodal error of a vector magnitude. The error is defined as OriginalVariable minus ModifiedVariable.
rInterfaceModelPart | interface modelpart in where the residual is computed |
rOriginalVariable | variable with the reference value |
rModifiedVariable | variable with the computed vvalue |
rErrorStorageVariable | variable to store the error nodal value |
|
inline |
Create a coupling condition-based skin object This method creates an condition-based skin model part by copying the conditions of a given skin model part.
rOriginInterfaceModelPart | Origin skin model part to copy the conditions from |
rDestinationInterfaceModelPart | Empty destination modelpart to create the skin nodes and conditions |
|
inline |
Pressure to positive face pressure interpolator This method interpolates the pressure in the embedded fluid background mesh to the positive face pressure of all nodes in the provided structure interface model part.
rFluidModelPart | embedded fluid background mesh model part |
rStructureSkinModelPart | structure interface model part |
|
inline |
This function returns the interface length in 2D or the interface area in 3D.
rInterfaceModelPart | interface modelpart in where the are is computed |
|
inline |
This function computes the interface residual size as the number of interface nodes times the problem domain size.
|
inlineprotectedvirtual |
Reimplemented in Kratos::TrilinosPartitionedFSIUtilities< TSpace, TValueType, TDim >.
|
inlineprotected |
Get the skin condition name Auxiliary method that returns the auxiliary embedded skin condition type name.
|
inlineprotected |
Get the skin element name Auxiliary method that returns the auxiliary embedded skin element type name.
|
inline |
This function fills the interface vector (length equal to the residual size).
rInterfaceModelPart | interface modelpart from where the interface vector is filled |
rOriginVariable | variable to retrieve the values from |
rInterfaceVector | reference to the itnerface vector to be filled |
Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::KRATOS_CLASS_POINTER_DEFINITION | ( | PartitionedFSIUtilities< TSpace, TValueType, TDim > | ) |
|
inlineprotectedvirtual |
Reimplemented in Kratos::TrilinosPartitionedFSIUtilities< TSpace, TValueType, TDim >.
|
inlinevirtual |
This function resizes and sets to zero an interface vector (length equal to the residual size).
rInterfaceModelPart | interface modelpart in where the residual is computed |
Reimplemented in Kratos::TrilinosPartitionedFSIUtilities< TSpace, TValueType, TDim >.
|
inline |
Corrected guess local array value update This auxiliar method helps to update the nodal database using the values from a corrected guess vector values. Note that it performs the array components loop.
rCorrectedGuess | Corrected coupling guess vector |
rValueToUpdate | Reference in where the updated value is to be stored |
AuxPosition | Position in rCorrectedGuess from where the value is retrieved |
|
inline |
Corrected guess local double value update This auxiliar method helps to update the nodal database using the values from a corrected guess vector values.
rCorrectedGuess | Corrected coupling guess vector |
rValueToUpdate | Reference in where the updated value is to be stored |
AuxPosition | Position in rCorrectedGuess from where the value is retrieved |
|
inlinevirtual |
Sets the values in the corrected guess vector inside the selected nodal array variable.
rInterfaceModelPart | interface modelpart in where the residual is computed |
rSolutionVariable | variable in where the corrected solution is to be stored |
rCorrectedGuess | vector containing the interface corrected values |