Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Classes | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
serac::LiquidCrystElastomerBrighenti Struct Reference

Brighenti's liquid crystal elastomer model. More...

#include <liquid_crystal_elastomer.hpp>

Classes

struct  State
 internal variables for the liquid crystal elastomer model More...
 

Public Member Functions

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

Public Attributes

double density
 mass density More...
 
double shear_modulus_
 shear modulus in stress-free configuration
 
double bulk_modulus_
 bulk modulus in stress-free configuration
 
double order_constant_
 Order constant.
 
double initial_order_parameter_
 initial value of order parameter
 
double transition_temperature_
 Transition temperature.
 
double N_b_squared_
 Kuhn segment parameters.
 

Static Public Attributes

static constexpr int dim = 3
 this model is only intended to be used in 3D
 

Detailed Description

Brighenti's liquid crystal elastomer model.

Definition at line 29 of file liquid_crystal_elastomer.hpp.

Constructor & Destructor Documentation

◆ 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
rhomass density of the material (in reference configuration)
shear_modulusShear modulus of the material
bulk_modulusBulk modulus of the material
order_constanttemperature-valued constant in exponential factor for order parameter
order_parameterInitial value of the order parameter
transition_temperatureCharacteristic temperature of the order-disorder transition
N_b_squaredNumber of Kunh segments/chain, times square of Kuhn segment length

Definition at line 51 of file liquid_crystal_elastomer.hpp.

Member Function Documentation

◆ 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
SType of the deformation gradient
TType of angle of the liquid crystals
UType of the unit length normal of liquid crystal alignment (first-order tensor)
Parameters
stateState variables for this material
F_hatDot product of the current and previous deformation gradients
thetaIn-plane angle of the liquid crystals in the elastomer matrix
normalUnit 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
DispGradTypenumber-like type for the displacement gradient tensor
OrderParamTypenumber-like type for the order parameter
GammaAngleTypenumber-like type for the orientation angle gamma
Parameters
[in,out]stateA state variable object for this material. The value is updated in place.
[in]displacement_graddisplacement gradient with respect to the reference configuration
[in]temperature_tuplethe temperature
[in]gamma_tuplethe 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.

Member Data Documentation

◆ density

double serac::LiquidCrystElastomerBrighenti::density

mass density


Definition at line 189 of file liquid_crystal_elastomer.hpp.


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