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