ply_filter_element_oddfract Subroutine

public subroutine ply_filter_element_oddfract(me, nDims, inLen, element_data, filtered_data)

Filter a polynomial representation in elements in one dimension according to its odd mode fraction.

Odd and even modes are weighed with their polynomial degree (so it's a little like actually using the derivative), squared and summed respectively. We then compute a spectral damping order with fraction of the odd mode energy in the total modal energy. The larger the fraction, the larger damping order and the weaker the filtering. uscale is given by (me%max_order-me%min_order+1). odd_fraction is given by odd / (odd+even) and the damping order is then computed by me%min_order + floor(uscale * odd_fraction**me%fract_ex) The fract_exponent provides a mechanism to take larger odd fractions stronger into account than smaller ones.

As we need to perform this operation in all dimensions, it would be good to shift the indices around. When doing this, we can stick to the same implementation for all directions, without the need to put any logic in here to decide on the current direction. In 3D we would end up with this chain: (x,y,z) -> filter_element for Z -> (z,x,y) -> filter_element for Y -> (y,z,x) -> filter_element for X -> (x,y,z) Thus, the logic is that we perform the filter on the last dimension, and cycle the indices in the output.

We can generalize this to arbitrary dimensions. In 2D it would look like this: (x,y) -> filter_element for Y -> (y,x) -> filter_element for X -> (x,y) And in 1D, we just need to perform one transformation: (x) -> filter_element for X -> (x)

We need: nDofs in the direction where the transformation is to be done and the nDofs for all normal directions.

Arguments

TypeIntentOptionalAttributesName
class(ply_filter_element_type), intent(in) :: me

Filter parameters.

integer, intent(in) :: nDims

Number of dimensions of the polynomial data.

integer, intent(in) :: inLen(nDims)

Number degrees of freedom for each direction in element_data.

The first index of element_data needs to have a length equal to the product of all inLen components. The splitting operation will be done in the last dimension.

real(kind=rk), intent(in) :: element_data(:,:)

Polynomial representation in the elements.

The first index are the degrees of freedom in elements, the second index are the elements. In the first index the shape of data has to be in the form (inLen(1), inLen(2), ... , inLen(nDims)). The filtering operation is performed on the last dimension in that data.

real(kind=rk), intent(out) :: filtered_data(:,:)

The filtered polynomial modes.

The ordering is rotated, such, that the filtered dimension becomes the first one, and all others are shifted by one to the right. Thus, the new data has the layout (inLen(nDims), inLen(1), inLen(2), ...)


Called by

proc~~ply_filter_element_oddfract~~CalledByGraph proc~ply_filter_element_oddfract ply_filter_element_oddfract proc~ply_filter_oddfract_1d ply_filter_oddfract_1D proc~ply_filter_oddfract_1d->proc~ply_filter_element_oddfract proc~ply_filter_oddfract_2d ply_filter_oddfract_2D proc~ply_filter_oddfract_2d->proc~ply_filter_element_oddfract proc~ply_filter_oddfract_3d ply_filter_oddfract_3D proc~ply_filter_oddfract_3d->proc~ply_filter_element_oddfract

Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private :: iDir
integer, private :: iParent
integer, private :: parentMode
integer, private :: maxcol
integer, private :: indep
integer, private :: nIndeps
integer, private :: nParents
integer, private :: parentpos
integer, private :: filterpos
integer, private :: upper_scale
real(kind=rk), private :: dof_fract
real(kind=rk), private :: damping
real(kind=rk), private, allocatable:: even(:,:)
real(kind=rk), private, allocatable:: odd(:,:)
real(kind=rk), private :: spenergy
integer, private, allocatable:: damp_ord(:,:)