-
-
Notifications
You must be signed in to change notification settings - Fork 25
/
Copy pathjsolve.jl
167 lines (141 loc) · 6.41 KB
/
jsolve.jl
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
#this file contains functions/wrappers related to the solve function in FEniCS
#https://fenicsproject.org/olddocs/dolfin/1.3.0/python/programmers-reference/fem/solving/solve.html
using FEniCS
function solve(A::Matrix, x, b::Matrix, solvers...)
fenics.solve(A.pyobject, x, b.pyobject, solvers...)
end
export solve
#lvsolve is the linear variational solver
function lvsolve(a, L, u; solver_parameters::Dict = Dict("linear_solver" => "default"),
form_compiler_parameters::Dict = Dict("optimize" => true))
fenics.solve(a.pyobject == L.pyobject, u.pyobject,
solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
function lvsolve(a, L, u, bcs = nothing;
solver_parameters::Dict = Dict("linear_solver" => "default"),
form_compiler_parameters::Dict = Dict("optimize" => true))
fenics.solve(a.pyobject == L.pyobject, u.pyobject, bcs = bcs.pyobject,
solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
#allows BoundaryCondition to be provided in an AbstractArray (of type BoundaryCondition)
function lvsolve(a, L, u, bcs::AbstractArray;
solver_parameters::Dict = Dict("linear_solver" => "default"),
form_compiler_parameters::Dict = Dict("optimize" => true))
bcs_py = [bc.pyobject for bc in bcs]
fenics.solve(a.pyobject == L.pyobject, u.pyobject, bcs = bcs_py,
solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
export lvsolve
#Dict("linear_solver"=>"default")
#Dict("optimize"=>true)
#nlvsolve is the non-linear variational solver
function nlvsolve(F, u; J = nothing,
solver_parameters::Dict = Dict("nonlinear_solver" => "newton"),
form_compiler_parameters::Dict = Dict("optimize" => true))
fenics.solve(F.pyobject == 0, u.pyobject, J = J, solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
function nlvsolve(F, u, bcs = nothing; J = nothing,
solver_parameters::Dict = Dict("nonlinear_solver" => "newton"),
form_compiler_parameters::Dict = Dict("optimize" => true))
fenics.solve(F.pyobject == 0, u.pyobject, J = J, solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
function nlvsolve(F, u, bcs::BoundaryCondition; J = nothing,
solver_parameters::Dict = Dict("nonlinear_solver" => "newton"),
form_compiler_parameters::Dict = Dict("optimize" => true))
fenics.solve(F.pyobject == 0, u.pyobject, bcs = bcs.pyobject, J = J,
solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
#allows BoundaryCondition to be provided in an AbstractArray (of type BoundaryCondition)
function nlvsolve(F, u, bcs::AbstractArray; J = nothing,
solver_parameters::Dict = Dict("nonlinear_solver" => "newton"),
form_compiler_parameters::Dict = Dict("optimize" => true))
bcs_py = [bc.pyobject for bc in bcs]
fenics.solve(F.pyobject == 0, u.pyobject, bcs = bcs, J = J,
solver_parameters = solver_parameters,
form_compiler_parameters = form_compiler_parameters)
end
export nlvsolve
#anlvsolve corresponds to the adapative non-linear solver.
function anlvsolve(F, a, u, bcs, tol, M)
fenics.solve(F.pyobject == a.pyobject, u.pyobject, bcs = bcs.pyobject, tol = tol, M = M)
end
#this function hasnt been tested yet, so isnt exported
function norm(u::FeFunction; normType = "L2", mesh::Union{Nothing, Mesh} = nothing)
if isa(mesh, Nothing)
return fenics.norm(u.pyobject, normType)
else
return fenics.norm(u.pyobject, normType, mesh.pyobject)
end
end
"""
errornorm is the function to calculate the error between our exact and calculated
solution. The norm kwarg defines the norm measure used (by default it is the L2 norm)
"""
errornorm(ans, sol; norm = "L2") = fenics.errornorm(ans.pyobject, sol.pyobject, norm)
export errornorm
File(path::StringOrSymbol) = fenics.File(path) #used to store the solution in various formats
function File(path::StringOrSymbol, object::FeFunction)
vtkfile = File(path)
vtkfile << object.pyobject
end
function File(path::StringOrSymbol, object::FeFunction, time::Number)
vtkfile = File(path)
vtkfile << (object.pyobject, time)
end
function File(path::StringOrSymbol, object::Mesh)
vtkfile = File(path)
vtkfile << object.pyobject
end
function File(path::StringOrSymbol, object::Mesh, time::Number)
vtkfile = File(path)
vtkfile << (object.pyobject, time)
end
export File
XDMFFile(path::StringOrSymbol) = fenics.XDMFFile(path)
export XDMFFile
TimeSeries(path::StringOrSymbol) = fenics.TimeSeries(path)
retrieve(timeseries, placeholder, time) = timeseries.retrieve(placeholder, time)
export TimeSeries, retrieve
write(path, solution, time::Number) = path.write(solution.pyobject, time)
write(path, solution::PyObject, time::Number) = path.write(solution, time)
store(path, solution, time::Number) = path.store(solution.pyobject, time)
store(path, solution::PyObject, time::Number) = path.store(solution, time)
export write, store
array(matrix) = fenicspycall(matrix, :gather_on_zero)
vector(solution::FeFunction) = fenicspycall(solution, :vector) #
interpolate(ex, V::FunctionSpace) = FeFunction(fenics.interpolate(ex.pyobject, V.pyobject))
export vector, interpolate, array
"the following function is used to extract the array from a solution/form"
function get_array(form::Expression)
assembled_form = assemble(form)
return array(assembled_form)
end
function get_array(solution::FeFunction)
generic_vector = vector(solution)
instantiated_vector = fenics.Vector(generic_vector)
return instantiated_vector.gather_on_zero()
end
"""
we return the array from an assembled form
"""
function get_array(assembled_form::Matrix)
return array(assembled_form)
end
export get_array
"""
Return projection of given expression *v* onto the finite element space *V*
*Example of usage*
v = Expression("sin(pi*x[0])")
V = FunctionSpace(mesh, "Lagrange", 1)
Pv = project(v, V)
"""
function project(v::Union{FeFunction, Expression}, V::FunctionSpace)
FeFunction(fenics.project(v.pyobject, V.pyobject))
end
export project