Using Method of Manufactured solutions and solve the Helmholtz equation.

Problem taken from http://www.dealii.org/8.4.1/doxygen/deal.II/step_7.html


In [ ]:
using JuAFEM
using KrylovMethods
const  = Tensors.gradient
const Δ = Tensors.hessian;

In [ ]:
grid = generate_grid(Quadrilateral, (150, 150))
addnodeset!(grid, "dirichlet_boundary", x -> x[1]  1 ||  x[2]  1);

In [ ]:
dim = 2
ip = Lagrange{dim, RefCube, 1}()
qr = QuadratureRule{dim, RefCube}(2)
qr_face = QuadratureRule{dim-1, RefCube}(2)
cellvalues = CellScalarValues(qr, ip);
facevalues = FaceScalarValues(qr_face, ip);

In [ ]:
dh = DofHandler(grid)
push!(dh, :u, 1)
close!(dh)


Out[ ]:
DofHandler
  Fields:
    u dim: 1
  Total dofs: 22801
  Dofs per cell: 4

In [ ]:
function u_ana{T}(x::Vec{2, T}) 
    xs = (Vec{2}((-0.5,  0.5)),
          Vec{2}((-0.5, -0.5)),
          Vec{2}(( 0.5,  -0.5)))
    σ = 1/8
    s = zero(eltype(x))
    for i in 1:3
        s += exp(- norm(x - xs[i])^2 / σ^2)
    end
    return max(1e-15 * one(T), s) # Denormals, be gone
end;

In [ ]:
dbc = DirichletBoundaryConditions(dh)
add!(dbc, :u, getnodeset(grid, "dirichlet_boundary"), (x,t) -> u_ana(x))
close!(dbc)
update!(dbc, 0.0)

In [ ]:
K = create_sparsity_pattern(dh);

In [ ]:
function doassemble{dim}(cellvalues::CellScalarValues{dim}, facevalues::FaceScalarValues{dim},
                         K::SparseMatrixCSC, dh::DofHandler)
    b = 1.0
    f = zeros(ndofs(dh))
    assembler = start_assemble(K, f)
    
    n_basefuncs = getnbasefunctions(cellvalues)
    global_dofs = zeros(Int, ndofs_per_cell(dh))

    fe = zeros(n_basefuncs) # Local force vector
    Ke = zeros(n_basefuncs, n_basefuncs) # Local stiffness mastrix

    @inbounds for (cellcount, cell) in enumerate(CellIterator(dh))
        fill!(Ke, 0)
        fill!(fe, 0)
        coords = getcoordinates(cell)

        reinit!(cellvalues, cell)
        for q_point in 1:getnquadpoints(cellvalues)
             = getdetJdV(cellvalues, q_point)
            
            coords_qp = spatial_coordinate(cellvalues, q_point, coords)
            f_true = trace(hessian(u_ana, coords_qp)) + u_ana(coords_qp)
            for i in 1:n_basefuncs
                δu = shape_value(cellvalues, q_point, i)
                ∇δu = shape_gradient(cellvalues, q_point, i)
                fe[i] += (δu * f_true) * 
                for j in 1:n_basefuncs
                    u = shape_value(cellvalues, q_point, j)
                    ∇u = shape_gradient(cellvalues, q_point, j)
                    Ke[i, j] += (∇δu  ∇u + δu * u) * 
                end
            end
        end
        
        for face in 1:nfaces(cell)
            if onboundary(cell, face) && 
                   ((cellcount, face)  getfaceset(grid, "left") || 
                    (cellcount, face)  getfaceset(grid, "bottom"))
                reinit!(facevalues, cell, face)
                for q_point in 1:getnquadpoints(facevalues)
                    coords_qp = spatial_coordinate(facevalues, q_point, coords)
                    n = getnormal(facevalues, q_point)
                    g = gradient(u_ana, coords_qp)  n
                     = getdetJdV(facevalues, q_point)
                    for i in 1:n_basefuncs
                        δu = shape_value(facevalues, q_point, i)
                        fe[i] += (δu * g) * 
                    end
                end
            end
        end
   
        celldofs!(global_dofs, cell)
        assemble!(assembler, global_dofs, fe, Ke)
    end
    return K, f
end;

In [ ]:
K, f = doassemble(cellvalues, facevalues, K, dh);
apply!(K, f, dbc)
u = cholfact(Symmetric(K)) \ f;

In [ ]:
vtkfile = vtk_grid("helmholtz", dh, u)
vtk_save(vtkfile)


Out[ ]:
1-element Array{String,1}:
 "helmholtz.vtu"

In [ ]:
Base.Test.@test maximum(u)  0.05637592090022005
println("Helmholtz successful")


Helmholtz successful