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.
Public Types | Public Member Functions | List of all members
Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim > Class Template Reference

#include <partitioned_fsi_utilities.hpp>

Inheritance diagram for Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >:
Collaboration diagram for Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >:

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
 

Detailed Description

template<class TSpace, class TValueType, unsigned int TDim>
class Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >

Short class definition. Detail class definition.

Member Typedef Documentation

◆ MatrixPointerType

template<class TSpace , class TValueType , unsigned int TDim>
typedef TSpace::MatrixPointerType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::MatrixPointerType

◆ MatrixType

template<class TSpace , class TValueType , unsigned int TDim>
typedef TSpace::MatrixType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::MatrixType

◆ VectorPointerType

template<class TSpace , class TValueType , unsigned int TDim>
typedef TSpace::VectorPointerType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::VectorPointerType

◆ VectorType

template<class TSpace , class TValueType , unsigned int TDim>
typedef TSpace::VectorType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::VectorType

Type Definitions

Constructor & Destructor Documentation

◆ PartitionedFSIUtilities() [1/2]

template<class TSpace , class TValueType , unsigned int TDim>
Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::PartitionedFSIUtilities ( )
inline

Constructor. Empty constructor

◆ PartitionedFSIUtilities() [2/2]

template<class TSpace , class TValueType , unsigned int TDim>
Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::PartitionedFSIUtilities ( const PartitionedFSIUtilities< TSpace, TValueType, TDim > &  Other)
delete

Copy constructor.

◆ ~PartitionedFSIUtilities()

template<class TSpace , class TValueType , unsigned int TDim>
virtual Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::~PartitionedFSIUtilities ( )
virtualdefault

Destructor.

Member Function Documentation

◆ AuxSetLocalValue() [1/2]

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::AuxSetLocalValue ( VectorType rInterfaceResidual,
const array_1d< double, 3 > &  rResidualValue,
const int  AuxPosition 
)
inlinevirtual

Auxiliar call to SetLocalValue for array type This method serves as an auxiliar call to the SetLocalValue function for array variables.

Parameters
rInterfaceResidualInterface residual vector in where the value is set
rResidualValueDouble residual value
AuxPositionPosition 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.

◆ AuxSetLocalValue() [2/2]

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::AuxSetLocalValue ( VectorType rInterfaceResidual,
const double rResidualValue,
const int  AuxPosition 
)
inlinevirtual

Auxiliar call to SetLocalValue for double type This method serves as an auxiliar call to the SetLocalValue function for double variables.

Parameters
rInterfaceResidualInterface residual vector in where the value is set
rResidualValueDouble residual value
AuxPositionPosition in rInterfaceResidual where the value is to be set

◆ CalculateTractionFromPressureValues() [1/2]

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::CalculateTractionFromPressureValues ( ModelPart rModelPart,
const Variable< double > &  rPositivePressureVariable,
const Variable< double > &  rNegativePressureVariable,
const Variable< array_1d< double, 3 >> &  rTractionVariable,
const bool  SwapTractionSign 
)
inline

◆ CalculateTractionFromPressureValues() [2/2]

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::CalculateTractionFromPressureValues ( ModelPart rModelPart,
const Variable< double > &  rPressureVariable,
const Variable< array_1d< double, 3 >> &  rTractionVariable,
const bool  SwapTractionSign 
)
inline

◆ CheckCurrentCoordinatesFluid()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::CheckCurrentCoordinatesFluid ( ModelPart rModelPart,
const double  tolerance 
)
inlinevirtual

Checks if X = X0 + deltaX

Parameters
rModelPartreference to the model part in where the check has to be performed

◆ CheckCurrentCoordinatesStructure()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::CheckCurrentCoordinatesStructure ( ModelPart rModelPart,
const double  tolerance 
)
inlinevirtual

Checks if X = X0 + deltaX

Parameters
rModelPartreference to the model part in where the check has to be performed

◆ ComputeAndPrintFluidInterfaceNorms()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::ComputeAndPrintFluidInterfaceNorms ( ModelPart rInterfaceModelPart)
inlinevirtual

Computes and prints the fluid interface residual norms for debugging purposes

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed

◆ ComputeAndPrintStructureInterfaceNorms()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::ComputeAndPrintStructureInterfaceNorms ( ModelPart rInterfaceModelPart)
inlinevirtual

Computes and prints the structure interface residual norms for debugging purposes

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed

◆ ComputeConsistentResidual()

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::ComputeConsistentResidual ( ModelPart rInterfaceModelPart,
const Variable< TValueType > &  rOriginalVariable,
const Variable< TValueType > &  rModifiedVariable,
const Variable< TValueType > &  rErrorStorageVariable 
)
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.

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
rOriginalVariablevariable with the reference value
rModifiedVariablevariable with the computed vvalue
rErrorStorageVariablevariable to store the error nodal value

◆ ComputeInterfaceResidualNorm()

template<class TSpace , class TValueType , unsigned int TDim>
double Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::ComputeInterfaceResidualNorm ( ModelPart rInterfaceModelPart,
const Variable< TValueType > &  rOriginalVariable,
const Variable< TValueType > &  rModifiedVariable,
const Variable< TValueType > &  rResidualVariable,
const std::string  ResidualType = "nodal" 
)
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.

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
rOriginalVariableorigin variable to compute the residual
rModifiedVariableend variable to compute the residual
rInterfaceResidualreference to the residual vector
ResidualTyperesidual computation type (nodal or consistent)
Returns
double interface residual norm

◆ ComputeInterfaceResidualVector()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::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 
)
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.

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
rOriginalVariableorigin variable to compute the residual
rModifiedVariableend variable to compute the residual
rInterfaceResidualreference to the residual vector
ResidualTyperesidual computation type (nodal or consistent)
rArrayResidualVariablevariable to save the residual values
rResidualNormVariablevariable to save the residual norm

◆ ComputeNodeByNodeResidual()

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::ComputeNodeByNodeResidual ( ModelPart rInterfaceModelPart,
const Variable< TValueType > &  rOriginalVariable,
const Variable< TValueType > &  rModifiedVariable,
const Variable< TValueType > &  rErrorStorageVariable 
)
inlineprotected

This function computes the nodal error of a vector magnitude. The error is defined as OriginalVariable minus ModifiedVariable.

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
rOriginalVariablevariable with the reference value
rModifiedVariablevariable with the computed vvalue
rErrorStorageVariablevariable to store the error nodal value

◆ CreateCouplingSkin()

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::CreateCouplingSkin ( const ModelPart rOriginInterfaceModelPart,
ModelPart rDestinationInterfaceModelPart 
)
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.

Parameters
rOriginInterfaceModelPartOrigin skin model part to copy the conditions from
rDestinationInterfaceModelPartEmpty destination modelpart to create the skin nodes and conditions

◆ EmbeddedPressureToPositiveFacePressureInterpolator()

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::EmbeddedPressureToPositiveFacePressureInterpolator ( ModelPart rFluidModelPart,
ModelPart rStructureSkinModelPart 
)
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.

Parameters
rFluidModelPartembedded fluid background mesh model part
rStructureSkinModelPartstructure interface model part

◆ GetInterfaceArea()

template<class TSpace , class TValueType , unsigned int TDim>
double Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::GetInterfaceArea ( ModelPart rInterfaceModelPart)
inline

This function returns the interface length in 2D or the interface area in 3D.

Parameters
rInterfaceModelPartinterface modelpart in where the are is computed
Returns
the given modelpart interface length

◆ GetInterfaceResidualSize()

template<class TSpace , class TValueType , unsigned int TDim>
int Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::GetInterfaceResidualSize ( ModelPart rInterfaceModelPart)
inline

This function computes the interface residual size as the number of interface nodes times the problem domain size.

Returns
the model part residual size

◆ GetLocalValue()

template<class TSpace , class TValueType , unsigned int TDim>
virtual double Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::GetLocalValue ( const VectorType rVector,
int  LocalRow 
) const
inlineprotectedvirtual

◆ GetSkinConditionName()

template<class TSpace , class TValueType , unsigned int TDim>
std::string Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::GetSkinConditionName ( )
inlineprotected

Get the skin condition name Auxiliary method that returns the auxiliary embedded skin condition type name.

Returns
std::string Condition type registering name

◆ GetSkinElementName()

template<class TSpace , class TValueType , unsigned int TDim>
std::string Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::GetSkinElementName ( )
inlineprotected

Get the skin element name Auxiliary method that returns the auxiliary embedded skin element type name.

Returns
std::string Element type registering name

◆ InitializeInterfaceVector()

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::InitializeInterfaceVector ( const ModelPart rInterfaceModelPart,
const Variable< TValueType > &  rOriginVariable,
VectorType rInterfaceVector 
)
inline

This function fills the interface vector (length equal to the residual size).

Parameters
rInterfaceModelPartinterface modelpart from where the interface vector is filled
rOriginVariablevariable to retrieve the values from
rInterfaceVectorreference to the itnerface vector to be filled

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<class TSpace , class TValueType , unsigned int TDim>
Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::KRATOS_CLASS_POINTER_DEFINITION ( PartitionedFSIUtilities< TSpace, TValueType, TDim >  )

◆ SetLocalValue()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::SetLocalValue ( VectorType rVector,
int  LocalRow,
double  Value 
) const
inlineprotectedvirtual

◆ SetUpInterfaceVector()

template<class TSpace , class TValueType , unsigned int TDim>
virtual VectorPointerType Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::SetUpInterfaceVector ( ModelPart rInterfaceModelPart)
inlinevirtual

This function resizes and sets to zero an interface vector (length equal to the residual size).

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
Returns
pointer to the vector that has been created

Reimplemented in Kratos::TrilinosPartitionedFSIUtilities< TSpace, TValueType, TDim >.

◆ UpdateInterfaceLocalValue() [1/2]

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::UpdateInterfaceLocalValue ( const VectorType rCorrectedGuess,
array_1d< double, 3 > &  rValueToUpdate,
const int  AuxPosition 
)
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.

Parameters
rCorrectedGuessCorrected coupling guess vector
rValueToUpdateReference in where the updated value is to be stored
AuxPositionPosition in rCorrectedGuess from where the value is retrieved

◆ UpdateInterfaceLocalValue() [2/2]

template<class TSpace , class TValueType , unsigned int TDim>
void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::UpdateInterfaceLocalValue ( const VectorType rCorrectedGuess,
double rValueToUpdate,
const int  AuxPosition 
)
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.

Parameters
rCorrectedGuessCorrected coupling guess vector
rValueToUpdateReference in where the updated value is to be stored
AuxPositionPosition in rCorrectedGuess from where the value is retrieved

◆ UpdateInterfaceValues()

template<class TSpace , class TValueType , unsigned int TDim>
virtual void Kratos::PartitionedFSIUtilities< TSpace, TValueType, TDim >::UpdateInterfaceValues ( ModelPart rInterfaceModelPart,
const Variable< TValueType > &  rSolutionVariable,
const VectorType rCorrectedGuess 
)
inlinevirtual

Sets the values in the corrected guess vector inside the selected nodal array variable.

Parameters
rInterfaceModelPartinterface modelpart in where the residual is computed
rSolutionVariablevariable in where the corrected solution is to be stored
rCorrectedGuessvector containing the interface corrected values

The documentation for this class was generated from the following file: