evalDiff_forElement Subroutine

private recursive subroutine evalDiff_forElement(fun, varSys, elempos, time, tree, nElems, nDofs, res)

Evaluate the function pointers of the dependent variables, and then calculate the difference between these two. ( scalar or vector ) In lua file, first define new variable with varType operation kind as "difference" and provide two dependent variable via input_varname. If input_varname variable is not part of predefined solver variables then add also that variable via variable table for example spacetime function variable. For example: Define a variable called difference, which depend on density and spacetime. one can get an error between simulation results and analytical solution.

\verbatim -- in lua file, one can define as following: variable = {{ name = 'dens_reference', ncomponents = 1, vartype = "st_fun", st_fun = luaFun }, { name = 'dens_difference', ncomponents = 1, vartype = "operation", operation = {kind='difference', input_varname={density, dens_reference}} } ... } tracking = { variable = {'dens_difference'}, folder = 'tracking/', shape = {kind = 'canoND', object = {origin = {3.0,3.1,3.0} } }, format = 'ascii', time = {min = 0, max = tmax, interval = 1}, } \endverbatim

The interface has to comply to the abstract interface tem_varSys_module#tem_varSys_proc_element.

Arguments

TypeIntentOptionalAttributesName
class(tem_varSys_op_type), intent(in) :: fun

Description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varSys

The variable system to obtain the variable from.

integer, intent(in) :: elempos(:)

Position of the TreeID of the element to get the variable for in the global treeID list.

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

type(treelmesh_type), intent(in) :: tree

global treelm mesh info

integer, intent(in) :: nElems

Number of values to obtain for this variable (vectorized access).

integer, intent(in) :: nDofs

Number of degrees of freedom within an element.

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

Resulting values for the requested variable.

Linearized array dimension: (n requested entries) x (nComponents of this variable) x (nDegrees of freedom) Access: (iElem-1)fun%nComponentsnDofs + (iDof-1)*fun%nComponents + iComp


Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private :: iDep
integer, private :: posDepVar
real(kind=rk), private, allocatable:: input_varRes(:,:)
integer, private :: nTotal