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_PLATEN

Platen boundary integral.

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.

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.

enum RatelPlatenType

Specify platen type.

Values:

enumerator RATEL_PLATEN_NITSCHE
enumerator RATEL_PLATEN_PENALTY
typedef PetscErrorCode RatelBCFunction(PetscInt, PetscReal, const PetscReal*, PetscInt, PetscScalar*, void*)

Function signature for BC function to set on DMPlex face.

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 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 that knots[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, returns s = (t - knots[prev_knot_index]) / (knots[prev_knot_index + 1]). If t is before the first knot, i.e. prev_knot_index == -1, implicitly inserts a knot at time 0 and returns s = t / knots[0]. If t is past the final knot, returns s = 1. If t > knots[num_knots - 1], returns points[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], returns points[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

void FrictionCoulomb(const RatelFrictionParams *friction, const CeedScalar force_T[3], const CeedScalar force_n, CeedScalar traction[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} \]
where \( F_C = -\mu_k \sigma_n \) is the Coulomb friction force, \( \mu_k \) is the kinetic friction coefficient, \( \bm v_t \) is the tangential velocity, and \( \bm \sigma_t \) is the tangential force. See 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), \]
where \( \gamma \) is the penalty parameter, \( g(\bm u) \) is the gap function, and \( \bm u \) is the displacement.

Note

If force_n is non-negative, set traction to zero. Otherwise, clamp magnitude of traction to be less than or equal to -f * force_n

Note

Uses RatelFrictionParams::kinetic_coefficient as \( \mu_k \).

Parameters:
  • friction[in] Friction parameters

  • force_T[in] Tangential force

  • force_n[in] Normal force, possibly positive if Nitsche method is used

  • traction[out] Computed traction force

Returns:

An error code: 0 - success, otherwise - failure

void FrictionCoulomb_fwd(const RatelFrictionParams *friction, 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

  • 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

void FrictionThrelfall(const RatelFrictionParams *friction, const CeedScalar tangent_velocity[3], const CeedScalar force_n, CeedScalar traction[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} \]
where \( F_C = -\mu_k (\hat\sigma_n) \) is the Coulomb friction force, \( \bm v_t \) is the tangential velocity, \( v_0 \) is the tolerance velocity, and \( \mu_k \) is the kinetic friction coefficient. See RatelVecSignum() for the definition of \( \sgn(\cdot) \).

Parameters:
  • friction[in] Friction parameters

  • tangent_velocity[in] Velocity in the tangent plane to the contact surface

  • force_n[in] Normal force, possibly positive if Nitsche method is used

  • traction[out] Computed friction traction force

Returns:

An error code: 0 - success, otherwise - failure

void FrictionThrelfall_fwd(const RatelFrictionParams *friction, const CeedScalar velocity_T[3], const CeedScalar dvelocity_T[3], const CeedScalar force_n, const CeedScalar dforce_n, CeedScalar dtraction[3])

Compute the derivative of the friction traction using Threlfall model.

Parameters:
  • friction[in] Friction parameters

  • velocity_T[in] Tangential velocity

  • dvelocity_T[in] Linearization of tangential velocity

  • 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 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 PlatenBCs(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 platen contact boundary condition.

Updated to support initial gap using https://hal.archives-ouvertes.fr/hal-00722114/document.

In the frictionless case, compute

(f1_gamma < 0) ? (f1_gamma * n_platen * (w detNb)) : 0
where
  • f1_gamma = f1_n + gamma*gap

  • f1_n = (f1·n_platen)·n_platen is the contact pressure between the surface and platen

  • gap = (X + u - C)·n_platen is the gap between the surface and the platen in the current configuration

For linear elasticity f1 = sigma, where sigma is Cauchy stress tensor. For hyperelastic in initial configuration f1 = P, where P is First Piola Kirchhoff stress. For hyperelastic in current configuration f1 = tau, where tau is Kirchhoff stress.

For the frictional case, see the theory section of the documentation

Parameters:
  • ctx[in] QFunction platen, holding RatelBCPlatenParams

  • 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 (or u in static case)

    • input_data_offset + 2 - quadrature point coordinates

  • out[out] Output array

    • output_data_offset - action of QFunction on u

    • output_data_offset + - stored f1_n + gamma*gap and f1_gamma_T - gamma*(du/dt)_T

Returns:

An error code: 0 - success, otherwise - failure

int PlatenBCs_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 incremental traction due to platen contact boundary condition.

In the frictionless case compute

(f1_gamma < 0) ? ((df1_N + gamma * du_N) * n_platen * (w detNb)) : 0
where
  • df1_N = (df1·n_platen)·n_platen is the incremental surface force away from platen.

  • du_N = du·n_platen is the incremental displacement away from platen.

For linear elasticity df1 = dsigma. For hyperelastic in initial configuration df1 = dP. For hyperelastic in current configuration df1 = dtau - tau * grad_du^T.

For the frictional case, see the theory section of the documentation

Parameters:
  • ctx[in] QFunction platen, holding RatelBCPlatenParams

  • 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 u

    • input_data_offset + 1 - stored values from residual

  • out[out] Output array

    • 0 - action of QFunction

Returns:

An error code: 0 - success, otherwise - failure

int SetupPenaltyPlatens(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 by x. 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 of dxdX_{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, and N

Returns:

An error code: 0 - success, otherwise - failure

int PlatenPenaltyBCs(void *ctx, const 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

(gap < 0) ? (w detNb) * gap / epsilon * N : 0
where epsilon = sqrt(detNb) / platen->gamma is the penalty parameter

Parameters:
  • ctx[in] QFunction context, holding RatelBCPlatenParams

  • Q[in] Number of quadrature points

  • in[in] Input arrays

    • 0 - quadrature data

    • 1 - u

    • 2 - du/dt (or u 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 PlatenPenaltyBCsJacobian(void *ctx, const 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 RatelBCPlatenParams

  • 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 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 label

  • name[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 context

  • dm[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 context

  • dm[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 context

  • i[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 context

  • dm[in] DM to use for Dirichlet boundaries

  • incremental[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 context

  • dm[in] DM to use for MMS boundaries

  • op_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 context

  • i[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 context

  • dm[in] DM to use for Dirichlet boundaries

  • incremental[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 context

  • i[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 context

  • dm[in] DM to use for Neumann boundaries

  • label_name[in] DMPlex label name for surface

  • label_value[in] DMPlex label value for surface

  • restriction[out] CeedElemRestriction for QData

  • q_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 context

  • dm[in] DM to use for Neumann boundaries

  • op_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 context

  • dm_energy[in] DM for strain energy computation

  • op_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 context

  • dm[in] DM with surface force faces

Returns:

An error code: 0 - success, otherwise - failure

static PetscErrorCode RatelMaterialPlatenContextRegisterFields(RatelMaterial material, CeedQFunctionContext ctx)

Register material specific fields for platen boundary conditions CeedQFunctionContext

Not collective across MPI processes.

Parameters:
  • material[in] RatelMaterial to setup platen context

  • ctx[inout] CeedQFunctionContext for platen

Returns:

An error code: 0 - success, otherwise - failure

static PetscErrorCode RatelCeedPlatenContextCreate(RatelMaterial material, const RatelBCPlatenParamsCommon *params_platen, CeedQFunctionContext *ctx)

Create libCEED CeedQFunctionContext for platen boundary conditions.

Not collective across MPI processes.

Parameters:
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 context

  • params_friction[in] RatelFrictionParams to print

  • viewer[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 context

  • params_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 context

  • params_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 context

  • option_prefix[in] Prefix string for command line options

  • params_friction[out] RatelFrictionParams to set

Returns:

An error code: 0 - success, otherwise - failure

PetscErrorCode RatelBoundaryPlatenParamsCommonView(Ratel ratel, const RatelBCPlatenParamsCommon *params_platen, PetscViewer viewer)

Print platen boundary condition parameters to a viewer.

Collective across MPI processes.

Parameters:
Returns:

An error code: 0 - success, otherwise - failure

PetscErrorCode RatelBoundaryPlatenParamsCommonFromOptions(Ratel ratel, const char options_prefix[], RatelBCPlatenParamsCommon *params_platen)

Read platen boundary condition parameters from options.

Collective across MPI processes.

Parameters:
  • ratel[in] Ratel context

  • options_prefix[in] Prefix string for command line options

  • params_platen[out] RatelBCPlatenParamsCommon to set from options, must be allocated

Returns:

An error code: 0 - success, otherwise - failure

static PetscErrorCode RatelMaterialPlatenContextCreate(RatelMaterial material, const char *platen_name, PetscInt face_id, PetscInt face_orientation, PetscInt name_index, CeedQFunctionContext *ctx)

Process command line options for platen boundary conditions.

Collective across MPI processes.

Parameters:
  • material[in] RatelMaterial to setup platen context

  • platen_name[in] const char * name of platen

  • face_id[in] DMPlex face domain id

  • face_orientation[in] DMPlex face domain stratum value

  • name_index[in] Index of platen in name array

  • ctx[out] CeedQFunctionContext for platen

Returns:

An error code: 0 - success, otherwise - failure

PetscErrorCode RatelMaterialCeedAddBoundariesPlatenNitsche(RatelMaterial material, CeedOperator op_residual, CeedVector u_dot_loc, CeedOperator op_jacobian)

Add platen contact boundary conditions to the residual u term operator.

Collective across MPI processes.

Parameters:
  • material[in] RatelMaterial context

  • op_residual[out] Composite residual u term CeedOperator to add RatelMaterial sub-operator

  • u_dot_loc[out] CeedVector for passive input velocity field

  • op_jacobian[out] Composite Jacobian CeedOperator to add RatelMaterial sub-operator

Returns:

An error code: 0 - success, otherwise - failure

PetscErrorCode RatelMaterialSetupPlatenNitscheJacobianMultigridLevel(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 platen boundary conditions Jacobian.

Collective across MPI processes.

Parameters:
  • material[in] RatelMaterial context

  • dm_level[in] DMPlex for multigrid level to setup

  • m_loc[in] CeedVector holding multiplicity

  • sub_op_jacobian_fine[in] Fine grid Jacobian CeedOperator

  • op_jacobian_coarse[inout] Composite Jacobian CeedOperator to add RatelMaterial suboperators

  • op_prolong[inout] Composite prolongation CeedOperator to add RatelMaterial suboperators, unused but needed for compatable function signature

  • op_restrict[inout] Composite restriction CeedOperator to add RatelMaterial suboperators, unused but needed for compatable function signature

Returns:

An error code: 0 - success, otherwise - failure

PetscErrorCode RatelCeedAddBoundariesPlatenPenalty(Ratel ratel, DM dm, CeedOperator op_residual, CeedVector u_dot_loc, CeedOperator op_jacobian)
PetscErrorCode RatelMaterialSetupPlatenPenaltyJacobianMultigridLevel(RatelMaterial material, DM dm_level, CeedVector m_loc, CeedOperator sub_op_jacobian_fine, CeedOperator op_jacobian_coarse, CeedOperator op_prolong, CeedOperator op_restrict)
PetscErrorCode RatelBoundaryPressureDataFromOptions(Ratel ratel, PetscInt i, RatelBCPressureParams *params_pressure)

Read pressure boundary conditions from options database.

Collective across MPI processes.

Parameters:
  • ratel[in] Ratel context

  • i[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 context

  • dm[in] DM to use for pressure boundaries

  • op_residual[out] Composite residual u term CeedOperator to add pressure sub-operators

  • op_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 context

  • dm_level[in] DMPlex for multigrid level to setup

  • m_loc[in] CeedVector holding multiplicity

  • op_jacobian_fine[in] Fine grid Jacobian CeedOperator

  • op_jacobian_coarse[inout] Composite Jacobian CeedOperator to add RatelMaterial 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 context

  • option_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 context

  • option_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.

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))
FLOPS_Platen_without_df1 (7 * FLOPS_Dot3 + 64)
FLOPS_Jacobian_PlatenPenaltyBC (3 * FLOPS_Dot3 + FLOPS_Norm3 + 57)
FLOPS_Jacobian_PressureBC 30
struct RatelBCClampParams
#include <boundary.h>

Clamp boundary condition context.

Public Members

CeedScalar translation[4 * RATEL_MAX_BC_INTERP_POINTS]

Translation vectors.

CeedScalar rotation_axis[4 * 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_u

Number of components in field.

CeedScalar time

Current solver time.

CeedScalar dt

Current solver time step.

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_u

Number of components in field.

CeedScalar time

Current solver time.

CeedScalar dt

Current solver time step.

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.

CeedScalar time

Current solver time.

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.

CeedScalar time

Current solver time.

struct RatelBoundingBoxParams
#include <bounding-box.h>
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.

struct RatelBCPlatenParamsCommon
#include <platen.h>

Common Platen BC data.

Public Members

CeedScalar normal[3]

Exterior normal for platen.

CeedScalar center[3]

Center of the platen.

CeedScalar distance[RATEL_MAX_BC_INTERP_POINTS]

Displacement of the platen along normal 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 platen name lookup.

struct RatelBCPlatenParams
#include <platen.h>

Platen BC data, including material data.

Public Members

RatelBCPlatenParamsCommon platen

Common platen data.

CeedScalar time

Current solver time.

CeedScalar material[RATEL_MAX_MATERIAL_SIZE]

Material data.