-
Notifications
You must be signed in to change notification settings - Fork 1
/
SimplePatchPowell.py
274 lines (236 loc) · 11.6 KB
/
SimplePatchPowell.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
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
import sys
import numpy as np
from scipy import optimize as opt
import gmsh
from mpi4py import MPI as nMPI
import meshio
import dolfin
x_old = np.array([0.0, 0.0, 0.0])
rho_old = 0.0
def Cost(xp):
global x_old, rho_old
comm = nMPI.COMM_WORLD
mpi_rank = comm.Get_rank()
x1, x2 = xp
rs = 8.0 # radiation boundary radius
l = x1 # Patch length
w = 4.5 # Patch width
s1 = x2 * x1 / 2.0 # Feed offset
h = 1.0 # Patch height
t = 0.05 # Metal thickness
lc = 1.0 # Coax length
rc = 0.25 # Coax shield radius
cc = 0.107 #Coax center conductor 50 ohm air diel
eps = 1.0e-4
tol = 1.0e-6
eta = 377.0
eps_c = 1.0
k0 = 2.45 * 2.0 * np.pi / 30.0 # Frequency in GHz
ls = 0.025
lm = 0.8
lw = 0.06
lp = 0.3
if mpi_rank == 0:
print("x[0] = {0:<f}, x[1] = {1:<f} ".format(xp[0], xp[1]))
print("length = {0:<f}, width = {1:<f}, feed offset = {2:<f}".format(l, w, s1))
gmsh.initialize()
gmsh.option.setNumber('General.Terminal', 1)
gmsh.model.add("SimplePatchOpt")
# Radiation sphere
gmsh.model.occ.addSphere(0.0, 0.0, 0.0, rs, 1)
gmsh.model.occ.addBox(0.0, -rs, 0.0, rs, 2*rs, rs, 2)
gmsh.model.occ.intersect([(3,1)],[(3,2)], 3, removeObject=True, removeTool=True)
# Patch
gmsh.model.occ.addBox(0.0, -l/2, h, w/2, l, t, 4)
# coax center
gmsh.model.occ.addCylinder(0.0, s1, -lc, 0.0, 0.0, lc+h, cc, 5, 2.0*np.pi)
# patch center gnd via
#gmsh.model.occ.addCylinder(0.0, 0.0, 0.0, 0.0, 0.0, h, cc, 6, 2.0*np.pi)
# coax shield
gmsh.model.occ.addCylinder(0.0, s1, -lc, 0.0, 0.0, lc, rc, 7)
gmsh.model.occ.addBox(0.0, s1-rc, -lc, rc, 2.0*rc, lc, 8)
gmsh.model.occ.intersect([(3,7)], [(3,8)], 9, removeObject=True, removeTool=True)
gmsh.model.occ.fuse([(3,3)], [(3,9)], 10, removeObject=True, removeTool=True)
# cutout internal boundaries
#gmsh.model.occ.cut([(3,10)], [(3,4),(3,5),(3,6)], 11, removeObject=True, removeTool=True)
gmsh.model.occ.cut([(3,10)], [(3,4),(3,5)], 11, removeObject=True, removeTool=True)
gmsh.option.setNumber('Mesh.MeshSizeMin', ls)
gmsh.option.setNumber('Mesh.MeshSizeMax', lm)
gmsh.option.setNumber('Mesh.Algorithm', 6)
gmsh.option.setNumber('Mesh.Algorithm3D', 1)
# gmsh.option.setNumber('Mesh.MshFileVersion', 2.2)
gmsh.option.setNumber('Mesh.MshFileVersion', 4.1)
gmsh.option.setNumber('Mesh.Format', 1)
gmsh.option.setNumber('Mesh.MinimumCirclePoints', 36)
gmsh.option.setNumber('Mesh.CharacteristicLengthFromCurvature', 1)
gmsh.model.occ.synchronize()
pts = gmsh.model.getEntities(0)
gmsh.model.mesh.setSize(pts, lm) #Set background mesh density
pts = gmsh.model.getEntitiesInBoundingBox(-eps, -l/2-eps, h-eps, w/2+eps, l/2+eps, h+t+eps)
gmsh.model.mesh.setSize(pts, ls)
pts = gmsh.model.getEntitiesInBoundingBox(-eps, s1-rc-eps, -lc-eps, rc+eps, s1+rc+eps, h+eps)
gmsh.model.mesh.setSize(pts, lw)
pts = gmsh.model.getEntitiesInBoundingBox(-eps, -rc-eps, -eps, rc+eps, rc+eps, eps)
gmsh.model.mesh.setSize(pts, lw)
# Embed points to reduce mesh density on patch faces
fce1 = gmsh.model.getEntitiesInBoundingBox(-eps, -l/2-eps, h+t-eps, w/2+eps, l/2+eps, h+t+eps, 2)
gmsh.model.occ.synchronize()
gmsh.model.geo.addPoint(w/4, -l/4, h+t, lp, 1000)
gmsh.model.geo.addPoint(w/4, 0.0, h+t, lp, 1001)
gmsh.model.geo.addPoint(w/4, l/4, h+t, lp, 1002)
gmsh.model.geo.synchronize()
# sf1 = gmsh.model.occ.fragment(fce1, [(0,1000),(0,1001), (0,1002)], -1, True, True)
gmsh.model.occ.synchronize()
print(fce1)
fce2 = gmsh.model.getEntitiesInBoundingBox(-eps, -l/2-eps, h-eps, w/2+eps, l/2+eps, h+eps, 2)
gmsh.model.geo.addPoint(w/4, -9*l/32, h, lp, 1003)
gmsh.model.geo.addPoint(w/4, 0.0, h, lp, 1004)
gmsh.model.geo.addPoint(w/4, 9*l/32, h, lp, 1005)
gmsh.model.geo.synchronize()
for tt in fce1:
gmsh.model.mesh.embed(0, [1000, 1001, 1002], 2, tt[1])
for tt in fce2:
gmsh.model.mesh.embed(0, [1003, 1004, 1005], 2, tt[1])
# sf2 = gmsh.model.occ.fragment(fce2, [(0,1003),(0,1004), (0,1005)], -1, True, True)
print(fce2)
gmsh.model.occ.remove(fce1)
gmsh.model.occ.remove(fce2)
gmsh.model.occ.synchronize()
gmsh.model.addPhysicalGroup(3, [11], 1)
gmsh.model.setPhysicalName(3, 1, "Air")
gmsh.model.mesh.optimize("Relocate3D", niter=5)
gmsh.model.mesh.generate(3)
gmsh.write("SimplePatch.msh")
# gmsh.fltk.run()
gmsh.finalize()
msh = meshio.read("SimplePatch.msh")
for cell in msh.cells:
if cell.type == "tetra":
tetra_cells = cell.data
for key in msh.cell_data_dict["gmsh:physical"].keys():
if key == "tetra":
tetra_data = msh.cell_data_dict["gmsh:physical"][key]
tetra_mesh = meshio.Mesh(points=msh.points, cells={"tetra": tetra_cells},
cell_data={"VolumeRegions":[tetra_data]})
meshio.write("mesh.xdmf", tetra_mesh)
mesh = dolfin.Mesh()
with dolfin.XDMFFile("mesh.xdmf") as infile:
infile.read(mesh)
mvc = dolfin.MeshValueCollection("size_t", mesh, 3)
with dolfin.XDMFFile("mesh.xdmf") as infile:
infile.read(mvc, "VolumeRegions")
cf = dolfin.cpp.mesh.MeshFunctionSizet(mesh, mvc)
class PEC(dolfin.SubDomain):
def inside(self, x, on_boundary):
return on_boundary
class InputBC(dolfin.SubDomain):
def inside(self, x, on_boundary):
return on_boundary and dolfin.near(x[2], -lc, tol)
class OutputBC(dolfin.SubDomain):
def inside(self, x, on_boundary):
rr = np.sqrt(x[0]*x[0]+x[1]*x[1]+x[2]*x[2])
return on_boundary and dolfin.near(rr, 8.0, 1.0e-1)
class PMC(dolfin.SubDomain):
def inside(self, x, on_boundary):
return on_boundary and dolfin.near(x[0], 0.0, tol)
# Volume domains
dolfin.File("VolSubDomains.pvd").write(cf)
dolfin.File("Mesh.pvd").write(mesh)
# Mark boundaries
sub_domains = dolfin.MeshFunction("size_t", mesh, mesh.topology().dim() - 1)
sub_domains.set_all(4)
pec = PEC()
pec.mark(sub_domains, 0)
in_port = InputBC()
in_port.mark(sub_domains, 1)
out_port = OutputBC()
out_port.mark(sub_domains, 2)
pmc = PMC()
pmc.mark(sub_domains, 3)
dolfin.File("BoxSubDomains.pvd").write(sub_domains)
# Set up function spaces
cell = dolfin.tetrahedron
ele_type = dolfin.FiniteElement('N1curl', cell, 2, variant="integral") # H(curl) element for EM
V2 = dolfin.FunctionSpace(mesh, ele_type * ele_type)
V = dolfin.FunctionSpace(mesh, ele_type)
(u_r, u_i) = dolfin.TrialFunctions(V2)
(v_r, v_i) = dolfin.TestFunctions(V2)
dolfin.info(mesh)
#surface integral definitions from boundaries
ds = dolfin.Measure('ds', domain = mesh, subdomain_data = sub_domains)
#volume regions
dx_air = dolfin.Measure('dx', domain = mesh, subdomain_data = cf, subdomain_id = 1)
dx_subst = dolfin.Measure('dx', domain = mesh, subdomain_data = cf, subdomain_id = 2)
# with source and sink terms
u0 = dolfin.Constant((0.0, 0.0, 0.0)) #PEC definition
h_src = dolfin.Expression(('-(x[1] - s) / (2.0 * pi * (pow(x[0], 2.0) + pow(x[1] - s,2.0)))', 'x[0] / (2.0 * pi *(pow(x[0],2.0) + pow(x[1] - s,2.0)))', 0.0), degree = 2, s = s1)
e_src = dolfin.Expression(('x[0] / (2.0 * pi * (pow(x[0], 2.0) + pow(x[1] - s,2.0)))', 'x[1] / (2.0 * pi *(pow(x[0],2.0) + pow(x[1] - s,2.0)))', 0.0), degree = 2, s = s1)
Rrad = dolfin.Expression(('sqrt(x[0] * x[0] + x[1] * x[1] + x[2] * x[2])'), degree = 2)
#Boundary condition dictionary
boundary_conditions = {0: {'PEC' : u0},
1: {'InputBC': (h_src)},
2: {'OutputBC': Rrad}}
n = dolfin.FacetNormal(mesh)
#Build PEC boundary conditions for real and imaginary parts
bcs = []
for i in boundary_conditions:
if 'PEC' in boundary_conditions[i]:
bc = dolfin.DirichletBC(V2.sub(0), boundary_conditions[i]['PEC'], sub_domains, i)
bcs.append(bc)
bc = dolfin.DirichletBC(V2.sub(1), boundary_conditions[i]['PEC'], sub_domains, i)
bcs.append(bc)
# Build input BC source term and loading term
integral_source = []
integrals_load =[]
for i in boundary_conditions:
if 'InputBC' in boundary_conditions[i]:
r = boundary_conditions[i]['InputBC']
bb1 = 2.0 * (k0 * eta) * dolfin.inner(v_i, dolfin.cross(n, r)) * ds(i) #Factor of two from field equivalence principle
integral_source.append(bb1)
bb2 = dolfin.inner(dolfin.cross(n, v_i), dolfin.cross(n, u_r)) * k0 * np.sqrt(eps_c) * ds(i)
integrals_load.append(bb2)
bb2 = dolfin.inner(-dolfin.cross(n, v_r), dolfin.cross(n, u_i)) * k0 * np.sqrt(eps_c) * ds(i)
integrals_load.append(bb2)
for i in boundary_conditions:
if 'OutputBC' in boundary_conditions[i]:
r = boundary_conditions[i]['OutputBC']
bb2 = (dolfin.inner(dolfin.cross(n, v_i), dolfin.cross(n, u_r)) * k0 + 1.0 * dolfin.inner(dolfin.cross(n, v_i), dolfin.cross(n, u_i)) / r)* ds(i)
integrals_load.append(bb2)
bb2 = (dolfin.inner(-dolfin.cross(n, v_r), dolfin.cross(n, u_i)) * k0 + 1.0 * dolfin.inner(dolfin.cross(n, v_r), dolfin.cross(n, u_r)) / r)* ds(i)
integrals_load.append(bb2)
# for PMC, do nothing. Natural BC.
a = (dolfin.inner(dolfin.curl(v_r), dolfin.curl(u_r)) + dolfin.inner(dolfin.curl(v_i), dolfin.curl(u_i)) - eps_c * k0 * k0 * (dolfin.inner(v_r, u_r) + dolfin.inner(v_i, u_i))) * dx_subst + (dolfin.inner(dolfin.curl(v_r), dolfin.curl(u_r)) + dolfin.inner(dolfin.curl(v_i), dolfin.curl(u_i)) - k0 * k0 * (dolfin.inner(v_r, u_r) + dolfin.inner(v_i, u_i))) * dx_air + sum(integrals_load)
L = sum(integral_source)
u1 = dolfin.Function(V2)
vdim = u1.vector().size()
print("Solution vector size =", vdim)
dolfin.solve(a == L, u1, bcs, solver_parameters = {'linear_solver' : 'mumps'})
u1_r, u1_i = u1.split(True)
fp = dolfin.File("EField_r.pvd")
fp << u1_r
fp = dolfin.File("EField_i.pvd")
fp << u1_i
#fp = File('WaveFile.pvd')
#ut = u1_r.copy(deepcopy=True)
#for i in range(50):
# ut.vector().zero()
# ut.vector().axpy(cos(pi * i / 25.0 + pi / 2.0), u1_i.vector())
# ut.vector().axpy(cos(pi * i / 25.0), u1_r.vector())
# fp << (ut, i)
H = dolfin.interpolate(h_src, V) # Get input field
P = dolfin.assemble((-dolfin.dot(u1_r,dolfin.cross(dolfin.curl(u1_i),n))+dolfin.dot(u1_i,dolfin.cross(dolfin.curl(u1_r),n))) * ds(2))
P_refl = dolfin.assemble((-dolfin.dot(u1_i,dolfin.cross(dolfin.curl(u1_r), n)) + dolfin.dot(u1_r, dolfin.cross(dolfin.curl(u1_i), n))) * ds(1))
P_inc = dolfin.assemble((dolfin.dot(H, H) * eta / (2.0 * np.sqrt(eps_c))) * ds(1))
print("Integrated power on port 2:", P/(2.0 * k0 * eta))
print("Incident power at port 1:", P_inc)
print("Integrated reflected power on port 1:", P_inc - P_refl / (2.0 * k0 * eta))
rho_old = (P_inc - P_refl / (2.0 * k0 * eta)) / P_inc #Fraction of incident power reflected as objective function
return rho_old
#Optimization
x0 = np.array([5.0, 0.675])
bnds = [(5.0, 7.0), (0.5, 0.875)]
#res = opt.brute(Cost, bnds, Ns=1, full_output=True)
#res = opt.minimize(Cost, x0, method='Nelder-Mead', options={'maxiter':100, 'disp':True, 'fatol':0.003})
res = opt.minimize(Cost, x0, method='Powell', bounds=bnds, options={'maxiter':100, 'disp':True, 'ftol':0.01})
print(res)
sys.exit(0)