5 #ifndef COMPRESSIBLE_MOONEY_TRANSVERSE_ISOTROPIC_HPP 6 #define COMPRESSIBLE_MOONEY_TRANSVERSE_ISOTROPIC_HPP 22 tiMooney(Epetra_Comm & comm, Teuchos::ParameterList & Parameters){
24 std::string mesh_file = Teuchos::getParameter<std::string>(Parameters.sublist(
"Mesh"),
"mesh_file");
25 Mesh =
new mesh(comm, mesh_file, 1.0);
28 StandardMap =
new Epetra_Map(-1,3*Mesh->n_local_nodes_without_ghosts,&Mesh->local_dof_without_ghosts[0],0,*
Comm);
29 OverlapMap =
new Epetra_Map(-1,3*Mesh->n_local_nodes,&Mesh->local_dof[0],0,*
Comm);
36 for (
unsigned int e=0;
e<Mesh->n_cells/32; ++
e){
37 for (
unsigned int j=0;
j<32; ++
j){
55 mu = 2.0*mu1 + 4.0*mu2 + 2.0*
mu3;
57 ptrmbeta4 = std::pow(trm,beta4);
58 ptrmbeta5 = std::pow(trm,beta5);
67 if (
Comm->MyPID()==0){
74 Comm->Broadcast(&topcoord,1,0);
133 F.Update(-1.0,rhs_dir,1.0);
154 if (phase[e_gid] % 2){
155 cos_plyagl = std::cos(plyagl);
156 sin_plyagl = std::sin(plyagl);
159 cos_plyagl = std::cos(plyagl);
160 sin_plyagl = -std::sin(plyagl);
164 void get_constitutive_tensors(Epetra_SerialDenseMatrix & deformation_gradient, Epetra_SerialDenseVector & piola_stress, Epetra_SerialDenseMatrix & tangent_piola){
166 double det = deformation_gradient(0,0)*deformation_gradient(1,1)*deformation_gradient(2,2)-deformation_gradient(0,0)*deformation_gradient(1,2)*deformation_gradient(2,1)-deformation_gradient(0,1)*deformation_gradient(1,0)*deformation_gradient(2,2)+deformation_gradient(0,1)*deformation_gradient(1,2)*deformation_gradient(2,0)+deformation_gradient(0,2)*deformation_gradient(1,0)*deformation_gradient(2,1)-deformation_gradient(0,2)*deformation_gradient(1,1)*deformation_gradient(2,0);
168 Epetra_SerialDenseMatrix C(3,3);
169 Epetra_SerialDenseMatrix CC(3,3);
170 Epetra_SerialDenseMatrix ddJ5(6,6);
171 Epetra_SerialDenseVector LML(6);
172 Epetra_SerialDenseVector eye(6);
173 Epetra_SerialDenseVector dJ5(6);
174 Epetra_SerialDenseVector M(6);
175 Epetra_SerialDenseVector
L(6);
176 Epetra_SerialDenseVector c(6);
178 C.Multiply(
'T',
'N',1.0,deformation_gradient,deformation_gradient,0.0);
179 CC.Multiply(
'N',
'N',1.0,C,C,0.0);
181 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);
183 L(0) = (1.0/(det*det))*(C(1,1)*C(2,2)-C(1,2)*C(2,1));
184 L(1) = (1.0/(det*det))*(C(0,0)*C(2,2)-C(0,2)*C(2,0));
185 L(2) = (1.0/(det*det))*(C(0,0)*C(1,1)-C(0,1)*C(1,0));
186 L(3) = (1.0/(det*det))*(C(0,2)*C(1,0)-C(0,0)*C(1,2));
187 L(4) = (1.0/(det*det))*(C(0,1)*C(1,2)-C(0,2)*C(1,1));
188 L(5) = (1.0/(det*det))*(C(0,2)*C(2,1)-C(0,1)*C(2,2));
190 M(0) = mu4*sin_plyagl*sin_plyagl+mu5*cos_plyagl*
cos_plyagl;
191 M(1) = mu4*cos_plyagl*cos_plyagl+mu5*sin_plyagl*
sin_plyagl;
195 M(5) = (mu4-
mu5)*cos_plyagl*sin_plyagl;
197 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);
198 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);
199 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);
200 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));
201 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));
202 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));
204 eye(0) = 1.0; eye(1) = 1.0; eye(2) = 1.0; eye(3) = 0.0; eye(4) = 0.0; eye(5) = 0.0;
206 double I1 = C(0,0) + C(1,1) + C(2,2);
207 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);
208 double I2 = (1.0/2.0)*(I1*I1-II1);
210 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);
211 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);
212 double J5 = I5 - I1*I4 + I2*
trm;
213 double pI3 = std::pow(I3,-beta3);
214 double pI4 = std::pow(I4, beta4);
215 double pJ5 = std::pow(J5, beta5);
217 for (
unsigned int i=0;
i<6; ++
i){
218 dJ5(
i) = J5*
L(
i) - I3*LML(
i);
219 piola_stress(
i) = 2.0*mu1*eye(
i) + 2.0*mu2*(I1*eye(
i)-c(
i)) + (2.0*mu3*det*det-mu)*
L(
i)
239 tensor_product(4.0*beta4*pI4/(I4*ptrmbeta4),M,M,tangent_piola,1.0);
240 tensor_product(4.0*beta5*pJ5/(J5*ptrmbeta5),dJ5,dJ5,tangent_piola,1.0);
245 ddJ5.Scale(4.0*pJ5/ptrmbeta5);
246 tangent_piola += ddJ5;
249 Epetra_SerialDenseVector
get_neumannBc(Epetra_SerialDenseMatrix & matrix_X, Epetra_SerialDenseMatrix & xg,
unsigned int & gp){
250 Epetra_SerialDenseVector h(3);
251 std::cout <<
"**Err: Not using that method in this example!\n";
255 Epetra_SerialDenseVector
get_forcing(
double & x1,
double & x2,
double & x3,
unsigned int & e_lid,
unsigned int & gp){
256 Epetra_SerialDenseVector f(3);
257 std::cout <<
"**Err: Not using that method in this example!\n";
264 void get_stress_for_recover(Epetra_SerialDenseMatrix & deformation_gradient,
double & det, Epetra_SerialDenseMatrix & piola_stress){
void setup_dirichlet_conditions()
std::vector< int > local_nodes
void tensor_product(double scalarAB, Epetra_SerialDenseVector &A, Epetra_SerialDenseVector &B, Epetra_SerialDenseMatrix &AoB, double scalarThis)
std::vector< double > nodes_coord
void set_parameters(Epetra_SerialDenseVector &x)
void get_stress_for_recover(Epetra_SerialDenseMatrix &deformation_gradient, double &det, Epetra_SerialDenseMatrix &piola_stress)
Epetra_SerialDenseVector get_forcing(double &x1, double &x2, double &x3, unsigned int &e_lid, unsigned int &gp)
tiMooney(Epetra_Comm &comm, Teuchos::ParameterList &Parameters)
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
Epetra_Import * ImportToOverlapMap
void apply_dirichlet_conditions(Epetra_FECrsMatrix &K, Epetra_FEVector &F, double &displacement)
Epetra_SerialDenseVector get_neumannBc(Epetra_SerialDenseMatrix &matrix_X, Epetra_SerialDenseMatrix &xg, unsigned int &gp)
void get_material_parameters_for_recover(unsigned int &e_lid)
void sym_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_constitutive_tensors(Epetra_SerialDenseMatrix &deformation_gradient, Epetra_SerialDenseVector &piola_stress, Epetra_SerialDenseMatrix &tangent_piola)
void get_matrix_and_rhs(Epetra_Vector &x, Epetra_FECrsMatrix &K, Epetra_FEVector &F)
void set_plyagl(double &Plyagl)