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) && !isdefined(@__MODULE__, :CUDA) 43 error( 44 string( 45 "In order to use user Q-functions with a CUDA backend, the CUDA.jl package ", 46 "must be loaded", 47 ), 48 ) 49 end 50 # COV_EXCL_STOP 51 QFunction(ref, f) 52end 53 54""" 55 create_interior_qfunction(ceed::Ceed, name::AbstractString) 56 57Create a [`QFunction`](@ref) from the Q-function gallery, using the provided name. 58 59# Examples 60 61- Build and apply the 3D mass operator 62``` 63build_mass_qf = create_interior_qfunction(c, "Mass3DBuild") 64apply_mass_qf = create_interior_qfunction(c, "MassApply") 65``` 66- Build and apply the 3D Poisson operator 67``` 68build_poi_qf = create_interior_qfunction(c, "Poisson3DBuild") 69apply_poi_qf = create_interior_qfunction(c, "Poisson3DApply") 70``` 71""" 72function create_interior_qfunction(c::Ceed, name::AbstractString) 73 ref = Ref{C.CeedQFunction}() 74 C.CeedQFunctionCreateInteriorByName(c.ref[], name, ref) 75 QFunction(ref) 76end 77 78""" 79 create_identity_qfunction(c::Ceed, size, inmode::EvalMode, outmode::EvalMode) 80 81Create an identity [`QFunction`](@ref). Inputs are written into outputs in the order given. 82This is useful for [`Operators`](@ref Operator) that can be represented with only the action 83of a [`ElemRestriction`](@ref) and [`Basis`](@ref), such as restriction and prolongation 84operators for p-multigrid. Backends may optimize `CeedOperators` with this Q-function to 85avoid the copy of input data to output fields by using the same memory location for both. 86""" 87function create_identity_qfunction(c::Ceed, size, inmode::EvalMode, outmode::EvalMode) 88 ref = Ref{C.CeedQFunction}() 89 C.CeedQFunctionCreateIdentity(c[], size, inmode, outmode, ref) 90 QFunction(ref) 91end 92 93function add_input!(qf::AbstractQFunction, name::AbstractString, size, emode) 94 C.CeedQFunctionAddInput(qf[], name, size, emode) 95end 96 97function add_output!(qf::AbstractQFunction, name::AbstractString, size, emode) 98 C.CeedQFunctionAddOutput(qf[], name, size, emode) 99end 100 101""" 102 set_context!(qf::QFunction, ctx::Context) 103 104Associate a [`Context`](@ref) object `ctx` with the given Q-function `qf`. 105""" 106function set_context!(qf::QFunction, ctx) 107 # Preserve the context data from the GC by storing a reference 108 qf.ctx = ctx 109 C.CeedQFunctionSetContext(qf[], ctx[]) 110end 111 112function get_field_sizes(qf::AbstractQFunction) 113 ninputs = Ref{CeedInt}() 114 noutputs = Ref{CeedInt}() 115 inputs = Ref{Ptr{C.CeedQFunctionField}}() 116 outputs = Ref{Ptr{C.CeedQFunctionField}}() 117 C.CeedQFunctionGetFields(qf[], ninputs, inputs, noutputs, outputs) 118 119 input_sizes = zeros(CeedInt, ninputs[]) 120 output_sizes = zeros(CeedInt, noutputs[]) 121 122 for i = 1:ninputs[] 123 field = unsafe_load(inputs[], i) 124 C.CeedQFunctionFieldGetSize(field, pointer(input_sizes, i)) 125 end 126 127 for i = 1:noutputs[] 128 field = unsafe_load(outputs[], i) 129 C.CeedQFunctionFieldGetSize(field, pointer(output_sizes, i)) 130 end 131 132 input_sizes, output_sizes 133end 134 135""" 136 apply!(qf::QFunction, Q, vin, vout) 137 138Apply the action of a [`QFunction`](@ref) to an array of input vectors, and store the result 139in an array of output vectors. 140""" 141function apply!(qf::QFunction, Q, vin, vout) 142 vins = map(x -> x[], vin) 143 vouts = map(x -> x[], vout) 144 GC.@preserve vin vout begin 145 C.CeedQFunctionApply(qf[], Q, pointer(vins), pointer(vouts)) 146 end 147end 148