5 #ifndef COMPRESSIBLEHYPERELASTICITY_LINEARPATCHTEST_HPP 6 #define COMPRESSIBLEHYPERELASTICITY_LINEARPATCHTEST_HPP 23 std::string mesh_file = Teuchos::getParameter<std::string>(Parameters.sublist(
"Mesh"),
"mesh_file");
24 Mesh =
new mesh(comm, mesh_file, 1000.0);
27 StandardMap =
new Epetra_Map(-1,3*Mesh->n_local_nodes_without_ghosts,&Mesh->local_dof_without_ghosts[0],0,*
Comm);
28 OverlapMap =
new Epetra_Map(-1,3*Mesh->n_local_nodes,&Mesh->local_dof[0],0,*
Comm);
47 mu = 2.0*mu1 + 4.0*mu2 + 2.0*
mu3;
49 ptrmbeta4 = std::pow(trm,beta4);
50 ptrmbeta5 = std::pow(trm,beta5);
52 cos_plyagl = std::cos(plyagl);
53 sin_plyagl = std::sin(plyagl);
56 Epetra_SerialDenseVector
exactSolution(
double & x1,
double & x2,
double & x3){
57 Epetra_SerialDenseVector
u(3);
64 double errorL2(Epetra_Vector & uStandardMap){
75 Epetra_SerialDenseVector uExact(3), vH(3);
77 Epetra_SerialDenseMatrix u_G(3,n_gauss_points), x_G(3,n_gauss_points);
81 for (
unsigned int inode=0; inode<
Mesh->
el_type; ++inode){
92 for (
unsigned int gp=0; gp<n_gauss_points; ++gp){
95 vH(0) = uExact(0) - u_G(0,gp);
96 vH(1) = uExact(1) - u_G(1,gp);
97 vH(2) = uExact(2) - u_G(2,gp);
102 Comm->SumAll(&error,&totalError,1);
103 totalError = std::sqrt(totalError);
120 if(x==0||y==0||z==0||x==1.0||y==1.0||z==1.0){
132 if(x==0||y==0||z==0||x==1.0||y==1.0||z==1.0){
152 if (x==0||y==0||z==0||x==1.0||y==1.0||z==1.0){
153 Epetra_SerialDenseVector
u(3);
163 F.Update(-1.0,rhs_dir,1.0);
170 if (x==0||y==0||z==0||x==1.0||y==1.0||z==1.0){
182 void get_constitutive_tensors(Epetra_SerialDenseMatrix & deformation_gradient, Epetra_SerialDenseVector & piola_stress, Epetra_SerialDenseMatrix & tangent_piola){
184 double det = deformation_gradient(0,0)*deformation_gradient(1,1)*deformation_gradient(2,2)
185 - deformation_gradient(0,0)*deformation_gradient(1,2)*deformation_gradient(2,1)
186 - deformation_gradient(0,1)*deformation_gradient(1,0)*deformation_gradient(2,2)
187 + deformation_gradient(0,1)*deformation_gradient(1,2)*deformation_gradient(2,0)
188 + deformation_gradient(0,2)*deformation_gradient(1,0)*deformation_gradient(2,1)
189 - deformation_gradient(0,2)*deformation_gradient(1,1)*deformation_gradient(2,0);
191 Epetra_SerialDenseMatrix C(3,3);
192 Epetra_SerialDenseMatrix CC(3,3);
193 Epetra_SerialDenseMatrix ddJ5(6,6);
194 Epetra_SerialDenseVector LML(6);
195 Epetra_SerialDenseVector eye(6);
196 Epetra_SerialDenseVector dJ5(6);
197 Epetra_SerialDenseVector M(6);
198 Epetra_SerialDenseVector
L(6);
199 Epetra_SerialDenseVector c(6);
201 C.Multiply(
'T',
'N',1.0,deformation_gradient,deformation_gradient,0.0);
202 CC.Multiply(
'N',
'N',1.0,C,C,0.0);
204 c(0) = C(0,0); c(1) = C(1,1); c(2) = C(2,2); c(3) = C(1,2); c(4) = C(0,2); c(5) = C(0,1);
206 L(0) = (1.0/(det*det))*(C(1,1)*C(2,2)-C(1,2)*C(2,1));
207 L(1) = (1.0/(det*det))*(C(0,0)*C(2,2)-C(0,2)*C(2,0));
208 L(2) = (1.0/(det*det))*(C(0,0)*C(1,1)-C(0,1)*C(1,0));
209 L(3) = (1.0/(det*det))*(C(0,2)*C(1,0)-C(0,0)*C(1,2));
210 L(4) = (1.0/(det*det))*(C(0,1)*C(1,2)-C(0,2)*C(1,1));
211 L(5) = (1.0/(det*det))*(C(0,2)*C(2,1)-C(0,1)*C(2,2));
213 M(0) = mu4*sin_plyagl*sin_plyagl+mu5*cos_plyagl*
cos_plyagl;
214 M(1) = mu4*cos_plyagl*cos_plyagl+mu5*sin_plyagl*
sin_plyagl;
218 M(5) = (mu4-
mu5)*cos_plyagl*sin_plyagl;
220 LML(0) = M(0)*
L(0)*
L(0)+2.0*M(4)*
L(0)*
L(4)+2.0*M(5)*
L(0)*
L(5)+M(2)*
L(4)*
L(4)+2.0*M(3)*
L(4)*
L(5)+M(1)*
L(5)*
L(5);
221 LML(1) = M(1)*
L(1)*
L(1)+2.0*M(4)*
L(1)*
L(3)+2.0*M(5)*
L(1)*
L(5)+M(2)*
L(3)*
L(3)+2.0*M(4)*
L(3)*
L(5)+M(0)*
L(5)*
L(5);
222 LML(2) = M(2)*
L(2)*
L(2)+2.0*M(3)*
L(2)*
L(3)+2.0*M(4)*
L(2)*
L(4)+M(1)*
L(3)*
L(3)+2.0*M(5)*
L(3)*
L(4)+M(0)*
L(4)*
L(4);
223 LML(3) =
L(2)*(
L(1)*M(3)+
L(3)*M(2)+
L(5)*M(4))+
L(3)*(
L(1)*M(1)+
L(3)*M(3)+
L(5)*M(5))+
L(4)*(
L(5)*M(0)+
L(1)*M(5)+
L(3)*M(4));
224 LML(4) =
L(2)*(
L(0)*M(4)+
L(4)*M(2)+
L(5)*M(3))+
L(3)*(
L(0)*M(5)+
L(5)*M(1)+
L(4)*M(3))+
L(4)*(
L(0)*M(0)+
L(4)*M(4)+
L(5)*M(5));
225 LML(5) =
L(1)*(
L(0)*M(5)+
L(5)*M(1)+
L(4)*M(3))+
L(3)*(
L(0)*M(4)+
L(4)*M(2)+
L(5)*M(3))+
L(5)*(
L(0)*M(0)+
L(4)*M(4)+
L(5)*M(5));
227 eye(0) = 1.0; eye(1) = 1.0; eye(2) = 1.0; eye(3) = 0.0; eye(4) = 0.0; eye(5) = 0.0;
229 double I1 = C(0,0) + C(1,1) + C(2,2);
230 double II1 = C(0,0)*C(0,0) + C(1,1)*C(1,1) + C(2,2)*C(2,2) + 2.0*C(1,2)*C(1,2) + 2.0*C(0,2)*C(0,2) + 2.0*C(0,1)*C(0,1);
231 double I2 = (1.0/2.0)*(I1*I1-II1);
233 double I4 = C(0,0)*M(0) + C(1,1)*M(1) + C(2,2)*M(2) + 2.0*C(0,1)*M(5) + 2.0*C(0,2)*M(4) + 2.0*C(1,2)*M(3);
234 double I5 = CC(0,0)*M(0) + CC(1,1)*M(1) + CC(2,2)*M(2) + 2.0*CC(0,1)*M(5) + 2.0*CC(0,2)*M(4) + 2.0*CC(1,2)*M(3);
235 double J5 = I5 - I1*I4 + I2*
trm;
236 double pI3 = std::pow(I3,-beta3);
237 double pI4 = std::pow(I4, beta4);
238 double pJ5 = std::pow(J5, beta5);
240 for (
unsigned int i=0;
i<6; ++
i){
241 dJ5(
i) = J5*
L(
i) - I3*LML(
i);
242 piola_stress(
i) = 2.0*mu1*eye(
i) + 2.0*mu2*(I1*eye(
i)-c(
i)) + (2.0*mu3*det*det-mu)*
L(
i)
262 tensor_product(4.0*beta4*pI4/(I4*ptrmbeta4),M,M,tangent_piola,1.0);
263 tensor_product(4.0*beta5*pJ5/(J5*ptrmbeta5),dJ5,dJ5,tangent_piola,1.0);
268 ddJ5.Scale(4.0*pJ5/ptrmbeta5);
269 tangent_piola += ddJ5;
272 Epetra_SerialDenseVector
get_neumannBc(Epetra_SerialDenseMatrix & matrix_X, Epetra_SerialDenseMatrix & xg,
unsigned int & gp){
273 Epetra_SerialDenseVector h(3);
274 std::cout <<
"**Err: Not using that method in this example!\n";
278 Epetra_SerialDenseVector
get_forcing(
double & x1,
double & x2,
double & x3,
unsigned int & e_lid,
unsigned int & gp){
279 Epetra_SerialDenseVector f(3);
280 std::cout <<
"**Err: Not using that method in this example!\n";
287 void get_stress_for_recover(Epetra_SerialDenseMatrix & deformation_gradient,
double & det, Epetra_SerialDenseMatrix & piola_stress){
void get_constitutive_tensors(Epetra_SerialDenseMatrix &deformation_gradient, Epetra_SerialDenseVector &piola_stress, Epetra_SerialDenseMatrix &tangent_piola)
std::vector< int > local_nodes
void tensor_product(double scalarAB, Epetra_SerialDenseVector &A, Epetra_SerialDenseVector &B, Epetra_SerialDenseMatrix &AoB, double scalarThis)
void get_material_parameters(unsigned int &e_lid, unsigned int &gp)
void get_stress_for_recover(Epetra_SerialDenseMatrix &deformation_gradient, double &det, Epetra_SerialDenseMatrix &piola_stress)
std::vector< double > nodes_coord
~compressibleHyperelasticity_linearPatchTest()
Epetra_SerialDenseVector N_cells
Epetra_SerialDenseVector exactSolution(double &x1, double &x2, double &x3)
Epetra_SerialDenseVector detJac_cells
void get_matrix_and_rhs(Epetra_Vector &x, Epetra_FECrsMatrix &K, Epetra_FEVector &F)
std::vector< int > local_cells
clc clear all close all mesh
void assemblePureDirichlet_homogeneousForcing(Epetra_Vector &x, Epetra_FECrsMatrix &K, Epetra_FEVector &F)
int n_local_nodes_without_ghosts
std::vector< int > cells_nodes
void set_parameters(Epetra_SerialDenseVector &x, double &angle)
Epetra_Import * ImportToOverlapMap
void apply_dirichlet_conditions(Epetra_FECrsMatrix &K, Epetra_FEVector &F, double &displacement)
unsigned int n_gauss_cells
Epetra_SerialDenseVector get_forcing(double &x1, double &x2, double &x3, unsigned int &e_lid, unsigned int &gp)
void sym_tensor_product(double scalarAB, Epetra_SerialDenseVector &A, Epetra_SerialDenseVector &B, Epetra_SerialDenseMatrix &AoB, double scalarThis)
Epetra_SerialDenseVector get_neumannBc(Epetra_SerialDenseMatrix &matrix_X, Epetra_SerialDenseMatrix &xg, unsigned int &gp)
void setup_dirichlet_conditions()
compressibleHyperelasticity_linearPatchTest(Epetra_Comm &comm, Teuchos::ParameterList &Parameters)
Epetra_SerialDenseVector gauss_weight_cells
double errorL2(Epetra_Vector &uStandardMap)
void get_material_parameters_for_recover(unsigned int &e_lid)