hvs_asciiSpatial_dump_elem_data Subroutine

public subroutine hvs_asciiSpatial_dump_elem_data(asciiSpatial, varpos, varSys, bary, mesh, subtree, time, nDofs)

Write a spatial representation for elements into an ascii tracking file

Each time this routine is called, a new file is written Filename: {tracking_folder}{tracking_label}spatial_{timestamp}.res e.g.: tracking/lineProbe_spatial_00001_01_01378.1.res Each process writes its own files

Arguments

Type IntentOptional Attributes Name
type(hvs_asciiSpatial_type), intent(inout) :: asciiSpatial

The file description to open

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

Positions of the variables to write

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

solver-provided variable systems

real(kind=rk), intent(in) :: bary(:,:)

Barycenter of elements

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

Mesh to write the data on.

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

Optional restriction of the elements to output.

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

current global time

integer, intent(in) :: nDofs

The number of dofs for each scalar variable of the equation system


Calls

proc~~hvs_asciispatial_dump_elem_data~~CallsGraph proc~hvs_asciispatial_dump_elem_data hvs_asciiSpatial_dump_elem_data proc~tem_get_element_chunk tem_get_element_chunk proc~hvs_asciispatial_dump_elem_data->proc~tem_get_element_chunk

Called by

proc~~hvs_asciispatial_dump_elem_data~~CalledByGraph proc~hvs_asciispatial_dump_elem_data hvs_asciiSpatial_dump_elem_data proc~hvs_output_write hvs_output_write proc~hvs_output_write->proc~hvs_asciispatial_dump_elem_data proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_write

Contents


Source Code

  subroutine hvs_asciiSpatial_dump_elem_data( asciiSpatial, varPos, varSys,    &
    &                        bary, mesh, subTree, time, nDofs )
    ! ---------------------------------------------------------------------------
    !> The file description to open
    type(hvs_asciiSpatial_type), intent(inout) :: asciiSpatial

    !> solver-provided variable systems
    type(tem_varSys_type), intent(in)       :: varSys

    !> Positions of the variables to write
    integer, intent(in) :: varpos(:)

    !> Barycenter of elements
    real(kind=rk), intent(in) :: bary(:,:)

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

    !> Optional restriction of the elements to output.
    type(tem_subtree_type), optional, intent(in) :: subtree

    !> current global time
    type(tem_time_type ), intent(in)  :: time

    !> The number of dofs for each scalar variable of the equation system
    integer, intent(in) :: nDofs
    ! ---------------------------------------------------------------------------
    integer :: nVars, nElems, nScalars, elemOff, nChunkElems, elemSize
    integer :: iElem, iChunk, iScalar, iDof
    integer :: buf_start, buf_end
    real(kind=rk), allocatable :: res(:)
    integer, allocatable :: elemPos(:)
    character(len=1024)  :: buffer
    ! ---------------------------------------------------------------------------
    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)

    ! Size of a single element
    elemsize = nScalars*nDofs

    if (present(subTree)) then
      nElems = subTree%nElems
    else
      nElems = mesh%nElems
    end if

    ! allocate elemPos to size of chunkSize
    allocate(elemPos(asciiSpatial%chunkSize))

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

      ! number of elements written to THIS chunk
      nChunkElems = min(asciiSpatial%chunkSize, nElems-elemOff)

      ! Compute the element lower and upper bound for the current chunk
      buf_start = elemOff + 1
      buf_end = elemOff + nChunkElems

      if (present(subTree)) then
        elemPos(1:nChunkElems) = subTree%map2Global(buf_start:buf_end)
      else
        elemPos(1:nChunkElems) = (/ (iElem, iElem=buf_start, buf_end) /)
      end if

      ! evaluate all variables on current chunk
      call tem_get_element_chunk(varSys  = varSys,                 &
        &                        varPos  = varPos,                 &
        &                        elemPos = elemPos(1:nChunkElems), &
        &                        time    = time,                   &
        &                        tree    = mesh,                   &
        &                        nElems  = nChunkElems,            &
        &                        nDofs   = nDofs,                  &
        &                        res     = res                     )

      ! Then gather contents into buffer, and write buffer to file
      buffer = ''
      do iElem = 1, nChunkElems
        ! write coordinates to buffer
        write( buffer, '(3(1x,e24.16e3))' ) bary(elemOff+iElem, 1:3)

        ! append values in chuck to buffer
        do iDof = 1, nDofs
          do iScalar = 1, nScalars
            write( buffer, '(a,1x,e24.16e3)' ) trim(buffer),           &
              &  res( (iElem-1)*elemSize + (iDof-1)*nScalars + iScalar )
          end do
        end do

        ! write buffer to file
        write ( asciiSpatial%outUnit , '(a)' ) trim(buffer)

      end do !nChunkElems

    end do ! iChunk

    deallocate(elemPos)
    deallocate(res)

  end subroutine hvs_asciiSpatial_dump_elem_data