Boundary Conditions¶
These functions add boundary terms to the residual evaluation.
-
enum RatelBoundaryType¶
Specify boundary condition.
Values:
-
enumerator RATEL_BOUNDARY_CLAMP¶
Dirichlet clamped boundary.
-
enumerator RATEL_BOUNDARY_MMS¶
Dirichlet MMS boundary.
-
enumerator RATEL_BOUNDARY_TRACTION¶
Neumann traction boundary,.
-
enumerator RATEL_BOUNDARY_CONTACT¶
Contact boundary condition.
-
enumerator RATEL_BOUNDARY_CLAMP¶
-
enum RatelContactShapeType¶
Specify rigid contact shape.
Values:
-
enumerator RATEL_CONTACT_SHAPE_PLATEN¶
Platen contact shape.
-
enumerator RATEL_CONTACT_SHAPE_PLATEN¶
-
enum RatelFrictionType¶
Specify friction model.
Values:
-
enumerator RATEL_FRICTION_NONE¶
No friction.
-
enumerator RATEL_FRICTION_COULOMB¶
Coulomb friction model.
-
enumerator RATEL_FRICTION_THRELFALL¶
Threlfall friction model.
-
enumerator RATEL_FRICTION_NONE¶
-
enum RatelBCInterpolationType¶
Specify interpolation type for flexible boundary condition.
Values:
-
enumerator RATEL_BC_INTERP_NONE¶
No interpolation, piecewise constant.
-
enumerator RATEL_BC_INTERP_LINEAR¶
Linear interpolation.
-
enumerator RATEL_BC_INTERP_NONE¶
-
enum RatelContactEnforcementType¶
Specify contact enforcement type.
Values:
-
enumerator RATEL_CONTACT_ENFORCEMENT_NITSCHE¶
-
enumerator RATEL_CONTACT_ENFORCEMENT_PENALTY¶
-
enumerator RATEL_CONTACT_ENFORCEMENT_NITSCHE¶
-
typedef PetscErrorCode RatelBCFunction(PetscInt, PetscReal, const PetscReal*, PetscInt, PetscScalar*, void*)¶
Function signature for BC function to set on DMPlex face.
-
typedef int (*RatelBCContactProject)(const RatelBCContactParams *context, const void *model, const CeedScalar x[3], CeedScalar y[3], CeedScalar n_y[3])¶
Function pointer type for computing the projection of a point onto the contact surface.
-
typedef int (*RatelBCContactUpdateModel)(const RatelBCContactParams *context, CeedScalar time, void *updated_model)¶
Function pointer type for updating the contact model based on the current time.
-
typedef int (*RatelBCContactFriction)(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar force_T[3], const CeedScalar force_n, CeedScalar traction[3], CeedScalar stored_values[3])¶
-
typedef int (*RatelBCContactFriction_fwd)(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar dvelocity_T[3], const CeedScalar force_T[3], const CeedScalar dforce_T[3], const CeedScalar force_n, const CeedScalar dforce_n, CeedScalar dtraction[3])¶
-
int MMSBCs_CEED_ScalarBPs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for CEED scalar BPs.
- Parameters:
ctx – [in] QFunction context, holding
RatelMMSCEEDBPsParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSBCs_CEED_VectorBPs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for CEED vector BPs.
- Parameters:
ctx – [in] QFunction context, holding
RatelMMSCEEDBPsParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int ClampBCs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet (clamp) boundary values.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCClampParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelContactTranslate(const RatelBCContactParamsCommon *contact, CeedScalar time, const CeedScalar point[3], CeedScalar translated_point[3])¶
Translate a point in the contact surface.
- Parameters:
contact – [in]
RatelBCContactParamsCommon
contexttime – [in] Current solver time
point – [in] Point to be translated
translated_point – [out] Translated point
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelBCContactGap(const RatelBCContactParams *context, const CeedScalar x[3], const CeedScalar y[3], const CeedScalar n_y[3], CeedScalar *gap_n)¶
Compute the gap between the contact surfaces.
- Parameters:
context – [in]
RatelBCContactParams
contextx – [in] Current position
y – [in] Nearest point on contact surface
n_y – [in] Normal vector to contact surface at y
gap_n – [out] Normal gap
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelBCContactGap_fwd(const RatelBCContactParams *context, const CeedScalar x[3], const CeedScalar y[3], const CeedScalar n_y[3], const CeedScalar dx[3], CeedScalar *gap_n, CeedScalar *dgap_n)¶
Compute the gap and its derivative with respect to x at the contact surface.
- Parameters:
context – [in]
RatelBCContactParams
context, which holds the contact parametersx – [in] Position in the current configuration
y – [in] Nearest point on contact surface at the current time step
n_y – [in] Normal to the contact surface at y
dx – [in] Variation in x, typically just du
gap_n – [out] Computed normal gap at the current time step
dgap_n – [out] Computed variation of the normal gap with respect to x
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelBCContactSlipVelocity(const RatelBCContactParams *context, const void *model, const void *model_prev, RatelBCContactProject fn_project, const CeedScalar v_x[3], const CeedScalar n_y[3], const CeedScalar y[3], const CeedScalar gap_n, CeedScalar v[3])¶
Compute the slip velocity at the contact surface.
- Parameters:
context – [in]
RatelBCContactParams
context, which holds the contact parametersmodel – [in] Model specific data for the current timestep
model_prev – [in] Model specific data for the previous timestep
fn_project – [in] Model projection function to compute the nearest point on the contact surface and normal vector
v_x – [in] Velocity at the mesh point
n_y – [in] Normal vector to the contact surface at the nearest point
y – [in] Nearest point on the contact surface at the current time step
gap_n – [in] Computed normal gap at the current time step
v – [out] Compute slip velocity at the contact surface
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionCoulomb(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar force_T[3], const CeedScalar force_n, CeedScalar traction[3], CeedScalar stored_values[3])¶
Compute friction traction using Coulomb model.
The Coulomb friction model is given by
\[ F = \begin{cases} F_C \ \sgn(\bm v_t), & \text{if } \|\bm v_t\| > 0, \\ \min(F_C, \|\bm \sigma_t\|) \sgn(\bm \sigma_t), & \text{if } \|\bm v_t\| = 0, \end{cases} \]RatelVecSignum()
for the definition of \( \sgn(\cdot) \).The Coulomb model is replaced by the following non-smooth operator for the Nitsche method:
\[ F = -\min(\max(0, -\mu_k (\sigma_n + \gamma g(\bm u)), \|\bm \sigma_t - \gamma \bm v_t\|) \sgn(\bm \sigma_t - \gamma \bm v_t), \]Note
If
force_n
is non-negative, settraction
to zero. Otherwise, clamp magnitude oftraction
to be less than or equal to-f * force_n
Note
Uses
RatelFrictionParams::kinetic_coefficient
as \( \mu_k \).- Parameters:
friction – [in] Friction parameters
velocity_T – [in] Tangential velocity
force_T – [in] Tangential force
force_n – [in] Normal force, possibly positive if Nitsche method is used
traction – [out] Computed traction force
stored_values – [out] Stored values for Jacobian evaluation
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionCoulomb_fwd(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar dvelocity_T[3], const CeedScalar force_T[3], const CeedScalar dforce_T[3], const CeedScalar force_n, const CeedScalar dforce_n, CeedScalar dtraction[3])¶
Compute the linearization of the friction traction using Coulomb model.
Note
Uses
RatelFrictionParams::kinetic_coefficient
as \( \mu_k \).- Parameters:
friction – [in] Friction parameters
velocity_T – [in] Tangential velocity on contact surface (unused)
dvelocity_T – [in] Linearization of tangential velocity on contact surface
force_T – [in] Tangential force on contact surface
dforce_T – [in] Linearization of tangential force on contact surface
force_n – [in] Normal force on contact surface
dforce_n – [in] Linearization of normal force on contact surface
dtraction – [out] Linearization of the computed friction traction
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionThrelfall(const RatelFrictionParams *friction, const CeedScalar tangent_velocity[3], const CeedScalar force_T[3], const CeedScalar force_n, CeedScalar traction[3], CeedScalar stored_values[3])¶
Compute friction traction using Threlfall model.
The Threlfall model is given by
\[ F = \begin{cases} F_C \left(\frac{1 - \exp\left(\frac{-3 ||\bm v_t||}{v_0}\right)}{1 - \exp(-3)}\right) \sgn(\bm v_t), & \text{if } ||\bm v_t|| \leq v_0, \\ F_C \ \sgn(\bm v_t), & \text{otherwise}, \end{cases} \]RatelVecSignum()
for the definition of \( \sgn(\cdot) \).Note
Uses
RatelFrictionParams::kinetic_coefficient
as \( \mu_k \) andRatelFrictionParams::tolerance_velocity
as \( v_0 \).- Parameters:
friction – [in] Friction parameters
tangent_velocity – [in] Velocity in the tangent plane to the contact surface
force_T – [in] Tangential force (unused)
force_n – [in] Normal force, possibly positive if Nitsche method is used
traction – [out] Computed friction traction force
stored_values – [out] Stored values for Jacobian evaluation
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionThrelfall_fwd(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar dvelocity_T[3], const CeedScalar force_T[3], const CeedScalar dforce_T[3], const CeedScalar force_n, const CeedScalar dforce_n, CeedScalar dtraction[3])¶
Compute the derivative of the friction traction using Threlfall model.
Note
Uses
RatelFrictionParams::kinetic_coefficient
as \( \mu_k \) andRatelFrictionParams::tolerance_velocity
as \( v_0 \).- Parameters:
friction – [in] Friction parameters
velocity_T – [in] Tangential velocity
dvelocity_T – [in] Linearization of tangential velocity
force_T – [in] Tangential force on contact surface
dforce_T – [in] Linearization of tangential force on contact surface
force_n – [in] Normal force
dforce_n – [in] Linearization of normal force
dtraction – [out] Linearization of the computed friction traction
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionNone(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar force_T[3], const CeedScalar force_n, CeedScalar traction[3], CeedScalar stored_values[3])¶
No-op for frictionless case.
- Parameters:
friction – [in] Friction parameters (unused)
velocity_T – [in] Tangential velocity (unused)
force_T – [in] Tangential force (unused)
force_n – [in] Normal force, possibly positive if Nitsche method is used (unused)
traction – [out] Computed traction force (set to zero)
stored_values – [out] Stored values for Jacobian evaluation (set to zero)
- Returns:
An error code: 0 - success, otherwise - failure
-
int FrictionNone_fwd(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar dvelocity_T[3], const CeedScalar force_T[3], const CeedScalar dforce_T[3], const CeedScalar force_n, const CeedScalar dforce_n, CeedScalar dtraction[3])¶
No-op for frictionless case.
- Parameters:
friction – [in] Friction parameters (unused)
velocity_T – [in] Tangential velocity (unused)
dvelocity_T – [in] Linearization of tangential velocity (unused)
force_T – [in] Tangential force on contact surface (unused)
dforce_T – [in] Linearization of tangential force on contact surface (unused)
force_n – [in] Normal force (unused)
dforce_n – [in] Linearization of normal force (unused)
dtraction – [out] Linearization of the computed friction traction (set to zero)
- Returns:
An error code: 0 - success, otherwise - failure
-
int ContactBCs_Nitsche(void *ctx, CeedInt Q, RatelComputef1 compute_f1, bool has_state_values, bool has_stored_values, CeedInt num_active_field_eval_modes, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the surface integral of the Nitsche contact boundary condition.
Updated to support initial gap using https://hal.archives-ouvertes.fr/hal-00722114/document.
In the frictionless case, compute
where(f1_gamma < 0) ? (f1_gamma * n_contact * (w detNb)) : 0
f1_gamma = f1_n + gamma*gap
f1_n = (f1·n_contact)·n_contact
is the contact pressure between the target surface and contact surfacegap = (X + u - C)·n_contact
is the gap between the target surface and the contact surface in the current configurationn_contact
is the outward facing normal to the analytic contact surface
For linear elasticity
f1 = sigma
, wheresigma
is Cauchy stress tensor. For hyperelastic in initial configurationf1 = P
, whereP
is First Piola Kirchhoff stress. For hyperelastic in current configurationf1 = tau
, where tau is Kirchhoff stress.For the frictional case, see the theory section of the documentation
- Parameters:
ctx – [in] QFunction context, holding
RatelBCContactParams
Q – [in] Number of quadrature points
compute_f1 – [in] Function to compute f1
has_state_values – [in] Boolean flag indicating model state values in residual evaluation
has_stored_values – [in] Boolean flag indicating model stores values in residual evaluation
num_active_field_eval_modes – [in] Number of active field evaluation modes
in – [in] Input arrays
0 - qdata
input_data_offset
-u
input_data_offset + 1
-du/dt
(oru
in static case)input_data_offset + 2
- quadrature point coordinates
out – [out] Output array
output_data_offset
- action of QFunction onu
output_data_offset +
- storedf1_n + gamma*gap
andf1_gamma_T - gamma*(du/dt)_T
- Returns:
An error code: 0 - success, otherwise - failure
-
int ContactBCs_Nitsche_Jacobian(void *ctx, CeedInt Q, RatelComputef1_fwd compute_df1, bool has_state_values, bool has_stored_values, CeedInt num_active_field_eval_modes, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the Jacobian of the surface integral for Nitsche method based contact on the constrained faces.
In the frictionless case compute
where(f1_gamma < 0) ? ((df1_N + gamma * du_N) * n_contact * (w detNb)) : 0
df1_N = (df1·n_contact)·n_contact
is the incremental surface force away from contact surface.du_N = du·n_contact
is the incremental displacement away from contact surface.
For linear elasticity
df1 = dsigma
. For hyperelastic in initial configurationdf1 = dP
. For hyperelastic in current configurationdf1 = dtau - tau * grad_du^T
.For the frictional case, see the theory section of the documentation
- Parameters:
ctx – [in] QFunction context, holding
RatelBCContactParams
Q – [in] Number of quadrature points
compute_df1 – [in] Function to compute f1
has_state_values – [in] Boolean flag indicating model state values in residual evaluation
has_stored_values – [in] Boolean flag indicating model stores values in residual evaluation
num_active_field_eval_modes – [in] Number of active fields
in – [in] Input arrays
0 - qdata
input_data_offset
-du
, incremental change to uinput_data_offset + 1
- stored values from residual
out – [out] Output array
0 - action of QFunction
- Returns:
An error code: 0 - success, otherwise - failure
-
int SetupSurfaceGeometrySplit(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute geometric factors for integration, gradient transformations, and coordinate transformations on element faces.
Reference (parent) 2D coordinates are given by
X
and physical (current) 3D coordinates are given byx
. The change of coordinate matrix is given bydxdX_{i,j} = dx_i/dX_j (indicial notation) [3 * 2]
.(N_1, N_2, N_3)
is given by the cross product of the columns ofdxdX_{i,j}
.detNb
is the magnitude of(N_1, N_2, N_3)
.- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - Jacobian of face coordinates
1 - quadrature weights
out – [out] Output array
0 - qdata,
w
,detNb
, andN
- Returns:
An error code: 0 - success, otherwise - failure
-
int ContactBCs_Penalty(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the surface integral for penalty method based contact on the constrained faces.
In the frictionless case, compute
where(gap < 0) ? (w detNb) * gap / epsilon * N : 0
epsilon = sqrt(detNb) / context->gamma
is the penalty parameter- Parameters:
ctx – [in] QFunction context, holding
RatelBCContactParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - quadrature data
1 -
u
2 -
du/dt
(oru
in static case)3 - quadrature point coordinates
out – [out] Output array
0 -
v
1 - stored values from residual
- Returns:
An error code: 0 - success, otherwise - failure
-
int ContactBCs_Penalty_Jacobian(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the Jacobian of the surface integral for penalty method based contact on the constrained faces.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCContactParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - quadrature data
1 - incremental change in
u
2 - stored values from residual
out – [out] Output array
0 - action of QFunction on
du
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelBCContactProject_Platen(const RatelBCContactParams *context, const void *model, const CeedScalar x[3], CeedScalar y[3], CeedScalar n_y[3])¶
Compute the nearest point on the contact surface.
- Parameters:
context – [in]
RatelBCContactParams
contextmodel – [in] Model specific data for the current timestep,
RatelBCContactPlatenParams
x – [in] Current position
y – [out] Nearest point on contact surface
n_y – [out] Normal vector to contact surface at y
- Returns:
An error code: 0 - success, otherwise - failure
-
int RatelBCContactUpdateModel_Platen(const RatelBCContactParams *context, CeedScalar time, void *updated_model)¶
Update the model for the given time.
- Parameters:
context – [in]
RatelBCContactParams
contexttime – [in] Current solver time
updated_model – [out] Model with updated parameters
- Returns:
An error code: 0 - success, otherwise - failure
-
int SetupDirichletBCs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Setup coordinate and scaling data for Dirichlet boundary.
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - nodal coordinates
1 - nodal multiplicity
out – [out] Output array
0 - stored coordinates
1 - stored multiplicity scaling factor
- Returns:
An error code: 0 - success, otherwise - failure
-
CeedInt RatelBCInterpGetPreviousKnotIndex(CeedScalar t, CeedInt num_knots, const CeedScalar *knots)¶
Compute the index
j \in [0, num_knots]
such thatknots[j] <= t < knots[j+1]
.If
t < knots[0]
, returns -1.- Parameters:
t – [in] Current time
num_knots – [in] Number of knots
knots – [in] Values of knots
- Returns:
A scalar value: the index of the previous knot
-
CeedScalar RatelBCInterpScaleTime(CeedScalar t, CeedInt prev_knot_index, CeedInt num_knots, const CeedScalar *knots)¶
Compute rescaled time between knots, \(s \in [0,1)\).
If
0 <= prev_knot_index < num_knots-1
, returnss = (t - knots[prev_knot_index]) / (knots[prev_knot_index + 1])
. Ift
is before the first knot, i.e. prev_knot_index == -1, implicitly inserts a knot at time 0 and returnss = t / knots[0]
. Ift
is past the final knot, returnss = 1
. Ift > knots[num_knots - 1]
, returnspoints[num_knots - 1]
.- Parameters:
t – [in] Current time
prev_knot_index – [in] Index of last knot
num_knots – [in] Number of knots
knots – [in] Values of knots
- Returns:
A scalar value: the rescaled time between knots
-
int RatelBCInterpolate(RatelBCInterpolationType type, CeedScalar t, CeedInt prev_knot_index, CeedInt num_knots, const CeedScalar *knots, CeedInt point_dim, const CeedScalar *points, CeedScalar *out_point)¶
Compute the interpolated value between the given points and knots.
For
t > knots[num_knots-1]
, returnspoints[num_knots-1]
.- Parameters:
type – [in] If type is
RATEL_BC_INTERP_NONE
- returns the point corresponding to the first knot after time t.RATEL_BC_INTERP_LINEAR
- returns the result of linearly interpolating between adjacent points.
t – [in] Current time
prev_knot_index – [in] Index of last knot
num_knots – [in] Number of knots
knots – [in] Values of knots
point_dim – [in] Dimension of coordinates
points – [in] Coordinates of knots
out_point – [out] Interpolated value
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSBCs_Linear(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for linear elasticity MMS.
Uses an MMS adapted from doi:10.1016/j.compstruc.2019.106175
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSBCs_MixedLinear(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for linear elasticity MMS.
Uses an MMS adapted from doi:10.1016/j.compstruc.2019.106175
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSBCs_PoroElasticityLinear_u(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for linear poroelasticity MMS.
Uses an MMS adapted from doi:10.1016/j.compstruc.2019.106175
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values (size of 3)
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSBCs_PoroElasticityLinear_p(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet MMS boundary values for linear poroelasticity MMS.
Uses an MMS adapted from doi:10.1016/j.compstruc.2019.106175
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored coordinates
1 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values (size of 1)
- Returns:
An error code: 0 - success, otherwise - failure
-
int PressureBCs(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the surface integral for pressure on the constrained faces.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCPressureParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - gradient of coordinates
1 - quadrature weights
2 - gradient of u with respect to reference coordinates
out – [out] Output array
0 -
\int v^T . (pn) ds
- Returns:
An error code: 0 - success, otherwise - failure
-
int PressureBCsJacobian(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute the Jacobian of the surface integral for pressure on the constrained faces.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCPressureParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - gradient of coordinates
1 - quadrature weights
2 - gradient of incremental change in
u
out – [out] Output array
0 - action of QFunction on
dug
- Returns:
An error code: 0 - success, otherwise - failure
-
int SetupSlipBCs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Setup scaling data for Dirichlet (clamp) boundary.
- Parameters:
ctx – [in] QFunction context, unused
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - nodal multiplicity
out – [out] Output array
0 - stored multiplicity scaling factor
- Returns:
An error code: 0 - success, otherwise - failure
-
int SlipBCs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute Dirichlet (slip) boundary values.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCSlipParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - stored multiplicity scaling factor
out – [out] Output array
0 - Dirichlet values
- Returns:
An error code: 0 - success, otherwise - failure
-
int TractionBCs(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute surface integral of the traction vector on the constrained faces.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCTractionParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - surface qdata
out – [out] Output array
0 -
v = t * (w detNb)
- Returns:
An error code: 0 - success, otherwise - failure
-
int TractionEnergy(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute external energy cause by the traction vector on the constrained faces.
- Parameters:
ctx – [in] QFunction context, holding
RatelBCTractionParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - surface qdata
1 - displacement u
out – [out] Output array
0 - traction energy
u t * (w detNb)
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSError_MixedLinear(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute error from true solution for mixed linear elasticity MMS.
- Parameters:
ctx – [in] QFunction context, holding
RatelLinearElasticityParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - qdata
1 - quadrature point coordinates
2 - computed solution, displacement field
3 - computed solution, pressure field
out – [out] Output array
0 - error from true solution, displacement field
0 - error from true solution, pressure field
- Returns:
An error code: 0 - success, otherwise - failure
-
int MMSError_PoroElasticityLinear(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out)¶
Compute error from true solution for mixed linear poroelasticity MMS.
- Parameters:
ctx – [in] QFunction context, holding
RatelMMSPoroElasticityLinearParams
Q – [in] Number of quadrature points
in – [in] Input arrays
0 - qdata
1 - quadrature point coordinates
2 - computed solution, displacement field
3 - computed solution, pressure field
out – [out] Output array
0 - error from true solution, displacement field
1 - error from true solution, pressure field
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelCreateBCLabel(DM dm, const char name[])¶
Create boundary label.
Not collective across MPI processes.
- Parameters:
dm – [inout]
DM
to add boundary labelname – [in] Name for new boundary label
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelDMAddBoundariesDirichlet(Ratel ratel, DM dm)¶
Add Dirichlet boundary conditions to
DM
.Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [out]
DM
to update with Dirichlet boundaries
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelDMAddBoundariesSlip(Ratel ratel, DM dm)¶
Add Dirichlet slip boundaries to
DM
.Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [out]
DM
to update with dirichlet slip boundaries
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundarySlipDataFromOptions(Ratel ratel, PetscInt i, PetscInt j, RatelBCSlipParams *params_slip)¶
Read slip boundary conditions from options database.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contexti – [in] Field index
j – [in] Boundary index
params_slip – [out] Data structure to store slip data
-
PetscErrorCode RatelCeedAddBoundariesDirichletSlip(Ratel ratel, DM dm, PetscBool incremental, CeedOperator op_dirichlet)¶
Add Dirichlet slip boundary conditions to Dirichlet boundary condition operator.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for Dirichlet boundariesincremental – [in] Set up operator to incrementally apply slip boundary conditions
op_dirichlet – [out] Composite Dirichlet
CeedOperator
to add sub-operators to
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelCeedAddBoundariesDirichletMMS(Ratel ratel, DM dm, CeedOperator op_dirichlet)¶
Add Dirichlet MMS boundary conditions to Dirichlet boundary condition operator.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for MMS boundariesop_dirichlet – [out] Composite Dirichlet
CeedOperator
to add sub-operators to
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryClampDataFromOptions(Ratel ratel, PetscInt i, PetscInt j, RatelBCClampParams *params_clamp)¶
Read Dirichlet boundary conditions from options database.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contexti – [in] Field index
j – [in] Boundary index
params_clamp – [out] Data structure to store clamp data
-
PetscErrorCode RatelCeedAddBoundariesDirichletClamp(Ratel ratel, DM dm, PetscBool incremental, CeedOperator op_dirichlet)¶
Add Dirichlet clamp boundary conditions to Dirichlet boundary condition operator.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for Dirichlet boundariesincremental – [in] Set up operator to incrementally apply slip boundary conditions
op_dirichlet – [out] Composite Dirichlet
CeedOperator
to add sub-operators to
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryTractionDataFromOptions(Ratel ratel, PetscInt i, RatelBCTractionParams *params_traction)¶
Read traction boundary conditions from options database.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contexti – [in] Boundary index
params_traction – [out] Data structure to store traction data
-
PetscErrorCode RatelCeedSetupSurfaceQData(Ratel ratel, DM dm, const char *label_name, PetscInt label_value, CeedElemRestriction *restriction, CeedVector *q_data)¶
Compute
CeedOperator
surface QData.Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for Neumann boundarieslabel_name – [in]
DMPlex
label name for surfacelabel_value – [in]
DMPlex
label value for surfacerestriction – [out]
CeedElemRestriction
for QDataq_data – [out]
CeedVector
holding QData
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelCeedAddBoundariesNeumann(Ratel ratel, DM dm, CeedOperator op_residual)¶
Add Neumann boundary conditions to residual operator.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for Neumann boundariesop_residual – [out] Composite residual u term
CeedOperator
to add Neumann sub-operators
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelSetupTractionEnergySuboperator(Ratel ratel, DM dm_energy, CeedOperator op_external_energy)¶
Add traction energy to energy operator.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm_energy – [in]
DM
for strain energy computationop_external_energy – [out] Composite external energy
CeedOperator
to add traction energy sub-operators
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelSetupSurfaceForceCentroids(Ratel ratel, DM dm)¶
Setup surface force face centroid computation.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
with surface force faces
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundarySetupSurfaceDisplacementSuboperators(Ratel ratel, DM dm_surface_displacement, CeedOperator ops_surface_displacement[])¶
Setup average surface displacement
CeedOperators
- Parameters:
ratel – [in]
Ratel
contextdm_surface_displacement – [in]
DM
for surface displacement computationops_surface_displacement – [out] Composite surface displacement
CeedOperator
objects to addRatelMaterial
suboperator to
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelMaterialContactContextRegisterFields(RatelMaterial material, CeedQFunctionContext ctx)¶
Register material specific fields for contact boundary conditions
CeedQFunctionContext
Not collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
to setup contact contextctx – [inout]
CeedQFunctionContext
for contact
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelContactCommonRegisterFields(Ratel ratel, const RatelBCContactParamsCommon *params_contact, CeedQFunctionContext ctx)¶
Register common fields for contact boundary conditions
CeedQFunctionContext
Not collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams_contact – [in]
RatelBCContactParamsCommon
for contactctx – [out] Context to register fields in
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelBoundaryContactCeedContextCreate(RatelMaterial material, RatelBCContactParams *params_contact, CeedQFunctionContext *ctx)¶
Create libCEED
CeedQFunctionContext
for contact boundary conditions.Not collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
to setup contact contextparams_contact – [in]
RatelBCContactParams
contextctx – [out]
CeedQFunctionContext
for contact
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelFrictionParamsView(Ratel ratel, const RatelFrictionParams *params_friction, PetscViewer viewer)¶
Print friction parameters to a viewer.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams_friction – [in]
RatelFrictionParams
to printviewer – [inout]
PetscViewer
to print to
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelFrictionParamsCheck_Coulomb(Ratel ratel, RatelFrictionParams *params_friction)¶
Check Coulomb friction parameters.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams_friction – [in] Friction parameters to check
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelFrictionParamsCheck_Threlfall(Ratel ratel, RatelFrictionParams *params_friction)¶
Check Threlfall friction parameters.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams_friction – [in] Friction parameters to check
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelFrictionParamsFromOptions(Ratel ratel, const char option_prefix[], RatelFrictionParams *params_friction)¶
Set friction parameters from command line options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoption_prefix – [in] Prefix string for command line options
params_friction – [out]
RatelFrictionParams
to set
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryContactPlatenParamsView(Ratel ratel, const RatelBCContactPlatenParams *params_platen, PetscViewer viewer)¶
Print platen boundary condition parameters to a viewer.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams_platen – [in]
RatelBCContactParams
to printviewer – [inout]
PetscViewer
to print to
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryContactParamsView(Ratel ratel, const RatelBCContactParams *params, PetscViewer viewer)¶
Print contact boundary condition parameters to a viewer.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextparams – [in]
RatelBCContactParams
to printviewer – [inout]
PetscViewer
to print to
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryPlatenParamsFromOptions(Ratel ratel, const char options_prefix[], RatelBCContactPlatenParams *params_platen)¶
Read platen boundary condition parameters from options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoptions_prefix – [in] Prefix string for command line options
params_platen – [out]
RatelBCContactParams
to set from options, must be allocated
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryContactParamsCommonFromOptions(Ratel ratel, const char options_prefix[], RatelContactShapeType shape, RatelBCContactParamsCommon *params_contact)¶
Read contact boundary condition parameters from options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoptions_prefix – [in] Prefix string for command line options
shape – [in] Model shape type
params_contact – [out]
RatelBCContactParamsCommon
to set from options, must be allocated
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryContactParamsFromOptions(Ratel ratel, const char options_prefix[], RatelBCContactParams *params_contact)¶
Read contact boundary condition parameters from options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoptions_prefix – [in] Prefix string for command line options
params_contact – [out]
RatelBCContactParams
to set from options, must be allocated
- Returns:
An error code: 0 - success, otherwise - failure
-
static PetscErrorCode RatelMaterialContactContextCreate(RatelMaterial material, const char *contact_name, PetscInt face_id, PetscInt face_orientation, PetscInt name_index, CeedQFunctionContext *ctx, CeedSize *jacobian_flops_estimate)¶
Process command line options for contact boundary conditions.
Collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
to setup contact contextcontact_name – [in]
const char *
name of contact faceface_id – [in]
DMPlex
face domain idface_orientation – [in]
DMPlex
face domain stratum valuename_index – [in] Index of contact face in name array
ctx – [out]
CeedQFunctionContext
for contact BCjacobian_flops_estimate – [out] Estimate of Jacobian FLOPS for the contact model, or NULL if not needed.
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelMaterialCeedAddBoundariesContactNitsche(RatelMaterial material, CeedVector u_dot_loc, CeedOperator op_residual, CeedOperator op_jacobian)¶
Add contact boundary conditions to the residual u term and Jacobian operators.
Collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
contextu_dot_loc – [in]
CeedVector
for passive input velocity fieldop_residual – [out] Composite residual u term
CeedOperator
to add RatelMaterial sub-operatorop_jacobian – [out] Composite Jacobian
CeedOperator
to addRatelMaterial
sub-operator
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelMaterialSetupContactNitscheJacobianMultigridLevel(RatelMaterial material, DM dm_level, CeedVector m_loc, CeedOperator sub_op_jacobian_fine, CeedOperator op_jacobian_coarse, CeedOperator op_prolong, CeedOperator op_restrict)¶
Setup multigrid level suboperator for
RatelMaterial
contact boundary conditions Jacobian.Collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
contextdm_level – [in]
DMPlex
for multigrid level to setupm_loc – [in]
CeedVector
holding multiplicitysub_op_jacobian_fine – [in] Fine grid Jacobian
CeedOperator
op_jacobian_coarse – [inout] Composite Jacobian
CeedOperator
to addRatelMaterial
suboperatorsop_prolong – [inout] Composite prolongation
CeedOperator
to addRatelMaterial
suboperators, unused but needed for compatable function signatureop_restrict – [inout] Composite restriction
CeedOperator
to addRatelMaterial
suboperators, unused but needed for compatable function signature
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelCeedAddBoundariesContactPenalty(Ratel ratel, DM dm, CeedVector u_dot_loc, CeedOperator op_residual, CeedOperator op_jacobian)¶
Add contact penalty boundary conditions to the residual u term and Jacobian operators.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for pressure boundariesu_dot_loc – [in] Local CeedVector for time derivative of u
op_residual – [out] Composite residual u term
CeedOperator
to add contact sub-operatorsop_jacobian – [out] Composite Jacobian
CeedOperator
to add contact sub-operators
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelMaterialSetupContactPenaltyJacobianMultigridLevel(RatelMaterial material, DM dm_level, CeedVector m_loc, CeedOperator sub_op_jacobian_fine, CeedOperator op_jacobian_coarse, CeedOperator op_prolong, CeedOperator op_restrict)¶
Setup multigrid level suboperator for
RatelMaterial
contact penalty boundary conditions Jacobian.Collective across MPI processes.
- Parameters:
material – [in]
RatelMaterial
contextdm_level – [in]
DMPlex
for multigrid level to setupm_loc – [in]
CeedVector
holding multiplicitysub_op_jacobian_fine – [in] Fine grid Jacobian
CeedOperator
op_jacobian_coarse – [inout] Composite Jacobian
CeedOperator
to addRatelMaterial
suboperatorsop_prolong – [inout] Composite prolongation
CeedOperator
to addRatelMaterial
suboperators, unused but needed for compatable function signatureop_restrict – [inout] Composite restriction
CeedOperator
to addRatelMaterial
suboperators, unused but needed for compatable function signature
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundaryPressureDataFromOptions(Ratel ratel, PetscInt i, RatelBCPressureParams *params_pressure)¶
Read pressure boundary conditions from options database.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contexti – [in] Boundary index
params_pressure – [out] Data structure to store pressure data
-
PetscErrorCode RatelCeedAddBoundariesPressure(Ratel ratel, DM dm, CeedOperator op_residual, CeedOperator op_jacobian)¶
Add pressure boundary conditions to residual and Jacobian operators.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm – [in]
DM
to use for pressure boundariesop_residual – [out] Composite residual u term
CeedOperator
to add pressure sub-operatorsop_jacobian – [out] Composite Jacobian
CeedOperator
to add pressure sub-operators
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelSetupPressureJacobianMultigridLevel(Ratel ratel, DM dm_level, CeedVector m_loc, CeedOperator op_jacobian_fine, CeedOperator op_jacobian_coarse)¶
Setup multigrid level suboperators for
RatelMaterial
pressure boundary conditions Jacobian.Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextdm_level – [in]
DMPlex
for multigrid level to setupm_loc – [in]
CeedVector
holding multiplicityop_jacobian_fine – [in] Fine grid Jacobian
CeedOperator
op_jacobian_coarse – [inout] Composite Jacobian
CeedOperator
to addRatelMaterial
suboperators
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelBoundingBoxParamsFromOptions(Ratel ratel, const char option_prefix[], const char name[], RatelBoundingBoxParams *bounding_box)¶
Read bounding box from options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoption_prefix – [in] Prefix for options
name – [in] Name of face
bounding_box – [out] Output bounding box
- Returns:
An error code: 0 - success, otherwise - failure
-
PetscErrorCode RatelFaceLabelValueFromOptions(Ratel ratel, const char option_prefix[], const char name[], PetscInt *label_value)¶
Read label value from options.
Collective across MPI processes.
- Parameters:
ratel – [in]
Ratel
contextoption_prefix – [in] Prefix for options
name – [in] Name of face
label_value – [out] Label value read from options
- Returns:
An error code: 0 - success, otherwise - failure
-
RATEL_MAX_BC_INTERP_POINTS 256¶
Maximum interpolation nodes for boundary conditions.
-
RATEL_MAX_CONTACT_MODEL_SIZE 128¶
Maximum size of contact model.
-
FLOPS_ContactGap (3 + FLOPS_Dot3)¶
-
FLOPS_ContactGap_fwd (FLOPS_ContactGap + FLOPS_Dot3)¶
-
FLOPS_FrictionCoulomb (1 + FLOPS_VecSignum + 3 * (2))¶
-
FLOPS_FrictionCoulomb_fwd (2 + 3 * (5) + FLOPS_VecSignum_fwd)¶
-
FLOPS_FrictionThrelfall (1 + FLOPS_VecSignum + 7 + 2 + 3 * (7))¶
-
FLOPS_FrictionThrelfall_fwd (2 + FLOPS_VecSignum_fwd + FLOPS_Dot3 + 7 + 9 + 3 * (3 + 4))¶
-
NUM_COMPONENTS_STORED_Contact_Nitsche 7¶
-
FLOPS_Contact_Nitsche_Jacobian_without_df1 (FLOPS_ContactGap_fwd + 3 * (1 + FLOPS_Dot3) + FLOPS_Dot3 + 3 * 2 + 3 + 3 + 3 + FLOPS_Dot3 + 3 * 4 + 3 * 3)¶
-
NUM_COMPONENTS_STORED_Contact_Penalty 7¶
-
FLOPS_Contact_Penalty_Jacobian (3 + FLOPS_Dot3 + 2 + 3 * 2 + 3 * 4 + 3 * 5)¶
-
FLOPS_ContactProject_Platen (3 + FLOPS_Dot3 + 3 * 2)¶
-
FLOPS_Jacobian_PressureBC 30¶
-
struct RatelBCClampParams¶
- #include <boundary.h>
Clamp boundary condition context.
Public Members
-
CeedScalar translation[3 * RATEL_MAX_BC_INTERP_POINTS]¶
Translation vectors.
-
CeedScalar rotation_axis[3 * RATEL_MAX_BC_INTERP_POINTS]¶
Rotation axes.
-
CeedScalar rotation_polynomial[2 * RATEL_MAX_BC_INTERP_POINTS]¶
Rotation polynomial coefficients.
-
CeedScalar times[RATEL_MAX_BC_INTERP_POINTS]¶
Transition times.
-
CeedInt num_times¶
Number of transition times.
-
RatelBCInterpolationType interpolation_type¶
Type of interpolation between points.
-
CeedInt num_comp¶
Number of components in field.
-
CeedScalar time¶
Current solver time.
-
CeedScalar dt¶
Current solver time step.
-
CeedScalar translation[3 * RATEL_MAX_BC_INTERP_POINTS]¶
-
struct RatelBCSlipParams¶
- #include <boundary.h>
Slip boundary condition context.
Public Members
-
CeedScalar translation[4 * RATEL_MAX_BC_INTERP_POINTS]¶
Translation vectors.
-
CeedScalar times[RATEL_MAX_BC_INTERP_POINTS]¶
Transition times.
-
CeedInt num_times¶
Number of transition times.
-
RatelBCInterpolationType interpolation_type¶
Type of interpolation between points.
-
CeedInt components[4]¶
Components to constrain.
-
CeedInt num_comp_slip¶
Number of constrained components.
-
CeedInt num_comp¶
Number of components in field.
-
CeedScalar time¶
Current solver time.
-
CeedScalar dt¶
Current solver time step.
-
CeedScalar translation[4 * RATEL_MAX_BC_INTERP_POINTS]¶
-
struct RatelBCTractionParams¶
- #include <boundary.h>
Traction boundary condition context.
Public Members
-
CeedScalar direction[3 * RATEL_MAX_BC_INTERP_POINTS]¶
Traction vector.
-
CeedScalar times[RATEL_MAX_BC_INTERP_POINTS]¶
Transition times.
-
CeedInt num_times¶
Number of transition times.
-
RatelBCInterpolationType interpolation_type¶
Type of interpolation between points.
-
CeedInt num_comp¶
Number of components in field.
-
CeedScalar time¶
Current solver time.
-
CeedScalar dt¶
Current solver time step.
-
CeedScalar direction[3 * RATEL_MAX_BC_INTERP_POINTS]¶
-
struct RatelBCPressureParams¶
- #include <boundary.h>
Pressure boundary condition context.
Public Members
-
CeedScalar pressure[RATEL_MAX_BC_INTERP_POINTS]¶
Pressure.
-
CeedScalar times[RATEL_MAX_BC_INTERP_POINTS]¶
Transition times.
-
CeedInt num_times¶
Number of transition times.
-
RatelBCInterpolationType interpolation_type¶
Type of interpolation between points.
-
CeedInt num_comp¶
Number of components in field.
-
CeedScalar time¶
Current solver time.
-
CeedScalar dt¶
Current solver time step.
-
CeedScalar pressure[RATEL_MAX_BC_INTERP_POINTS]¶
-
struct RatelBoundingBoxParams¶
- #include <bounding-box.h>
-
struct RatelBCContactParamsCommon¶
- #include <contact-common.h>
Common Contact BC data.
Public Members
-
CeedScalar translation[3 * RATEL_MAX_BC_INTERP_POINTS]¶
Translation vector.
-
CeedScalar times[RATEL_MAX_BC_INTERP_POINTS]¶
Transition times.
-
CeedInt num_times¶
Number of transition times.
-
RatelBCInterpolationType interpolation_type¶
Type of interpolation between points.
-
CeedScalar gamma¶
Nitsche’s method parameter.
-
RatelFrictionParams friction¶
Friction parameters.
-
CeedInt face_id¶
DM face id.
-
CeedInt face_domain_value¶
DM face orientation.
-
CeedInt name_index¶
Index for face name lookup.
-
CeedScalar translation[3 * RATEL_MAX_BC_INTERP_POINTS]¶
-
struct RatelBCContactParams¶
- #include <contact-common.h>
Contact BC data, including material data.
Public Members
-
RatelBCContactParamsCommon contact¶
Common contact BC data.
-
CeedScalar model[RATEL_MAX_CONTACT_MODEL_SIZE]¶
Specific contact model data.
-
RatelContactShapeType shape¶
Model shape.
-
CeedScalar time¶
Current solver time.
-
CeedScalar dt¶
Current time step.
-
CeedScalar material[RATEL_MAX_MATERIAL_SIZE]¶
Material data.
-
RatelBCContactParamsCommon contact¶
-
struct RatelFrictionParams¶
- #include <friction.h>
Friction parameters for all friction models, some of which may be unused in a specific model.
Public Members
-
RatelFrictionType model¶
Friction model.
-
CeedScalar static_coefficient¶
Static coefficient of friction.
-
CeedScalar kinetic_coefficient¶
Kinetic coefficient of friction, assumed to equal the static coefficient for some models.
-
CeedScalar viscous_coefficient¶
Viscous damping coefficient, Ns/m.
-
CeedScalar tolerance_velocity¶
Tolerance velocity, m/s.
-
RatelFrictionType model¶
-
struct RatelBCContactPlatenParams¶
- #include <platen.h>
Common Platen BC data.