tem_update_surfPos Subroutine

public subroutine tem_update_surfPos(me, levelDesc, globTree, movement, time, iLevel, IBMUnit, useInitPos, movPredef)

This subroutine updates the surface points and the parentIDs array as well as sets the correct property bits.

Arguments

Type IntentOptional Attributes Name
type(tem_surfData_type), intent(inout) :: me

datatype to store the surface information

type(tem_levelDesc_type), intent(inout) :: levelDesc

the level descriptor incl. ghost and halo elements as well as the communicator information on the level iLevel

type(treelmesh_type), intent(inout) :: globTree

global Tree information

type(tem_spacetime_fun_type) :: movement

spacetime function to define the motion of the surface points

type(tem_time_type) :: time

timing information

integer, intent(inout) :: iLevel

the current level

integer, intent(in), optional :: IBMUnit(0:tem_last_lu)

optional output log unit other than the global logUnit

logical, intent(in), optional :: useInitPos

shall the initial points be stored and used for updating the points later on ???

logical, intent(in) :: movPredef

logical to define wether the motion is predefined or not if not: initialize the values differently


Calls

proc~~tem_update_surfpos~~CallsGraph proc~tem_update_surfpos tem_update_surfPos interface~tem_spacetime_for tem_spacetime_for proc~tem_update_surfpos->interface~tem_spacetime_for proc~tem_init_surfdata tem_init_surfData proc~tem_update_surfpos->proc~tem_init_surfdata proc~tem_spacetime_for_coord tem_spacetime_for_coord interface~tem_spacetime_for->proc~tem_spacetime_for_coord proc~tem_spacetime_for_stcoord tem_spacetime_for_stcoord interface~tem_spacetime_for->proc~tem_spacetime_for_stcoord proc~tem_spacetime_for_treeids tem_spacetime_for_treeIDs interface~tem_spacetime_for->proc~tem_spacetime_for_treeids proc~tem_spacetime_scalar_for_index tem_spacetime_scalar_for_index interface~tem_spacetime_for->proc~tem_spacetime_scalar_for_index proc~tem_spacetime_vector_for_coord tem_spacetime_vector_for_coord interface~tem_spacetime_for->proc~tem_spacetime_vector_for_coord proc~tem_spacetime_vector_for_index tem_spacetime_vector_for_index interface~tem_spacetime_for->proc~tem_spacetime_vector_for_index proc~tem_spacetime_vector_for_treeids tem_spacetime_vector_for_treeIDs interface~tem_spacetime_for->proc~tem_spacetime_vector_for_treeids interface~append~24 append proc~tem_init_surfdata->interface~append~24 proc~tem_coordofreal tem_CoordOfReal proc~tem_init_surfdata->proc~tem_coordofreal proc~tem_idofcoord tem_IdOfCoord proc~tem_init_surfdata->proc~tem_idofcoord proc~tem_pathcomparison tem_PathComparison proc~tem_init_surfdata->proc~tem_pathcomparison proc~tem_pathof tem_PathOf proc~tem_init_surfdata->proc~tem_pathof proc~tem_treeidintotal tem_treeIDinTotal proc~tem_init_surfdata->proc~tem_treeidintotal

Source Code

  subroutine tem_update_surfPos( me, levelDesc, globTree, movement, time, &
    &                            iLevel, IBMUnit, useInitPos, movPredef )
    ! ---------------------------------------------------------------------------
    !> datatype to store the surface information
    type( tem_surfData_type ), intent(inout) :: me
    !> the level descriptor incl. ghost and halo elements as well as the
    !! communicator information on the level iLevel
    type( tem_levelDesc_type ), intent(inout) :: levelDesc
    !> global Tree information
    type( treelmesh_type ), intent(inout) :: globTree
    !> spacetime function to define the motion of the surface points
    type( tem_spacetime_fun_type ) :: movement
    !> timing information
    type(tem_time_type) :: time
    !> the current level
    integer, intent(inout) :: iLevel
    !> optional output log unit other than the global logUnit
    integer, optional, intent(in) :: IBMUnit(0:tem_last_lu)
    !> shall the initial points be stored and used for updating the points
    !! later on ???
    logical, optional, intent(in) :: useInitPos
    !> logical to define wether the motion is predefined or not
    !! if not: initialize the values differently
    logical, intent(in) :: movPredef
    ! ---------------------------------------------------------------------------
    ! counters
    integer :: iPoint
    real(kind=rk) :: pos(1,3)
    real(kind=rk) :: huge_real
    integer :: IBMUnit_loc(0:tem_last_lu)
    logical :: tmp_useInitPos
    ! ---------------------------------------------------------------------------

    if( present( useInitPos ))then
      tmp_useInitPos = useInitPos
    else
      tmp_useInitPos = .false.
    end if

    if( present( IBMUnit ))then
      IBMUnit_loc = IBMUnit
    else
      IBMUnit_loc = logUnit
    end if

    huge_real = huge( huge_real )

    ! loop over the surface points and ...
    do iPoint = 1, me%nUniquePoints_total
      if( movPredef )then
        ! only update points which have fluid elements as parents
        if( me%ParentIDs(iLevel)%ptrs( iPoint ) > 0 )then
          ! ... clean the IBM property bits from the element
          levelDesc%property( me%ParentIDs(iLevel)%ptrs( iPoint)) =              &
            &     ibclr( levelDesc%property( me%ParentIDs(iLevel)%ptrs( iPoint)),&
            &            prp_hasIBM )
        end if
        ! check wether the initial point coordinates shall be used ...
        if( tmp_useInitPos )then
          ! ... yes: use array backPointCoords and ...
          ! ... store the coordinates in a temporary variable
          pos(1,1:3) = me%backPointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 )
        else
          ! ... yes: use array pointCoords and ...
          ! ... store the coordinates in a temporary variable
          pos(1,1:3) = me%pointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 )
        end if
        ! ... apply the movement and store the new positions
        pos = tem_spacetime_for( me    = movement, &
          &                      coord = pos,      &
          &                      time  = time,     &
          &                      n     = 1,        &
          &                      nComp = 3         )
        ! ... copy back the new positions
        me%pointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 ) = pos(1,1:3)
      else ! .not. movPredef -> initialize Xk with non local fluid parents with
           ! huge
        ! only update points which have fluid elements as parents
        if( me%ParentIDs(iLevel)%ptrs( iPoint ) > 0 .and.                        &
          & me%ParentIDs(iLevel)%ptrs( iPoint ) <= levelDesc%elem%nElems( eT_fluid ) )then
          ! ... clean the IBM property bits from the element
          levelDesc%property( me%ParentIDs(iLevel)%ptrs( iPoint)) =              &
            &     ibclr( levelDesc%property( me%ParentIDs(iLevel)%ptrs( iPoint)),&
            &            prp_hasIBM )
          ! check wether the initial point coordinates shall be used ...
          if( tmp_useInitPos )then
            ! ... yes: use array backPointCoords and ...
            ! ... store the coordinates in a temporary variable
            pos(1,1:3) = me%backPointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 )
          else
            ! ... yes: use array pointCoords and ...
            ! ... store the coordinates in a temporary variable
            pos(1,1:3) = me%pointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 )
          end if
          ! ... apply the movement and store the new positions
          pos = tem_spacetime_for( me    = movement, &
            &                      coord = pos,      &
            &                      time  = time,     &
            &                      n     = 1,        &
            &                      nComp = 3         )
          ! ... copy back the new positions
          me%pointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 ) = pos(1,1:3)
        else
          ! set the coordinates to be infinity such that they are not recognized
          ! by this proc any more
          me%pointCoords( (iPoint-1)*3+1:(iPoint-1)*3+3 ) = (/huge_real,         &
            &                                                 huge_real,         &
            &                                                 huge_real/)
        end if
      end if
    end do ! iPoint

    ! update the parentIDs array and set the correct property bits
    call tem_init_surfData( me        = me,                                    &
      &                     levelDesc = levelDesc,                             &
      &                     globTree  = globTree,                              &
      &                     iLevel    = iLevel )

  end subroutine tem_update_surfPos