Extract component index of any vectorial variable. In lua file, first define new variable with varType operation kind as "extract" and provide name of the variable from which to extract component index via input_varname (it must be single variable) and index to extract via input_varIndex. If input_varname variable is not part of predefined solver variables then add also that variable via variable table.
\verbatim -- in lua file, one can define as following: variable = {{ name = 'vel_y', ncomponents = 1, vartype = "operation", operation = {kind='extract', input_varname={'velocity'}, input_varindex = {2} } }, } \endverbatim
The interface has to comply to the abstract interface tem_varSys_module#tem_varSys_proc_element.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
recursive subroutine extract_forElement( fun, varsys, elempos, time, tree, &
& nElems, nDofs, res )
! ---------------------------------------------------------------------- !
!> 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.
class(tem_varSys_op_type), intent(in) :: fun
!> The variable system to obtain the variable from.
type(tem_varSys_type), intent(in) :: varSys
!> Position of the TreeID of the element to get the variable for in the
!! global treeID list.
integer, intent(in) :: elempos(:)
!> Point in time at which to evaluate the variable.
type(tem_time_type), intent(in) :: time
!> global treelm mesh info
type(treelmesh_type), intent(in) :: tree
!> Number of values to obtain for this variable (vectorized access).
integer, intent(in) :: nElems
!> Number of degrees of freedom within an element.
integer, intent(in) :: nDofs
!> 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%nComponents*nDofs +
!! (iDof-1)*fun%nComponents + iComp
real(kind=rk), intent(out) :: res(:)
! ---------------------------------------------------------------------- !
integer :: iElem, depVar_pos, depVar_nComps, iDof, comp_index, iComp
real(kind=rk), allocatable :: input_varRes(:)
! ---------------------------------------------------------------------- !
! get the position of dependent variable, assuming only one
depVar_pos = fun%input_varPos(1)
! get the nComp of dependent variable
depVar_nComps = varSys%method%val( depVar_pos )%nComponents
! allocate array to save results
allocate(input_varRes( nElems * nDofs * depVar_nComps ))
! derive dependent variable
call varSys%method%val(depVar_pos)%get_element( &
& varSys = varSys, &
& elemPos = elemPos, &
& time = time, &
& tree = tree, &
& nElems = nElems, &
& nDofs = nDofs, &
& res = input_varRes )
! Extract component index from input_varRes
do iElem = 1, nElems
do iDof = 1, nDofs
do iComp = 1, fun%nComponents
comp_index = fun%input_varIndex(iComp)
res((( ielem-1)* fun%ncomponents* ndofs+( idof-1)* fun%ncomponents+icomp) ) &
& = input_varRes( &
& (( ielem-1)* depvar_ncomps* ndofs+( idof-1)* depvar_ncomps+comp_index))
end do
end do
end do
deallocate(input_varRes)
end subroutine extract_forElement