Brighenti's liquid crystal elastomer model.
More...
#include <liquid_crystal_elastomer.hpp>
|
struct | State |
| internal variables for the liquid crystal elastomer model More...
|
|
|
| LiquidCrystElastomerBrighenti (double rho, double shear_modulus, double bulk_modulus, double order_constant, double order_parameter, double transition_temperature, double N_b_squared) |
| Constructor. More...
|
|
template<typename DispGradType , typename OrderParamType , typename GammaAngleType > |
SERAC_HOST_DEVICE auto | operator() (State &state, const tensor< DispGradType, dim, dim > &displacement_grad, OrderParamType temperature_tuple, GammaAngleType gamma_tuple) const |
| Material response. More...
|
|
template<typename S , typename T , typename U > |
auto | calculateDistributionTensor (const State &state, const tensor< S, dim, dim > &F_hat, const T theta, const tensor< U, dim > &normal) const |
| Compute the distribution tensor using Brighenti's model. More...
|
|
|
static constexpr int | dim = 3 |
| this model is only intended to be used in 3D
|
|
Brighenti's liquid crystal elastomer model.
Definition at line 29 of file liquid_crystal_elastomer.hpp.
◆ LiquidCrystElastomerBrighenti()
serac::LiquidCrystElastomerBrighenti::LiquidCrystElastomerBrighenti |
( |
double |
rho, |
|
|
double |
shear_modulus, |
|
|
double |
bulk_modulus, |
|
|
double |
order_constant, |
|
|
double |
order_parameter, |
|
|
double |
transition_temperature, |
|
|
double |
N_b_squared |
|
) |
| |
|
inline |
Constructor.
- Parameters
-
rho | mass density of the material (in reference configuration) |
shear_modulus | Shear modulus of the material |
bulk_modulus | Bulk modulus of the material |
order_constant | temperature-valued constant in exponential factor for order parameter |
order_parameter | Initial value of the order parameter |
transition_temperature | Characteristic temperature of the order-disorder transition |
N_b_squared | Number of Kunh segments/chain, times square of Kuhn segment length |
Definition at line 51 of file liquid_crystal_elastomer.hpp.
◆ calculateDistributionTensor()
template<typename S , typename T , typename U >
auto serac::LiquidCrystElastomerBrighenti::calculateDistributionTensor |
( |
const State & |
state, |
|
|
const tensor< S, dim, dim > & |
F_hat, |
|
|
const T |
theta, |
|
|
const tensor< U, dim > & |
normal |
|
) |
| const |
|
inline |
Compute the distribution tensor using Brighenti's model.
- Template Parameters
-
S | Type of the deformation gradient |
T | Type of angle of the liquid crystals |
U | Type of the unit length normal of liquid crystal alignment (first-order tensor) |
- Parameters
-
state | State variables for this material |
F_hat | Dot product of the current and previous deformation gradients |
theta | In-plane angle of the liquid crystals in the elastomer matrix |
normal | Unit length normal of angle alignment |
- Returns
- distribution tensor
Definition at line 146 of file liquid_crystal_elastomer.hpp.
◆ operator()()
template<typename DispGradType , typename OrderParamType , typename GammaAngleType >
SERAC_HOST_DEVICE auto serac::LiquidCrystElastomerBrighenti::operator() |
( |
State & |
state, |
|
|
const tensor< DispGradType, dim, dim > & |
displacement_grad, |
|
|
OrderParamType |
temperature_tuple, |
|
|
GammaAngleType |
gamma_tuple |
|
) |
| const |
|
inline |
Material response.
- Template Parameters
-
DispGradType | number-like type for the displacement gradient tensor |
OrderParamType | number-like type for the order parameter |
GammaAngleType | number-like type for the orientation angle gamma |
- Parameters
-
[in,out] | state | A state variable object for this material. The value is updated in place. |
[in] | displacement_grad | displacement gradient with respect to the reference configuration |
[in] | temperature_tuple | the temperature |
[in] | gamma_tuple | the first polar angle used to define the liquid crystal orientation vector |
- Returns
- The calculated material response (Cauchy stress) for the material
Definition at line 83 of file liquid_crystal_elastomer.hpp.
◆ density
double serac::LiquidCrystElastomerBrighenti::density |
The documentation for this struct was generated from the following file: