-
Notifications
You must be signed in to change notification settings - Fork 112
Expand file tree
/
Copy pathConstantFEOperatorsTests.jl
More file actions
68 lines (51 loc) · 1.09 KB
/
ConstantFEOperatorsTests.jl
File metadata and controls
68 lines (51 loc) · 1.09 KB
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
module ConstantFEOperatorsTests
using Gridap
using ForwardDiff
using LinearAlgebra
using Test
using Gridap.FESpaces: get_algebraic_operator
using LineSearches: BackTracking
θ = 1.0
u(x,t) = (1.0-x[1])*x[1]*(1.0-x[2])*x[2]
u(t::Real) = x -> u(x,t)
∂tu = ∂t(u)
f(t) = x -> ∂t(u)(x,t)-Δ(u(t))(x)
domain = (0,1,0,1)
partition = (2,2)
model = CartesianDiscreteModel(domain,partition)
order = 2
reffe = ReferenceFE(lagrangian,Float64,order)
V0 = FESpace(
model,
reffe,
conformity=:H1,
dirichlet_tags="boundary"
)
U = TransientTrialFESpace(V0,u)
Ω = Triangulation(model)
degree = 2*order
dΩ = Measure(Ω,degree)
#
a(u,v) = ∫(∇(v)⋅∇(u))dΩ
b(v) = ∫(v*f(0.0))dΩ
m(ut,v) = ∫(ut*v)dΩ
op = TransientConstantFEOperator(m,a,b,U,V0)
t0 = 0.0
tF = 1.0
dt = 0.1
U0 = U(0.0)
uh0 = interpolate_everywhere(u(0.0),U0)
ls = LUSolver()
ode_solver = ThetaMethod(ls,dt,θ)
sol_t = solve(ode_solver,op,uh0,t0,tF)
l2(w) = w*w
tol = 1.0e-6
_t_n = t0
for (uh_tn, tn) in sol_t
global _t_n
_t_n += dt
e = u(tn) - uh_tn
el2 = sqrt(sum( ∫(l2(e))dΩ ))
@test el2 < tol
end
end #module