Initialize time reduction operation variable Loop over all variable in varSys and allocate redTrans%val for reduction_transient operation variable with nElems
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_varSys_type), | intent(in) | :: | varSys |
Global variable system |
||
type(treelmesh_type), | intent(in) | :: | tree |
treelmesh_type |
||
type(tem_varMap_type), | intent(out) | :: | redTransVarMap |
position of time reduction variable in varSys |
||
integer, | intent(in), | optional | :: | nDofs |
Solver nDegrees of freedom |
|
type(tem_time_type), | intent(in) | :: | time |
Current time |
subroutine tem_opVar_reduction_transient_init(varSys, tree, redTransVarMap,& & nDofs, time) ! ------------------------------------------------------------------------- !> Global variable system type(tem_varSys_type), intent(in) :: varSys !> treelmesh_type type(treelmesh_type), intent(in) :: tree !> position of time reduction variable in varSys type(tem_varMap_type), intent(out) :: redTransVarMap !> Solver nDegrees of freedom integer, intent(in), optional :: nDofs !> Current time type(tem_time_type), intent(in) :: time ! ------------------------------------------------------------------------- type(tem_varSys_op_data_type), pointer :: opData integer :: iVar, iElem, varPos, nDofs_loc, posDepVar, nCompMax, idxMax integer :: nRedVars type(grw_labelArray_type) :: redTransVarName integer :: elemPos(tree%nElems) real(kind=rk), allocatable :: input_varRes(:) ! ------------------------------------------------------------------------- if (present(nDofs)) then nDofs_loc = nDofs else nDofs_loc = 1 end if ! Gather list of variable names which has reduction_transient operation call init(redTransVarName) do iVar = 1, varSys%varName%nVals if (trim(varSys%method%val(iVar)%operType) == 'reduction_transient') then call append(me=redTransVarName, val=varSys%varName%val(iVar)) end if end do ! create varMap to store position of reduction_transient variable in varSys call tem_create_varMap( & & varName = redTransVarName%val(1:redTransVarName%nVals), & & varSys = varSys, & & varMap = redTransVarMap ) elemPos(1:tree%nElems) = (/ (iElem, iElem=1, tree%nElems) /) nRedVars = redTransVarMap%varPos%nVals nCompMax = maxval(varSys%method%val(redTransVarMap%varPos%val(1:nRedVars)) & & %nComponents) allocate(input_varRes(tree%nElems*nCompMax*nDofs_loc)) ! Initialize time reduction do iVar = 1, redTransVarMap%varPos%nVals varPos = redTransVarmap%varPos%val(iVar) call C_F_POINTER(varSys%method%val(varPos)%method_data, opData) call tem_reduction_transient_init( & & me = opData%redTrans, & & nElems = tree%nElems, & & nComponents = varSys%method%val(varPos)%nComponents, & & nDofs = nDofs_loc ) ! Fill last posDepVar = varSys%method%val(varPos)%input_varPos(1) call varSys%method%val(posDepVar)%get_element( & & varSys = varSys, & & elemPos = elemPos, & & time = time, & & tree = tree, & & nElems = tree%nElems, & & nDofs = nDofs_loc, & & res = input_varRes(:) ) idxMax = opData%redTrans%nEntries opData%redTrans%val(:, opData%redTrans%last) = input_varRes(1:idxMax) end do end subroutine tem_opVar_reduction_transient_init