$$ \gdef\pd#1#2{\frac{\partial #1}{\partial #2}} \gdef\d#1{\, \mathrm{d}#1} \gdef\dx{\d{x}} \gdef\tr#1{\operatorname{tr} (#1)} $$ $$ \gdef\norm#1{\left \lVert #1 \right\rVert} \gdef\seminorm#1{| #1 |} $$ $$ \gdef\vec#1{\mathbf{\boldsymbol{#1}}} \gdef\dvec#1{\bar{\vec #1}} $$
pub struct Quad9d2Element<T>where
    T: Scalar,{ /* private fields */ }
Expand description

A finite element representing quadratic basis functions on a quad, in two dimensions.

Implementations§

source§

impl<T> Quad9d2Element<T>where T: Scalar,

source

pub fn from_vertices(vertices: [Point2<T>; 9]) -> Self

source

pub fn vertices(&self) -> &[Point2<T>; 9]

source§

impl<T> Quad9d2Element<T>where T: Real,

source

pub fn reference() -> Self

Trait Implementations§

source§

impl<T> CanonicalMassQuadrature for Quad9d2Element<T>where T: Real,

source§

impl<T> CanonicalStiffnessQuadrature for Quad9d2Element<T>where T: Real,

source§

impl<T> Clone for Quad9d2Element<T>where T: Scalar + Clone,

source§

fn clone(&self) -> Quad9d2Element<T>

Returns a copy of the value. Read more
1.0.0 · source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
source§

impl<T> Debug for Quad9d2Element<T>where T: Scalar + Debug,

source§

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

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

impl<T> FiniteElement<T> for Quad9d2Element<T>where T: Real,

§

type GeometryDim = Const<2>

source§

fn reference_jacobian(&self, xi: &Point2<T>) -> Matrix2<T>

Compute the Jacobian of the transformation from the reference element to the given element at the given reference coordinates.
source§

fn map_reference_coords(&self, xi: &Point2<T>) -> Point2<T>

Maps reference coordinates to physical coordinates in the element.
source§

fn diameter(&self) -> T

The diameter of the finite element. Read more
source§

impl<T> FixedNodesReferenceFiniteElement<T> for Quad9d2Element<T>where T: Real,

§

type ReferenceDim = Const<2>

§

type NodalDim = Const<9>

source§

fn evaluate_basis(&self, xi: &Point2<T>) -> OMatrix<T, U1, U9>

Evaluates each basis function at the given reference coordinates. The result is given in a row vector where each entry is the value of the corresponding basis function.
source§

fn gradients(&self, xi: &Point2<T>) -> OMatrix<T, U2, U9>

Given nodal weights, construct a matrix whose columns are the gradients of each shape function in the element.
source§

impl<'a, T> From<&'a Quad4d2Element<T>> for Quad9d2Element<T>where T: Real,

source§

fn from(quad4: &'a Quad4d2Element<T>) -> Self

Converts to this type from the input type.
source§

impl<'a, T> From<Quad4d2Element<T>> for Quad9d2Element<T>where T: Real,

source§

fn from(quad4: Quad4d2Element<T>) -> Self

Converts to this type from the input type.
source§

impl<T> PartialEq<Quad9d2Element<T>> for Quad9d2Element<T>where T: Scalar + PartialEq,

source§

fn eq(&self, other: &Quad9d2Element<T>) -> bool

This method tests for self and other values to be equal, and is used by ==.
1.0.0 · source§

fn ne(&self, other: &Rhs) -> bool

This method tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
source§

impl<T> ReferenceFiniteElement<T> for Quad9d2Element<T>where T: Scalar, Quad9d2Element<T>: FixedNodesReferenceFiniteElement<T>, DefaultAllocator: BiDimAllocator<T, <Quad9d2Element<T> as FixedNodesReferenceFiniteElement<T>>::NodalDim, <Quad9d2Element<T> as FixedNodesReferenceFiniteElement<T>>::ReferenceDim>,

§

type ReferenceDim = <Quad9d2Element<T> as FixedNodesReferenceFiniteElement<T>>::ReferenceDim

source§

fn num_nodes(&self) -> usize

Returns the number of nodes in the element.
source§

fn populate_basis( &self, result: &mut [T], reference_coords: &OPoint<T, Self::ReferenceDim> )

Evaluates each basis function at the given reference coordinates. The result is given in a row vector where each entry is the value of the corresponding basis function. Read more
source§

fn populate_basis_gradients( &self, result: MatrixViewMut<'_, T, Self::ReferenceDim, Dyn>, reference_coords: &OPoint<T, Self::ReferenceDim> )

Given nodal weights, construct a matrix whose columns are the gradients of each shape function in the element.
source§

impl<T> TryFrom<Quad9d2Element<T>> for ConvexPolygon<T>where T: Real,

§

type Error = ConcavePolygonError

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

fn try_from(value: Quad9d2Element<T>) -> Result<Self, Self::Error>

Performs the conversion.
source§

impl<T> Copy for Quad9d2Element<T>where T: Scalar + Copy,

source§

impl<T> Eq for Quad9d2Element<T>where T: Scalar + Eq,

source§

impl<T> StructuralEq for Quad9d2Element<T>where T: Scalar,

source§

impl<T> StructuralPartialEq for Quad9d2Element<T>where T: Scalar,

Auto Trait Implementations§

§

impl<T> RefUnwindSafe for Quad9d2Element<T>where T: RefUnwindSafe,

§

impl<T> Send for Quad9d2Element<T>where T: Send,

§

impl<T> Sync for Quad9d2Element<T>where T: Sync,

§

impl<T> Unpin for Quad9d2Element<T>where T: Unpin,

§

impl<T> UnwindSafe for Quad9d2Element<T>where T: UnwindSafe,

Blanket Implementations§

source§

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

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

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

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for Twhere 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 Twhere 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.

§

impl<T> Pointable for T

§

const ALIGN: usize = _

The alignment of pointer.
§

type Init = T

The type for initializers.
§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
source§

impl<T> Same<T> for T

§

type Output = T

Should always be Self
§

impl<SS, SP> SupersetOf<SS> for SPwhere SS: SubsetOf<SP>,

§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
source§

impl<T> ToOwned for Twhere T: Clone,

§

type Owned = T

The resulting type after obtaining ownership.
source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
source§

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

§

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 Twhere U: TryFrom<T>,

§

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.
source§

impl<T> Scalar for Twhere T: 'static + Clone + PartialEq<T> + Debug,