![]() |
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.
|
Superconvergent patch recovery for linear meshes using quadratic polynomials. More...
#include <derivatives_recovery_utility.h>
Public Types | |
Type Definitions | |
typedef Node | NodeType |
Public Member Functions | |
Pointer Definitions | |
KRATOS_CLASS_POINTER_DEFINITION (DerivativesRecoveryUtility) | |
Static Public Member Functions | |
Operations | |
static void | Check (ModelPart &rModelPart) |
Check that all the required variables are present in the nodal database. More... | |
static void | CheckRequiredNeighborsPatch (ModelPart &rModelPart) |
Check that all the nodes have the required number of neighbors (6 in 2D and 10 in 3D) More... | |
static void | ExtendNeighborsPatch (ModelPart &rModelPart, const std::size_t RequiredNeighbors) |
Same as CheckRequiredNeighborsPatch with a custom number of required neighbors. More... | |
static void | CalculatePolynomialWeights (ModelPart &rModelPart) |
Fit a polynomial of degree=2 and store its derivative coefficients in FIRST_DERIVATIVE_WEIGHTS and SECOND_DERIVATIVE_WEIGHTS. More... | |
static void | RecoverDivergence (ModelPart &rModelPart, const Variable< array_1d< double, 3 >> &rOriginVariable, const Variable< double > &rDestinationVariable, const std::size_t BufferStep=0) |
Recover the divergence of a vector variable. More... | |
static void | RecoverGradient (ModelPart &rModelPart, const Variable< double > &rOriginVariable, const Variable< array_1d< double, 3 >> &rDestinationVariable, const std::size_t BufferStep=0) |
Recover the gradient of a scalar variable. More... | |
static void | RecoverLaplacian (ModelPart &rModelPart, const Variable< double > &rOriginVariable, const Variable< double > &rDestinationVariable, const std::size_t BufferStep=0) |
Recover the laplacian of a scalar variable. More... | |
static void | RecoverLaplacian (ModelPart &rModelPart, const Variable< array_1d< double, 3 >> &rOriginVariable, const Variable< array_1d< double, 3 >> &rDestinationVariable, const std::size_t BufferStep=0) |
Recover the laplacian of a vector variable as the gradient of the divergence. More... | |
Superconvergent patch recovery for linear meshes using quadratic polynomials.
Zhimin Zhangand and Ahmed Naga (2006) "A new finite element gradient recovery method: superconvergence property" SIAM J. Sci. Comput., 26(4), 1192–1213. (22 pages)
The polynomial follows the next convention: p_2(x, y, z, Node) = trans(P)*a trans(P) = (1, x, y, z, x^2, y^2, z^2, xy, xz, yz) trans(a) = a_i
typedef Node Kratos::DerivativesRecoveryUtility< TDim >::NodeType |
|
static |
Fit a polynomial of degree=2 and store its derivative coefficients in FIRST_DERIVATIVE_WEIGHTS and SECOND_DERIVATIVE_WEIGHTS.
It calls CheckRequiredNeighborsPatch
rModelPart |
|
static |
Check that all the required variables are present in the nodal database.
rModelPart |
|
static |
Check that all the nodes have the required number of neighbors (6 in 2D and 10 in 3D)
rModelPart |
|
static |
Same as CheckRequiredNeighborsPatch with a custom number of required neighbors.
rModelPart | |
RequiredNeighbors |
Kratos::DerivativesRecoveryUtility< TDim >::KRATOS_CLASS_POINTER_DEFINITION | ( | DerivativesRecoveryUtility< TDim > | ) |
|
static |
Recover the divergence of a vector variable.
rModelPart | |
rOriginVariable | |
rDestinationVariable | |
BufferStep |
|
static |
Recover the gradient of a scalar variable.
rModelPart | |
rOriginVariable | |
rDestinationVariable | |
BufferStep |
|
static |
Recover the laplacian of a vector variable as the gradient of the divergence.
rModelPart | |
rOriginVariable | |
rDestinationVariable | |
BufferStep |
|
static |
Recover the laplacian of a scalar variable.
rModelPart | |
rOriginVariable | |
rDestinationVariable | |
BufferStep |