Call the spacetime function for scalar 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 = 1, 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_scalarByTreeID( 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) :: st_res(1) integer(kind=long_k) :: treeID(1) ! -------------------------------------------------------------------------! !call C_F_POINTER( fun%method_Data, fPtr ) ! avoid compiler warning about unused varsys call C_F_POINTER( varsys%method%val(fun%mypos)%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 do iStFun = 1, fPtr%nVals if (fPtr%val(iStFun)%subTree%useGlobalMesh) then res(:) = res(:) & & + tem_spacetime_for( me = fPtr%val(iStFun), & & treeIDs = tree%treeID(elemPos(:)), & & tree = tree, & & time = time, & & n = nElems ) else 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 ) if (pos > 0) then oldPos = pos treeID = tree%treeID(elemPos(iElem)) st_res = tem_spacetime_for( me = fPtr%val(iStFun), & & treeIDs = treeID, & & tree = tree, & & time = time, & & n = 1 ) res(iElem) = res(iElem) + st_res(1) end if ! elemPos(iElem) part of subtree end do ! iElem end if ! global mesh end do ! nr. st_funs end subroutine evaluate_add_spacetime_scalarByTreeID