Struct QFunction

Source
pub struct QFunction<'a> { /* private fields */ }

Implementations§

Source§

impl<'a> QFunction<'a>

Source

pub fn create( ceed: &Ceed, vlength: usize, user_f: Box<QFunctionUserClosure>, ) -> Result<Self>

Source

pub fn apply(&self, Q: usize, u: &[Vector<'_>], v: &[Vector<'_>]) -> Result<i32>

Apply the action of a QFunction

  • Q - The number of quadrature points
  • input - Array of input Vectors
  • output - Array of output Vectors
let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let qf = ceed
    .q_function_interior(1, Box::new(user_f))?
    .input("u", 1, EvalMode::Interp)?
    .input("weights", 1, EvalMode::Weight)?
    .output("v", 1, EvalMode::Interp)?;

const Q: usize = 8;
let mut w = [0.; Q];
let mut u = [0.; Q];
let mut v = [0.; Q];

for i in 0..Q {
    let x = 2. * (i as Scalar) / ((Q as Scalar) - 1.) - 1.;
    u[i] = 2. + 3. * x + 5. * x * x;
    w[i] = 1. - x * x;
    v[i] = u[i] * w[i];
}

let uu = ceed.vector_from_slice(&u)?;
let ww = ceed.vector_from_slice(&w)?;
let mut vv = ceed.vector(Q)?;
vv.set_value(0.0);
{
    let input = vec![uu, ww];
    let mut output = vec![vv];
    qf.apply(Q, &input, &output)?;
    vv = output.remove(0);
}

vv.view()?
    .iter()
    .zip(v.iter())
    .for_each(|(computed, actual)| {
        assert_eq!(
            *computed, *actual,
            "Incorrect value in QFunction application"
        );
    });
Source

pub fn input( self, fieldname: &str, size: usize, emode: EvalMode, ) -> Result<Self>

Add a QFunction input

  • fieldname - Name of QFunction field
  • size - Size of QFunction field, (ncomp * dim) for Grad or (ncomp * 1) for None, Interp, and Weight
  • emode - EvalMode::None to use values directly, EvalMode::Interp to use interpolated values, EvalMode::Grad to use gradients, EvalMode::Weight to use quadrature weights
let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let mut qf = ceed.q_function_interior(1, Box::new(user_f))?;

qf = qf.input("u", 1, EvalMode::Interp)?;
qf = qf.input("weights", 1, EvalMode::Weight)?;
Source

pub fn output( self, fieldname: &str, size: usize, emode: EvalMode, ) -> Result<Self>

Add a QFunction output

  • fieldname - Name of QFunction field
  • size - Size of QFunction field, (ncomp * dim) for Grad or (ncomp * 1) for None and Interp
  • emode - EvalMode::None to use values directly, EvalMode::Interp to use interpolated values, EvalMode::Grad to use gradients
let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let mut qf = ceed.q_function_interior(1, Box::new(user_f))?;

qf.output("v", 1, EvalMode::Interp)?;
Source

pub fn inputs(&self) -> Result<&[QFunctionField<'_>]>

Get a slice of QFunction inputs

let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let mut qf = ceed
    .q_function_interior(1, Box::new(user_f))?
    .input("u", 1, EvalMode::Interp)?
    .input("weights", 1, EvalMode::Weight)?;

let inputs = qf.inputs()?;

assert_eq!(inputs.len(), 2, "Incorrect inputs array");
Source

pub fn outputs(&self) -> Result<&[QFunctionField<'_>]>

Get a slice of QFunction outputs

let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let mut qf = ceed
    .q_function_interior(1, Box::new(user_f))?
    .output("v", 1, EvalMode::Interp)?;

let outputs = qf.outputs()?;

assert_eq!(outputs.len(), 1, "Incorrect outputs array");

Trait Implementations§

Source§

impl<'a> Display for QFunction<'a>

View a QFunction

let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let qf = ceed
    .q_function_interior(1, Box::new(user_f))?
    .input("u", 1, EvalMode::Interp)?
    .input("weights", 1, EvalMode::Weight)?
    .output("v", 1, EvalMode::Interp)?;

println!("{}", qf);
Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<'a> Drop for QFunction<'a>

Source§

fn drop(&mut self)

Executes the destructor for this type. Read more
Source§

impl<'a> From<&'a QFunction<'_>> for QFunctionOpt<'a>

Construct a QFunctionOpt reference from a QFunction reference

Source§

fn from(qfunc: &'a QFunction<'_>) -> Self

Converts to this type from the input type.

Auto Trait Implementations§

§

impl<'a> Freeze for QFunction<'a>

§

impl<'a> !RefUnwindSafe for QFunction<'a>

§

impl<'a> !Send for QFunction<'a>

§

impl<'a> !Sync for QFunction<'a>

§

impl<'a> Unpin for QFunction<'a>

§

impl<'a> !UnwindSafe for QFunction<'a>

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> ToString for T
where T: Display + ?Sized,

Source§

fn to_string(&self) -> String

Converts the given value to a String. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.