hvs_ascii_dump_point_data Subroutine

public subroutine hvs_ascii_dump_point_data(ascii, outProc, varpos, varsys, mesh, time, subtree)

Write single log for the right scheme into its ascii file. This routine calls the get_point routine and dumps the exact point data for the point specified in the tracking table

Arguments

Type IntentOptional Attributes Name
type(hvs_ascii_type), intent(inout) :: ascii

ascii file to write data to.

type(tem_comm_env_type), intent(in) :: outProc

Parallel environment to use for the output.

integer, intent(in) :: varpos(:)

Position of the variable to write

type(tem_varSys_type), intent(in) :: varsys

Description of the available variable system to get the given varnames from.

type(treelmesh_type), intent(in) :: mesh

Mesh to write the data on.

type(tem_time_type), intent(in) :: time

Point in time to use for this data.

Can be important for space-time function evaluations.

type(tem_subTree_type), intent(in), optional :: subtree

Optional restriction of the points to output. Contains array of points passed in the config to output.


Calls

proc~~hvs_ascii_dump_point_data~~CallsGraph proc~hvs_ascii_dump_point_data hvs_ascii_dump_point_data proc~tem_reduction_spatial_append tem_reduction_spatial_append proc~hvs_ascii_dump_point_data->proc~tem_reduction_spatial_append proc~tem_get_point_chunk tem_get_point_chunk proc~hvs_ascii_dump_point_data->proc~tem_get_point_chunk proc~tem_reduction_spatial_close tem_reduction_spatial_close proc~hvs_ascii_dump_point_data->proc~tem_reduction_spatial_close proc~tem_baryofid tem_BaryOfId proc~hvs_ascii_dump_point_data->proc~tem_baryofid proc~tem_reduction_spatial_open tem_reduction_spatial_open proc~hvs_ascii_dump_point_data->proc~tem_reduction_spatial_open proc~tem_reduction_spatial_tochunk tem_reduction_spatial_toChunk proc~hvs_ascii_dump_point_data->proc~tem_reduction_spatial_tochunk proc~tem_elemsize tem_ElemSize proc~tem_reduction_spatial_append->proc~tem_elemsize proc~tem_levelof tem_LevelOf proc~tem_reduction_spatial_append->proc~tem_levelof mpi_reduce mpi_reduce proc~tem_reduction_spatial_close->mpi_reduce proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_baryofid->proc~tem_elemsizelevel proc~tem_coordofid tem_CoordOfId proc~tem_baryofid->proc~tem_coordofid proc~tem_elemsize->proc~tem_elemsizelevel proc~tem_elemsize->proc~tem_levelof proc~tem_coordofid->proc~tem_levelof

Called by

proc~~hvs_ascii_dump_point_data~~CalledByGraph proc~hvs_ascii_dump_point_data hvs_ascii_dump_point_data proc~hvs_output_write hvs_output_write proc~hvs_output_write->proc~hvs_ascii_dump_point_data proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_write

Contents


Source Code

  subroutine hvs_ascii_dump_point_data( ascii, outProc, varPos, varSys, mesh,  &
    &                                   time, subTree )
    ! ---------------------------------------------------------------------------
    !> ascii file to write data to.
    type(hvs_ascii_type), intent(inout) :: ascii

    !> Parallel environment to use for  the output.
    type(tem_comm_env_type ), intent(in)    :: outProc

    !> Position of the variable to write
    integer, intent(in) :: varpos(:)

    !> Description of the available variable system to get the given varnames
    !! from.
    type(tem_varSys_type), intent(in) :: varsys

    !> Mesh to write the data on.
    type(treelmesh_type), intent(in) :: mesh

    !> Point in time to use for this data.
    !!
    !! Can be important for space-time function evaluations.
    type(tem_time_type), intent(in) :: time

    !> Optional restriction of the points to output.
    !! Contains array of points passed in the config  to output.
    type(tem_subtree_type), optional, intent(in) :: subtree
    ! ---------------------------------------------------------------------------
    integer :: nVars, nPoints, nScalars, pointsOff, nChunkPoints
    integer :: iPoint, iChunk, iScalar, counter
    integer :: buf_start, buf_end
    real(kind=rk), allocatable :: res(:)
    character(len=4000) :: log_output ! output buffer
    real(kind=rk), allocatable :: points(:,:)
    ! ---------------------------------------------------------------------------
    if (present(subTree)) then
      nPoints = subTree%nPoints
    else
      nPoints = mesh%nElems
    end if

    allocate(res(io_buffer_size))

    ! Number of variables to dump
    nVars = size(varPos)

    ! Number of scalars in current output
    nScalars = sum(varSys%method%val(varPos(:))%nComponents)

    if (ascii%isReduce) then
      ! open spatial reduction
      call tem_reduction_spatial_open( me     = ascii%redSpatial, &
        &                              varSys = varSys,           &
        &                              varPos = varPos(:nVars)    )
    end if

    ! allocate points to size of chunkSize
    allocate(points(ascii%chunkSize,3))

    ! Process all chunks to derive the quantities defined in the tracking object
    do iChunk = 1, ascii%nChunks
      ! Number of points read so far in previous chunks.
      pointsOff = ((iChunk-1)*ascii%chunkSize)

      ! number of points written to THIS chunk
      nChunkPoints = min(ascii%chunkSize, nPoints-pointsOff)

      ! Compute the points lower and upper bound for the current chunk
      buf_start = pointsOff + 1
      buf_end = pointsOff + nChunkPoints

      if (present(subTree)) then
        points(1:nChunkPoints,:) = subTree%points(buf_start:buf_end,:)
      else
        counter = 0
        do iPoint = buf_start, buf_end
          counter = counter + 1
          points(counter, :) = tem_BaryOfId( mesh,               &
            &                                mesh%treeID(iPoint) )
        end do
      end if

      ! evaluate all variables on current chunk
      call tem_get_point_chunk(  varSys  = varSys,                   &
        &                        varPos  = varPos,                   &
        &                        point   = points(1:nChunkPoints,:), &
        &                        time    = time,                     &
        &                        tree    = mesh,                     &
        &                        nPnts   = nChunkPoints,             &
        &                        res     = res                       )

      ! preform spatial reduction
      if (ascii%isReduce) then
       ! NA : The reduction for the get points needs to be thought
       ! through as we would like at one point to have all the segment paoints
       ! stored and reduction applied. In that case, l2 norm and weighted sum
       ! requires the treeID information. Maybe we should write a new routine to
       ! do that tem_reduction_spatial_append_point
       ! KM: changed treeID argument optional because for point values
       ! no need to do volume reduction
       call tem_reduction_spatial_append( me     = ascii%redSpatial, &
         &                                chunk  = res,              &
         &                                nElems = nChunkPoints,     &
         &                                tree   = mesh,             &
         &                                varSys = varSys,           &
         &                                varPos = varPos            )
      end if
    end do ! iChunk

    ! If a reduction is present in the tracking, then the output is
    ! changed completely to only the reduced values
    ! For each variable, a reduction has to be assigned
    if( ascii%isReduce ) then

      ! Perform the global reduction on the data which was appended
      ! inside trackVariable by tem_reduction_spatial_append
      call tem_reduction_spatial_close( me   = ascii%redSpatial, &
        &                               proc = outProc           )

      ! Re-assign the chunk here to the stuff which was produced in the
      ! reduction operation
      call tem_reduction_spatial_toChunk(me          = ascii%redSpatial, &
        &                                chunk       = res,              &
        &                                nChunkElems = nChunkPoints      )
    end if

    ! If reduction is active only root of this output
    ! dumps data else all process writes their own result file
    if ( (ascii%isReduce .and. outProc%rank == 0) .or. &
      & (.not. ascii%isReduce) ) then

      ! First assemble the complete row consisting of the time ...
      write( log_output, '(e24.16e3)' ) time%sim

      ! add all the scalars entries of the variable systems for each elements

      do iPoint = 1, nChunkPoints
          do iScalar = 1, nScalars
            write( log_output, '(a,1x,e24.16e3)' ) trim(log_output),          &
              &    res( (iPoint-1)*nScalars +  iScalar )
        end do
      end do
      ! then write into the ascii file
      write ( ascii%outunit , '(a)' ) trim(log_output)
    end if

    deallocate(points)
    deallocate(res)

  end subroutine hvs_ascii_dump_point_data