Write a spatial representation for list of points 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
Type | Intent | Optional | 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 |
subroutine hvs_asciiSpatial_dump_point_data( asciiSpatial, varPos, varSys, &
& bary, mesh, subTree, time )
! ---------------------------------------------------------------------------
!> 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
! ---------------------------------------------------------------------------
integer :: nVars, nPoints, nScalars, pointsOff, nChunkPoints
integer :: iPoint, iChunk, iScalar, counter
integer :: buf_start, buf_end
real(kind=rk), allocatable :: res(:)
real(kind=rk), allocatable :: points(:,:)
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)
if (present(subTree)) then
nPoints = subTree%nPoints
else
nPoints = mesh%nElems
end if
! allocate points to size of chunkSize
allocate(points(asciiSpatial%chunkSize,3))
! Process all chunks to derive the quantities defined in the tracking object
do iChunk = 1, asciiSpatial%nChunks
! Number of points read so far in previous chunks.
pointsOff = ((iChunk-1)*asciiSpatial%chunkSize)
! number of points written to THIS chunk
nChunkPoints = min(asciiSpatial%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 )
! Then gather contents into buffer, and write buffer to file
buffer = ''
do iPoint = 1, nChunkPoints
! write coordinates to buffer
write( buffer, '(3(1x,e24.16e3))' ) bary(pointsOff+iPoint, 1:3)
! append values in chuck to buffer
do iScalar = 1, nScalars
write( buffer, '(a,1x,e24.16e3)' ) trim(buffer), &
& res( (iPoint-1)*nScalars + iScalar )
end do
! write buffer to file
write ( asciiSpatial%outUnit , '(a)' ) trim(buffer)
end do !nChunkElems
end do ! iChunk
deallocate(points)
deallocate(res)
end subroutine hvs_asciiSpatial_dump_point_data