-
Notifications
You must be signed in to change notification settings - Fork 0
/
CHiMaD3_DK_w_postprocess.i
322 lines (290 loc) · 6.22 KB
/
CHiMaD3_DK_w_postprocess.i
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
#
# CHiMaD benchmark problem #3, dendritic growth
# Dong-Uk Kim
# as of 06/20/2018 written
# as of 07/16/2018 revised
[Mesh]
type = GeneratedMesh
dim = 2
nx = 75
ny = 75
xmax = 960
ymax = 960
[]
[Variables]
[./p] #phase-field variable
order = FIRST
family = LAGRANGE
[../]
[./u] #dimensionless temperature variable
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./dpx]
order = FIRST
family = MONOMIAL
[../]
[./dpy]
order = FIRST
family = MONOMIAL
[../]
[./FEdensity]
order = FIRST
family = MONOMIAL
[../]
[]
[ICs]
[./p_circleIC]
type = SmoothCircleIC
variable = p
x1 = 0
y1 = 0
radius = 8
int_width = 1
invalue = 1.0
outvalue = -1.0
[../]
[./u_constantIC]
type = ConstantIC
variable = u
value = -0.3
[../]
[]
[BCs]
#Do nothing means no-flux boundary condition
[]
[AuxKernels]
[./get_dpx]
type = VariableGradientComponent
variable = dpx
gradient_variable = p
component = x
execute_on = LINEAR
[../]
[./get_dpy]
type = VariableGradientComponent
variable = dpy
gradient_variable = p
component = y
execute_on = LINEAR
[../]
[]
[Kernels]
#Phase-field
[./time_derivative_p]
type = AnisotropicTimeDerivative
variable = p
tau_name = tau_aniso
gradient_component_names = 'dpx dpy'
gradmag_threshold = 1e-4
[../]
[./InterfacialE]
type = AnisotropicGradEnergy
variable = p
mob_name = L
kappa_name = Wsq_aniso
gradient_component_names = 'dpx dpy'
gradmag_threshold = 1e-4
[../]
[./localFE]
type = AllenCahn
variable = p
f_name = f_doublewell
mob_name = L
[../]
[./DrivingForce]
type = ACSwitching
variable = p
Fj_names = 'lambda_U'
hj_names = 'h'
args = 'p u'
mob_name = L
[../]
#Dimensionless temperature field
[./time_derivative_u]
type = TimeDerivative
variable = u
# type = CoefTimeDerivative
# variable = u
# Coefficient = 0.1
[../]
[./div_grad_mu]
type = ACInterface
variable = u
mob_name = L
kappa_name = D
variable_L = false
# type = Diffusion
# variable = u
[../]
[./Partition_term]
type = CoupledSwitchingTimeDerivative
variable = u
v = p
Fj_names = 'oneovertwo'
hj_names = 'switching_U'
args = 'p'
# type = CoefCoupledTimeDerivative
# variable = u
# coef = -0.05
# v = p
[../]
[]
[Materials]
[./f_doublewell]
type = DerivativeParsedMaterial
f_name = f_doublewell
args = 'p'
function = '-0.5*p^2+0.25*p^4'
derivative_order = 2
#outputs = exodus
[../]
[./lambda_U]
type = DerivativeParsedMaterial
f_name = lambda_U
material_property_names = 'lambda'
args = 'u'
function = 'lambda*u'
[../]
[./h]
type = DerivativeParsedMaterial
f_name = h
args = 'p'
function = 'p*(1-2*p^2/3+p^4/5)'
derivative_order = 2
#outputs = exodus
[../]
[./switching_U]
type = DerivativeParsedMaterial
f_name = switching_U
args = 'p'
function = 'p'
derivative_order = 1
#outputs = exodus
[../]
[./constants]
type = GenericConstantMaterial
prop_names = 'tau0 W0 D eps4 L oneovertwo'
prop_values = '1 1 10 0.05 1 -0.5'
[../]
[./lambda]
type = ParsedMaterial
f_name = lambda
material_property_names = 'D tau0 W0'
function = 'D*tau0/0.6267/W0/W0'
[../]
[./Anisotropic_Wsquare]
type = DerivativeParsedMaterial
f_name = Wsq_aniso
material_property_names = 'W0 eps4'
args = 'dpx dpy'
function = 'if(sqrt(dpx^2 + dpy^2) > 1e-5, W0^2 * (1 + eps4 * (dpx^4 + dpy^4 - 6*dpx^2*dpy^2)/(dpx^2 + dpy^2)^2)^2, W0^2)'
derivative_order = 2
#outputs = exodus
[../]
[./Anisotropic_tau]
type = DerivativeParsedMaterial
f_name = tau_aniso
material_property_names = 'tau0 eps4'
args = 'dpx dpy'
function = 'if(sqrt(dpx^2 + dpy^2) > 1e-5, tau0 * (1 + eps4 * (dpx^4 + dpy^4 - 6*dpx^2*dpy^2)/(dpx^2 + dpy^2)^2)^2, tau0)'
derivative_order = 2
#outputs = exodus
[../]
# #To prevent divided by zero
# [./eps4]
# type = ParsedMaterial
# f_name = eps4
# material_property_names = 'eps0'
# #args = 'p'
# #function = 'if((0.99-p)*(0.99+p) >= 0, eps0, 0)'
# args = 'dpx dpy'
# function = 'if(sqrt(dpx^2 + dpy^2) > 1e-5, eps0, 0)'
# #outputs = exodus
# [../]
#For the postprocess
[./solid_volume]
type = ParsedMaterial
f_name = solid_fraction_per_element
args = 'p'
function = '(1 + p)/2/960/960'
[../]
[./FEdensity_material]
type = ParsedMaterial
f_name = f_density
material_property_names = 'Wsq_aniso f_doublewell lambda_U h'
args = 'p u dpx dpy'
function = '0.5*Wsq_aniso*(dpx^2+dpy^2) + f_doublewell + lambda_U*h'
#outputs = exodus
[../]
[]
[Preconditioning]
#active = ''
[./cw_coupling]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = PJFNK
scheme = bdf2
# petsc_options_iname = '-pc_type -sub_pc_type'
# petsc_options_value = 'asm lu'
# petsc_options_iname = '-pc_type -pc_asm_overlap'
# petsc_options_value = 'asm 1'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
l_max_its = 20
l_tol = 1e-4
nl_max_its = 20
nl_rel_tol = 1e-8
nl_abs_tol = 1e-11
[./TimeStepper]
type = IterationAdaptiveDT
dt = 0.003
growth_factor = 1.2
cutback_factor = 0.8
#optimal_iterations = 4
#iteration_window = 4
[../]
dtmax = 0.3
end_time = 1500
[./Adaptivity]
initial_adaptivity = 5
max_h_level = 5
refine_fraction = 0.95
coarsen_fraction = 0.10
weight_names = 'p u'
weight_values = '1.0 1.0'
[../]
[]
[Outputs]
exodus = true
csv = true
print_perf_log = true
[]
# [Debug]
# show_var_residual_norms = true
# []
[Postprocessors]
[./Total_solid_fraction]
type = ElementIntegralMaterialProperty
mat_prop = solid_fraction_per_element
[../]
[./Total_FE]
type = ElementIntegralMaterialProperty
mat_prop = f_density
[../]
[./Interface_location_along_x_axis]
type = FindValueOnLine
start_point = '0 0 0'
end_point = '960 0 0'
target = 0
depth = 15
tol = 1e-1
v = p
[../]
[]