-
Notifications
You must be signed in to change notification settings - Fork 1
/
binary.py
executable file
·325 lines (265 loc) · 9.03 KB
/
binary.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
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
# Importing relevant packages
# random is needed for noise generation in initialisation of system
# dolfin is the finite element solver package within the fenics environment
import random
from dolfin import *
import csv
import os
import time
import sys
from FHTaylor import taylorapprox_fullFH, taylorapprox_logonlyFH, heatofmixing, symquarticdoublewell, polynomialfit_fullFH
from Thermo.Thermo import RK, ThermoMix
from parameters.params import (
A_RAW,
NOISE_MAGNITUDE,
TIME_MAX,
DT,
N_CELLS,
DOMAIN_LENGTH,
theta_ch,
MESH_TYPE,
TIME_STRIDE,
SPECIES,
chi_AB,
N_A,
N_B,
GIBBS,
TEMPERATURE,
FINITE_ELEMENT,
FINITE_ELEMENT_ORDER,
SOLVER_CONFIG,
MOBILITY_MODEL,
SIZE_DISPARITY,
A_SYM
)
class CahnHilliardEquation(NonlinearProblem):
def __init__(self, a, L):
NonlinearProblem.__init__(self)
self.L = L
self.a = a
# self.reset_sparsity = True
def F(self, b, x):
assemble(self.L, tensor=b)
def J(self, A, x):
assemble(self.a, tensor=A)#, reset_sparsity=self.reset_sparsity)
# self.reset_sparsity = False
parameters["form_compiler"]["optimize"] = True
parameters["form_compiler"]["cpp_optimize"] = True
parameters["form_compiler"]["quadrature_degree"] = 2
# INITIAL AND BOUNDARY CONDITIONS OF THE PROBLEM #
# Initial conditions which include random noise as a perturbation
class InitialConditions(UserExpression):
def __init__(self, **kwargs):
random.seed(2 + MPI.rank(MPI.comm_world))
super().__init__(**kwargs)
def eval(self, values, x):
# [0] corresponds to the concentration field for species A
# [1] coresponds to the \mu_{AB} field
values[0] = A_RAW + 2.0 * NOISE_MAGNITUDE * (0.5 - random.random())
values[1] = 0.0
def value_shape(self):
return (2,)
# Setting up Mesh and Periodic Boundary Conditions #
## Setting up periodic boundary conditions.
## We are creating classes for defining the various parts of the boundaries (top, bottom, left and right)
class LeftBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and near(x[0], 0.0)
class RightBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and near(x[0], DOMAIN_LENGTH)
class BottomBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and near(x[1], 0.0)
class TopBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and near(x[1], DOMAIN_LENGTH)
class PeriodicBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and (near(x[0], 0.0))
# Map RightBoundary to LeftBoundary
def map(self, x, y):
y[0] = x[0] - DOMAIN_LENGTH
y[1] = x[1]
N = int(N_CELLS)
if MESH_TYPE == "structured":
mesh = RectangleMesh(
Point(0.0, 0.0), Point(DOMAIN_LENGTH, DOMAIN_LENGTH), N, N
)
else:
domain_vertices = [
Point(0.0, 0.0),
Point(DOMAIN_LENGTH, 0.0),
Point(DOMAIN_LENGTH, DOMAIN_LENGTH),
Point(0.0, DOMAIN_LENGTH),
]
domain = Polygon(domain_vertices)
mesh = generate_mesh(domain, N)
# CG stands for continuous galerkin can is a lagrange type element.
# The "1" corresponds to the order. So this is a linear lagrange element.
# CH = FunctionSpace(mesh, "CG", 1, constrained_domain=PeriodicBoundary())
# P1 = FiniteElement("Lagrange", interval, 1)
# # # This function has been deprecated for 2016.1.0, but still works.
# # # The requirement is to created a combined concentration and chemical potential function space
# # # Since we are solving for two variables.
P1 = FiniteElement(FINITE_ELEMENT, mesh.ufl_cell(), FINITE_ELEMENT_ORDER)
CH = FunctionSpace(mesh, P1*P1)
dch = TrialFunction(CH)
h_1, j_1 = TestFunctions(CH)
ch = Function(CH)
ch0 = Function(CH)
# Split mixed functions
da, dmu_AB = split(dch)
x_a, mu_AB = split (ch)
x_a0, mu0_AB = split(ch0)
x_a = variable(x_a)
ch_init = InitialConditions(degree=1)
ch.interpolate(ch_init)
ch0.interpolate(ch_init)
if GIBBS == "FH":
# r = RK("FH",["PB","PS"],[N_A,N_B])
g = ( x_a * ln(x_a) )/N_A + ((1.0-x_a)*ln(1-x_a)/ N_B) + x_a*(1.0-x_a)*chi_AB
# g = r.G_RK(x_a)
print("Vanilla FH")
elif GIBBS == "UNIFAC":
r = RK("UNIFAC",SPECIES,[N_A,N_B],[TEMPERATURE])
g = r.G_RK(x_a)
print("Redlich-Kister UNIFAC")
elif GIBBS == "PCSAFT_CR":
r = RK("PCSAFT",SPECIES,[N_A,N_B],[TEMPERATURE],CR="On")
g = r.G_RK(x_a)
print("Redlich-Kister PCSAFT")
elif GIBBS == "PCSAFT_Fit":
r = RK("PCSAFT",SPECIES,[N_A,N_B],[TEMPERATURE],CR="Off")
g = r.G_RK(x_a)
print("Redlich-Kister PCSAFT")
elif GIBBS == "TaylorApproxFullFH":
g = taylorapprox_fullFH(N_A, N_B, chi_AB, x_a)
print("full taylor approx of FH")
elif GIBBS == "TaylorApproxLogOnlyFH":
g = taylorapprox_logonlyFH(N_A, N_B, chi_AB, x_a)
print ("Taylor approx of log term in FH only")
elif GIBBS == "polynomialfit_fullFH":
g = polynomialfit_fullFH(N_A, N_B, chi_AB, x_a, 4)
print("Full curve fit ahaha")
elif GIBBS == "symquarticdoublewell_fullFH":
g = symquarticdoublewell(x_a, A_SYM)
print ("BS quartic polynomial")
elif GIBBS == "Heatofmixing":
g = heatofmixing(chi_AB, x_a)
else:
print ("work harder")
if SIZE_DISPARITY == "SMALL":
if GIBBS=="FH":
kappa = (2.0/3.0)*chi_AB
else:
chi_AB = r.chi()
kappa = (2.0/3.0)*chi_AB[0]
print ("about the same size")
elif SIZE_DISPARITY == "LARGE":
if GIBBS=="FH":
kappa = (1.0/3.0)*chi_AB
else:
chi_AB = r.chi()
kappa = (1.0/3.0)*chi_AB[0]
print ("big size difference")
print(kappa)
# print(2.0/3.0*r.RK(0)[0]/(N_A*N_B)**0.5)
# Using the fenics autodifferentiation toolkit
dgdx_a = diff(g,x_a)
mu_AB_mid = (1.0 - theta_ch) * mu0_AB + theta_ch * mu_AB
dt = DT
F_a = (
x_a * h_1 * dx
- x_a0 * h_1 * dx
+ dt * x_a * (1.0 - x_a) * dot(grad(mu_AB_mid), grad(h_1)) * dx
)
F_a_constant = (
x_a * h_1 * dx
- x_a0 * h_1 * dx
+ dt * dot(grad(mu_AB_mid), grad(h_1)) * dx
)
F_mu_AB = (
mu_AB * j_1 * dx
- dgdx_a * j_1 * dx
- kappa * dot(grad(x_a), grad(j_1)) * dx
)
if MOBILITY_MODEL == "Variable":
F = F_a + F_mu_AB
print ("work hard model")
elif MOBILITY_MODEL == "Constant":
F = F_a_constant + F_mu_AB
print("less work hehe")
else:
print("wrong model implemented")
sys.exit()
#Compute directional derivative about u in the direction of du
a = derivative(F, ch, dch)
if SOLVER_CONFIG == "LU":
problem = CahnHilliardEquation(a, F)
# problem.set_bounds(x_a_min, x_a_max)
solver = NewtonSolver()
solver.parameters["linear_solver"] = "lu"
# solver.parameters["linear_solver"] = "gmres"
# solver.parameters["preconditioner"] = "ilu"
solver.parameters["convergence_criterion"] = "residual"
solver.parameters["relative_tolerance"] = 1e-10
solver.parameters["absolute_tolerance"] = 1e-16
# solver.parameters["relaxation_parameter"] = 1.0
elif SOLVER_CONFIG == "KRYLOV":
class CustomSolver(NewtonSolver):
def __init__(self):
NewtonSolver.__init__(self, mesh.mpi_comm(),
PETScKrylovSolver(), PETScFactory.instance())
def solver_setup(self, A, P, problem, iteration):
self.linear_solver().set_operator(A)
PETScOptions.set("ksp_type", "gmres")
PETScOptions.set("ksp_monitor")
PETScOptions.set("pc_type", "hypre")
PETScOptions.set("pc_hypre_type", "euclid")
PETScOptions.set("ksp_rtol", "1.0e-8")
PETScOptions.set("ksp_atol", "1.0e-16")
PETScOptions.set('ksp_max_it', '1000')
self.linear_solver().set_from_options()
problem = CahnHilliardEquation(a, F)
solver = CustomSolver()
# Initialising the output files
gibbs_list = []
# file = XDMFFile("output.xdmf")
file = File("output.pvd", "compressed")
start = time.time()
t = 0.0
time_stride = TIME_STRIDE
timestep = 0
while (t < TIME_MAX):
t += dt
print (f"Time = {t}")
ch0.vector()[:] = ch.vector()
solver.solve(problem, ch.vector())
timestep += 1
# for i in ch.vector()[::2]:
# if i < 0.0:
# i == 0.0 + 1e-16
# elif i > 1.0:
# i == 1.0 - 1e-16
# Assembling the various terms of the Landau-Ginzburg free energy functional
homogenous_energy = assemble(g * dx())
gradient_energy = assemble(Constant(0.5)*kappa * dot(grad(x_a), grad(x_a)) * dx())
gibbs= homogenous_energy + gradient_energy
gibbs_list.append(gibbs)
fpath = "./output_gibbs.csv"
headers = ["time", "gibbs"]
end = time.time()
print(end-start)
# Write header row (for first timestep)
if not os.path.exists(fpath):
with open(fpath,"w") as f:
w = csv.DictWriter(f,headers)
w.writeheader()
# Appending case data
with open(fpath, "a") as f:
w = csv.DictWriter(f, headers)
w.writerow({"time": float(t), "gibbs": float(gibbs)})
if (timestep % time_stride ==0):
# file.write (ch.split()[0], t)
file << (ch.split()[0], t)