xref: /libCEED/julia/LibCEED.jl/src/QFunction.jl (revision 44554ea01e90fce366fc2a203c44be15754a38d6)
1abstract type AbstractQFunction end
2
3struct QFunctionNone <: AbstractQFunction end
4Base.getindex(::QFunctionNone) = C.CEED_QFUNCTION_NONE[]
5
6"""
7    QFunction
8
9A libCEED `CeedQFunction` object, typically created using the [`@interior_qf`](@ref) macro.
10
11A `QFunction` can also be created from the "Q-function gallery" using
12[`create_interior_qfunction`](@ref). The identity Q-function can be created using
13[`create_identity_qfunction`](@ref).
14"""
15mutable struct QFunction <: AbstractQFunction
16    ref::RefValue{C.CeedQFunction}
17    user_qf::Union{Nothing,UserQFunction}
18    ctx::Union{Nothing,Context}
19    function QFunction(ref, user_qf)
20        obj = new(ref, user_qf, nothing)
21        finalizer(obj) do x
22            # ccall(:jl_safe_printf, Cvoid, (Cstring, Cstring), "Finalizing %s.\n", repr(x))
23            destroy(x)
24        end
25        return obj
26    end
27end
28QFunction(ref::Ref{C.CeedQFunction}) = QFunction(ref, nothing)
29destroy(qf::QFunction) = C.CeedQFunctionDestroy(qf.ref) # COV_EXCL_LINE
30Base.getindex(qf::QFunction) = qf.ref[]
31Base.show(io::IO, ::MIME"text/plain", qf::QFunction) =
32    ceed_show(io, qf, C.CeedQFunctionView)
33
34function create_interior_qfunction(c::Ceed, f::UserQFunction; vlength=1)
35    ref = Ref{C.CeedQFunction}()
36    # Use empty string as source location to indicate to libCEED that there is
37    # no C source for this Q-function
38    C.CeedQFunctionCreateInterior(c[], vlength, f.fptr, "", ref)
39    # COV_EXCL_START
40    if !isnothing(f.cuf)
41        C.CeedQFunctionSetCUDAUserFunction(ref[], f.cuf)
42    elseif iscuda(c) && !cuda_is_loaded
43        error(string(
44            "In order to use user Q-functions with a CUDA backend, the CUDA.jl package ",
45            "must be loaded",
46        ))
47    end
48    # COV_EXCL_STOP
49    QFunction(ref, f)
50end
51
52"""
53    create_interior_qfunction(ceed::Ceed, name::AbstractString)
54
55Create a [`QFunction`](@ref) from the Q-function gallery, using the provided name.
56
57# Examples
58
59- Build and apply the 3D mass operator
60```
61build_mass_qf = create_interior_qfunction(c, "Mass3DBuild")
62apply_mass_qf = create_interior_qfunction(c, "MassApply")
63```
64- Build and apply the 3D Poisson operator
65```
66build_poi_qf = create_interior_qfunction(c, "Poisson3DBuild")
67apply_poi_qf = create_interior_qfunction(c, "Poisson3DApply")
68```
69"""
70function create_interior_qfunction(c::Ceed, name::AbstractString)
71    ref = Ref{C.CeedQFunction}()
72    C.CeedQFunctionCreateInteriorByName(c.ref[], name, ref)
73    QFunction(ref)
74end
75
76"""
77    create_identity_qfunction(c::Ceed, size, inmode::EvalMode, outmode::EvalMode)
78
79Create an identity [`QFunction`](@ref). Inputs are written into outputs in the order given.
80This is useful for [`Operators`](@ref Operator) that can be represented with only the action
81of a [`ElemRestriction`](@ref) and [`Basis`](@ref), such as restriction and prolongation
82operators for p-multigrid. Backends may optimize `CeedOperators` with this Q-function to
83avoid the copy of input data to output fields by using the same memory location for both.
84"""
85function create_identity_qfunction(c::Ceed, size, inmode::EvalMode, outmode::EvalMode)
86    ref = Ref{C.CeedQFunction}()
87    C.CeedQFunctionCreateIdentity(c[], size, inmode, outmode, ref)
88    QFunction(ref)
89end
90
91function add_input!(qf::AbstractQFunction, name::AbstractString, size, emode)
92    C.CeedQFunctionAddInput(qf[], name, size, emode)
93end
94
95function add_output!(qf::AbstractQFunction, name::AbstractString, size, emode)
96    C.CeedQFunctionAddOutput(qf[], name, size, emode)
97end
98
99"""
100    set_context!(qf::QFunction, ctx::Context)
101
102Associate a [`Context`](@ref) object `ctx` with the given Q-function `qf`.
103"""
104function set_context!(qf::QFunction, ctx)
105    # Preserve the context data from the GC by storing a reference
106    qf.ctx = ctx
107    C.CeedQFunctionSetContext(qf[], ctx[])
108end
109
110function get_field_sizes(qf::AbstractQFunction)
111    ninputs = Ref{CeedInt}()
112    noutputs = Ref{CeedInt}()
113
114    C.CeedQFunctionGetNumArgs(qf[], ninputs, noutputs)
115
116    inputs = Ref{Ptr{C.CeedQFunctionField}}()
117    outputs = Ref{Ptr{C.CeedQFunctionField}}()
118    C.CeedQFunctionGetFields(qf[], inputs, outputs)
119
120    input_sizes = zeros(CeedInt, ninputs[])
121    output_sizes = zeros(CeedInt, noutputs[])
122
123    for i = 1:ninputs[]
124        field = unsafe_load(inputs[], i)
125        C.CeedQFunctionFieldGetSize(field, pointer(input_sizes, i))
126    end
127
128    for i = 1:noutputs[]
129        field = unsafe_load(outputs[], i)
130        C.CeedQFunctionFieldGetSize(field, pointer(output_sizes, i))
131    end
132
133    input_sizes, output_sizes
134end
135
136"""
137    apply!(qf::QFunction, Q, vin, vout)
138
139Apply the action of a [`QFunction`](@ref) to an array of input vectors, and store the result
140in an array of output vectors.
141"""
142function apply!(qf::QFunction, Q, vin, vout)
143    vins = map(x -> x[], vin)
144    vouts = map(x -> x[], vout)
145    GC.@preserve vin vout begin
146        C.CeedQFunctionApply(qf[], Q, pointer(vins), pointer(vouts))
147    end
148end
149