Call the spacetime function for vector variable, which are stored in method_data which intern points to tem_st_fun_listElem_type. This spacetime function can be used as a analytical solution for reference. NOTE: nDofs is not used by this routine. So, to evaluate spacetime function for nDofs in an element, use tem_varSys_proc_point interface.
\verbatim -- in lua file, one can define it as following: variable = {{ name = 'reference', ncomponents = 3, vartype = st_fun, st_fun = luaFun }, ... } \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 evaluate_add_spacetime_vectorByTreeID( 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 :: iStFun, iElem, pos, oldPos
type(tem_st_fun_listElem_type), pointer :: fPtr
! element wise output
real(kind=rk), allocatable :: st_res(:,:)
real(kind=rk) :: st_res_elem(1, fun%nComponents)
integer(kind=long_k) :: treeID(1)
! -------------------------------------------------------------------------!
!write(*,*) 'derive variable label: ', trim(varSys%varname%val(fun%myPos))
call C_F_POINTER( fun%method_Data, fPtr )
res = 0.0_rk
! assuming nDofs = 1
if (nDofs > 1) then
write(logUnit(1),*) 'Spacetime function does not support nDofs>1'
call tem_abort()
end if
allocate(st_res(nElems, fun%nComponents))
st_res = 0.0_rk
do iStFun = 1, fPtr%nVals
if (fPtr%val(iStFun)%subTree%useGlobalMesh ) then
!write(*,*) 'global mesh'
st_res = st_res + &
& tem_spacetime_for( me = fPtr%val(iStFun), &
& treeIDs = tree%treeID(elemPos(:)), &
& tree = tree, &
& time = time, &
& nComp = fun%nComponents, &
& n = nElems )
else
!write(*,*) 'subTree Mesh'
oldPos = 1
do iElem = 1, nElems
! position of element matches with any of map2global then
! this element is part of this subTree
pos = tem_PositionInSorted( &
& me = fPtr%val(iStFun)%subTree%map2global, &
& val = elemPos(iElem), &
& lower = oldPos )
! pos>0 if elemPos(iElem) is part of subTree
if (pos > 0) then
oldPos = pos
treeID = tree%treeID(elemPos(iElem))
st_res_elem = &
& tem_spacetime_for( me = fPtr%val(iStFun), &
& treeIDs = treeID, &
& tree = tree, &
& time = time, &
& nComp = fun%nComponents, &
& n = 1 )
st_res(iElem,:) = st_res(iElem,:) + st_res_elem(1,:)
end if
end do
end if ! global mesh
end do ! iStFun
do iElem = 1, nElems
res( (iElem-1)*fun%nComponents + 1 : iElem*fun%nComponents ) &
& = st_res(iElem, :)
end do
end subroutine evaluate_add_spacetime_vectorByTreeID