tem_reduction_transient_module Module

Transient reduction enables the definition of quantities that are reduced over a given interval of iterations. It is mainly meant to provide temporally averaged values, but minimum, maximum and sum are also available. For this a double buffer is used, where the current inverval is first filled with the intermediate computation, while the resulting reduced value resides in the other storage place. This double buffer is provided for all elements in the domain. When the value of the transient reduction is requested (via the getElement function) the previously completed interval is returned. This means the data is only useful after the first completed interval and remains constant until the next interval is completed.

A variable providing the temporal average of a variable called pressure would for example be defined as follows:

  variable = {
    {
       name = 'press_timeavg',
       ncomponents = 1,
       vartype = "operation",
       operation = {
         kind='reduction_transient',
         input_varname={'pressure'},
         reduction_transient = {
           kind    = 'average',
           nrecord = 1000
         }
       }
    }
  }

Note that each transient variable definition increases the memory consumption by two reals for every degree of freedom of the input variable in the domain.



Derived Types

Contains transient reduction info loaded from variable table for reduction_transient operation kind

Components

Type Visibility Attributes Name Initial
character(len=labelLen), public :: reduceType

Type of time reduction operation

integer, public :: nRecord = 0

number of interations to record

type, public ::  tem_reduction_transient_type

all data needed for a transient reduction, operation to perform and necessary data from previous timesteps

Components

Type Visibility Attributes Name Initial
type(tem_reduction_transient_config_type), public :: config

reduction info loaded from config file

integer, public :: nTimes = 0

Number of "recorded" previous iterations

integer, public :: nEntries = 0

Number of values to store in the double buffer val

integer, public :: nComponents = 0

number of components

integer, public :: nDofs

Number of degrees of freedom

integer, public :: curr

Index of the storage for the currently filling (incomplete) position in val

integer, public :: last

Index of the storage location with the previously completed interval holding the reduced value over that interval.

real(kind=rk), public, allocatable :: val(:,:)

Double buffer to store data from previous timesteps size (nComponentsnDofstree%nElems,2) 2nd index is used to maintain last valid reduced value when nRecord is reached. It will be swapped to avoid copy operations


Subroutines

public subroutine tem_reduction_transient_load(me, conf, parent)

Read time reduction info from reduction_transient operation variable

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_transient_config_type), intent(out) :: me

time reduction

type(flu_State), intent(inout) :: conf

handle for lua file

integer, intent(in) :: parent

operation table handle

public subroutine tem_reduction_transient_init(me, nElems, nComponents, nDofs)

Initialize time reduction.

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_transient_type), intent(inout) :: me

current time reduction

integer, intent(in) :: nElems

Number of elements in tree

integer, intent(in) :: nComponents

nComponents of operation variable

integer, intent(in) :: nDofs

Number of degrees of freedom

public subroutine tem_reduction_transient_reset(me)

Reset time reduction.

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_transient_type), intent(inout) :: me

current time reduction

public subroutine tem_reduction_transient_update(me, res)

Apply time reduction on given res

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_transient_type), intent(inout) :: me

Time reduction data to update

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

Current values of the variable to reduce.

public subroutine tem_reduction_transient_getElement(me, elemPos, nElems, nDofs, res)

This routine returns the time reduced value for given elemPos

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_transient_type), intent(in) :: me

current time reduction

integer, intent(in) :: elemPos(:)

Position of elements in global tree is same as me%val

integer, intent(in) :: nElems

Number of elements to return

integer, intent(in) :: nDofs

Number of degrees of freedom to return

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

Result array