A routine to evaluate chunk of elements for given list of variables
If subTree is present, it will use map2Global from subTree else map2Global is created for current chunk of global mesh
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_varSys_type), | intent(in) | :: | varsys |
Variable system describing available data. |
||
integer, | intent(in) | :: | varPos(:) |
Position of variables to evaluate in varSys |
||
integer, | intent(in) | :: | elemPos(:) |
Position of treeID of the element to get the variable for. |
||
type(tem_time_type), | intent(in) | :: | time |
Time information for the current data. |
||
type(treelmesh_type), | intent(in) | :: | tree |
Mesh definition of the input data. |
||
integer, | intent(in) | :: | nElems |
number of elements to evaluate |
||
integer, | intent(in) | :: | nDofs |
Number of degrees of freedom. |
||
real(kind=rk), | intent(out) | :: | res(:) |
Output data size: io_buffer_size |
subroutine tem_get_element_chunk(varSys, varPos, elemPos, time, tree, nElems,&
& nDofs, res)
! --------------------------------------------------------------------------!
!> Variable system describing available data.
type(tem_varsys_type), intent(in) :: varsys
!> Position of variables to evaluate in varSys
integer, intent(in) :: varPos(:)
!> Mesh definition of the input data.
type(treelmesh_type), intent(in) :: tree
!> Time information for the current data.
type(tem_time_type), intent(in) :: time
!> Number of degrees of freedom.
integer, intent(in) :: nDofs
!> Position of treeID of the element to get the variable for.
integer, intent(in) :: elemPos(:)
!> number of elements to evaluate
integer, intent(in) :: nElems
!> Output data
!! size: io_buffer_size
real(kind=rk), intent(out) :: res(:)
! --------------------------------------------------------------------------!
integer :: maxComponents, nScalars, nComponents, elemsize, compOff
integer :: iElem, iVar, iDof, nVars, res_size
integer :: e_start, t_start, d_start
real(kind=rk), allocatable :: tmpdat(:)
! --------------------------------------------------------------------------!
! number of variables to evaluate
nVars = size(varPos)
! Number of scalars in current output
nScalars = sum(varSys%method%val(varPos(:))%nComponents)
! Size of a single element
elemsize = nScalars*nDofs
! Need to obtain the data variable for variable, and store it in an
! intermediate array, because all components should be put together in the
! res array.
! The temporary array therefore needs to be sufficiently large to store the
! maximal number of components.
maxComponents = maxval(varSys%method%val(varPos(:))%nComponents)
! Using a temporary array to store the variables and transfer them to res
! in the correct ordering afterwards.
allocate(tmpdat(nElems*maxComponents*nDofs))
compOff = 0
vars: do iVar=1, nVars
! get the number of components for variable iVar
nComponents = varSys%method%val(varPos(iVar))%nComponents
! get the size of the needed part of the res array
res_size = nElems * nDofs * nComponents
! derive the quantities for all the elements in the current chunk
call varSys%method%val(varpos(iVar))%get_element( &
& varSys = varSys, &
& elemPos = elemPos, &
& time = time, &
& tree = tree, &
& nElems = nElems, &
& nDofs = nDofs, &
& res = tmpdat(:res_size) )
! copy the information to the right positions in the result array
! res contains results for all variables,
! tmpdat is only for one variable
do iElem=1,nElems
e_start = (iElem-1)*elemsize
t_start = (iElem-1)*nComponents*nDofs
do iDof=1,nDofs
d_start = (iDof-1)*nScalars + compOff
res( (e_start+d_start+1) : (e_start+d_start+nComponents) ) &
& = tmpdat( t_start + (iDof-1)*nComponents + 1 &
& :t_start + iDof*nComponents )
end do
end do
! Increase the component offset for the next variables.
compOff = compOff + nComponents
end do vars
end subroutine tem_get_element_chunk