This subroutine call setup indices of input_variable
the interface has to comply to the abstract interface tem_varsys_module#tem_varsys_proc_setupIndices.
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. |
||
real(kind=rk), | intent(in) | :: | point(:,:) |
List of space coordinate points to store as growing array in method_data |
||
character(len=1), | intent(in), | optional | :: | offset_bit(:) |
Offset bit encoded as character for every point. Offset integer coord(3) is converted into a character with offset_bit = achar( (coord(1)+1) + (coord(2)+1)4 + (coord(3)+1)16 ) Backward transformation form character to 3 integer: coord(1) = mod(ichar(offset_bit),4) - 1 coord(2) = mod(ichar(offset_bit),16)/4 - 1 coord(3) = ichar(offset_bit)/16 - 1 If not present default is to center i.e offset_bit = achar(1+4+16) |
|
integer, | intent(in) | :: | iLevel |
Level to which input points belong to |
||
type(treelmesh_type), | intent(in) | :: | tree |
global treelm mesh info |
||
integer, | intent(in) | :: | nPnts |
Number of points to add in method_data of this variable |
||
integer, | intent(out) | :: | idx(:) |
Index of points in the growing array and variable val array. Size: nPoints This must be stored in boundary or source depends on who calls this routine. This index is required to return a value using getValOfIndex. |
recursive subroutine tem_opVar_setupIndices( fun, varSys, point, offset_bit, & & iLevel, tree, nPnts, idx ) ! ---------------------------------------------------------------------- ! !> 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 !> List of space coordinate points to store as growing array in !! method_data real(kind=rk), intent(in) :: point(:,:) !> Offset bit encoded as character for every point. !! !! Offset integer coord(3) is converted into a character with !! offset_bit = achar( (coord(1)+1) + (coord(2)+1)*4 + (coord(3)+1)*16 ) !! Backward transformation form character to 3 integer: !! coord(1) = mod(ichar(offset_bit),4) - 1 !! coord(2) = mod(ichar(offset_bit),16)/4 - 1 !! coord(3) = ichar(offset_bit)/16 - 1 !! !! If not present default is to center i.e offset_bit = achar(1+4+16) character, optional, intent(in) :: offset_bit(:) !> Level to which input points belong to integer, intent(in) :: iLevel !> global treelm mesh info type(treelmesh_type), intent(in) :: tree !> Number of points to add in method_data of this variable integer, intent(in) :: nPnts !> Index of points in the growing array and variable val array. !! Size: nPoints !! !! This must be stored in boundary or source depends on who !! calls this routine. !! This index is required to return a value using getValOfIndex. integer, intent(out) :: idx(:) ! ---------------------------------------------------------------------- ! type(tem_varSys_op_data_type), pointer :: opData integer :: iPnt, iDep type(grw_intArray_type), allocatable :: inputIndex_loc(:) integer, allocatable :: idxPerPnt(:) ! ---------------------------------------------------------------------- ! call C_F_POINTER( fun%method_Data, opData ) ! allcoate the index array for all inpits if (.not. allocated(opData%input_pntIndex)) then allocate( opData%input_pntIndex(fun%nInputs) ) end if ! allocate temporary inputIndex with size of nInputs and initialize ! growing array with length nPnts allocate(inputIndex_loc(fun%nInputs)) ! Now fill in the index arrays for the inputs call tem_opVar_fill_inputIndex( fun = fun, & & varSys = varSys, & & point = point, & & offset_bit = offset_bit, & & iLevel = iLevel, & & tree = tree, & & nPnts = nPnts, & & inputIndex = inputIndex_loc ) ! KM: Workaround for intel compiler in SuperMUC, the pointer ! gets corrupted after recursive call to depend variable call C_F_POINTER( fun%method_Data, opData ) ! initialize index with zero to identify points which does not ! belong to subTree allocate(idxPerPnt(fun%nInputs)) idx = 0 do iPnt = 1, nPnts do iDep = 1, fun%nInputs idxPerPnt(iDep) = inputIndex_loc(iDep)%val(iPnt) end do ! set index only when any of dependent variable has valid index if (any(idxPerPnt > 0)) then do iDep = 1, fun%nInputs call append(me = opData%input_pntIndex(iDep)%indexLvl(iLevel), & & val = inputIndex_loc(iDep)%val(iPnt) ) end do ! set index to last position in input_pntIndex of dep var 1 of ! indexLvl of iLevel idx(iPnt) = opData%input_pntIndex(1)%indexLvl(iLevel)%nVals end if end do do iDep = 1, fun%nInputs call truncate (opData%input_pntIndex(iDep)%indexLvl(iLevel) ) end do end subroutine tem_opVar_setupIndices