-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathElastic_complex.py
181 lines (138 loc) · 5.06 KB
/
Elastic_complex.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
from ufl import *
from dolfin import *
from multiphenics import *
import matplotlib.pyplot as plt
import numpy
a=0.1;
t=a/10;
####PML definition
############################################################
def compProd(a1,b1,a2,b2):
x1=a1*a2-b1*b2;
x2=a1*b2+a2*b1;
return (x1, x2)
class LR(UserExpression):
def eval(self, value, x):
a=0.1;t=a/10;
b1 = numpy.abs(x[0])-a/2; f1=10*(b1+numpy.abs(b1))/(2*t);
b2 = numpy.abs(x[1])-a/2; f2=10*(b2+numpy.abs(b2))/(2*t);
l1R=1+f1;l1I=-f1;
l2R=1+f2;l2I=-f2;
value[0] = l1R/(pow(l1R,2)+pow(l1I,2));
value[3] = l2R/(pow(l2R,2)+pow(l2I,2));
value[1] = 0;
value[2] = 0;
def value_shape(self):
return (2,2)
class LI(UserExpression):
def eval(self, value, x):
a=0.1;t=a/10;
b1 = numpy.abs(x[0])-a/2; f1=10*(b1+numpy.abs(b1))/(2*t);
b2 = numpy.abs(x[1])-a/2; f2=10*(b2+numpy.abs(b2))/(2*t);
l1R=1+f1;l1I=-f1;
l2R=1+f2;l2I=-f2;
value[0] = -l1I/(pow(l1R,2)+pow(l1I,2));
value[3] = -l2I/(pow(l2R,2)+pow(l2I,2));
value[1] = 0;
value[2] = 0;
def value_shape(self):
return (2,2)
class lR(UserExpression):
def eval(self, value, x):
a=0.1;t=a/10;
b1 = numpy.abs(x[0])-a/2; f1=10*(b1+numpy.abs(b1))/(2*t);
b2 = numpy.abs(x[1])-a/2; f2=10*(b2+numpy.abs(b2))/(2*t);
l1R=1+f1;l1I=-f1;
l2R=1+f2;l2I=-f2;
value[0], temp = compProd(l1R,l1I,l2R,l2I);
# def value_shape(self):
# return (1,)
class lI(UserExpression):
def eval(self, value, x):
a=0.1;t=a/10;
b1 = numpy.abs(x[0])-a/2; f1=10*(b1+numpy.abs(b1))/(2*t);
b2 = numpy.abs(x[1])-a/2; f2=10*(b2+numpy.abs(b2))/(2*t);
l1R=1+f1;l1I=-f1;
l2R=1+f2;l2I=-f2;
temp, value[0] = compProd(l1R,l1I,l2R,l2I);
# def value_shape(self):
# return (1,)
class Middle(SubDomain):
def inside(self, x, on_boundary):
return (near(x[1], 0) and between(x[0], (-1/1500, 1/1500)))
def OuterBoundaries(x, on_boundary):
return on_boundary
LambdaR = LR()
LambdaI = LI()
LLR = lR()
LLI = lI()
####Geometry and mesh
############################################################
nx, ny = 200, 200
mesh = RectangleMesh(Point(-a/2-t, -a/2-t), Point(a/2+t, a/2+t), nx, ny)
boundaries = MeshFunction("size_t", mesh, 1)
boundaries.set_all(0)
Middle().mark(boundaries, 1)
####Material properties
############################################################
E, rho, nu = 300E+8, 8000, 0.3
frequency=50000
omega=2*numpy.pi*frequency
lam = E*nu/((1+nu)*(1-2*nu))
mu = E/(2*(1+nu))
i,j,k,l = indices(4)
delta = Identity(2)
C=as_tensor((lam*(delta[i,j]*delta[k,l])+mu*(delta[i,k]*delta[j,l]+delta[i,l]*delta[j,k])),(i,j,k,l))
####Space of functions
############################################################
V = VectorFunctionSpace(mesh, "CG", 1)
Vcomplex = BlockFunctionSpace([V, V], restrict=[None, None])
#Applying the Boundary conditions
###########################################################
BC1 = DirichletBC(Vcomplex.sub(0), Constant((1.0, 0.0)), boundaries, 1)
BC2 = DirichletBC(Vcomplex.sub(1), Constant((0.0, 0.0)), boundaries, 1)
#BC3 = DirichletBC(Vcomplex.sub(0), Constant((0.0, 0.0)), OuterBoundaries)
#BC4 = DirichletBC(Vcomplex.sub(1), Constant((0.0, 0.0)), OuterBoundaries)
bcs = BlockDirichletBC([BC1, BC2])
#Real and Imaginary parts of the trial and test functions
###########################################################
(uR, uI) = BlockTrialFunction(Vcomplex)
(wR, wI) = BlockTestFunction(Vcomplex)
#Weak form and solve
###########################################################
strainR=0.5*(as_tensor((Dx(uR[i],k)*LambdaR[k,j]),(i,j)) + as_tensor((Dx(uR[k],i)*LambdaR[j,k]),(i,j))) - 0.5*(as_tensor((Dx(uI[i],k)*LambdaI[k,j]),(i,j)) + as_tensor((Dx(uI[k],i)*LambdaI[j,k]),(i,j)))
strainI=0.5*(as_tensor((Dx(uR[i],k)*LambdaI[k,j]),(i,j)) + as_tensor((Dx(uR[k],i)*LambdaI[j,k]),(i,j))) + 0.5*(as_tensor((Dx(uI[i],k)*LambdaR[k,j]),(i,j)) + as_tensor((Dx(uI[k],i)*LambdaR[j,k]),(i,j)))
tempR=LLR*LambdaR-LLI*LambdaI
tempI=LLR*LambdaI+LLI*LambdaR
sR=as_tensor((C[i,j,k,l]*strainR[k,l]),[i,j])
sI=as_tensor((C[i,j,k,l]*strainI[k,l]),[i,j])
stressR=as_tensor((sR[i,j]*tempR[j,k]),(i,k)) - as_tensor((sI[i,j]*tempI[j,k]),(i,k))
stressI=as_tensor((sR[i,j]*tempI[j,k]),(i,k)) + as_tensor((sI[i,j]*tempR[j,k]),(i,k))
a_11 = omega**2*rho*dot(LLR*uR,wR)*dx
a_12 = -omega**2*rho*dot(LLI*uI,wR)*dx
a_13 = -inner(stressR, grad(wR))*dx
a_14 = 0
a_21 = omega**2*rho*dot(LLI*uR,wI)*dx
a_22 = -omega**2*rho*dot(LLR*uI,wI)*dx
a_23 = 0
a_24 = -inner(stressI, grad(wI))*dx
a = [[a_11, a_12, a_13, a_14],
[a_21, a_22, a_23, a_24]]
L = [0, 0]
L = BlockForm(L, block_function_space=Vcomplex, block_form_rank=1)
A = block_assemble(a)
F = block_assemble(L)
bcs.apply(A)
bcs.apply(F)
sol = BlockFunction(Vcomplex)
block_solve(A, sol.block_vector(), F, "mumps")
uR, uI = block_split(sol)
plt.figure()
plot(uR)
plt.figure()
plot(uI)
plt.show()
file_uR = File("Elastic/uR.pvd")
file_uI = File("Elastic/uI.pvd")
file_uR << uR
file_uI << uI