VTK_GEO_XML_STRG_R8 Function

private function VTK_GEO_XML_STRG_R8(nx1, nx2, ny1, ny2, nz1, nz2, NN, X, Y, Z) result(E_IO)

Function for saving mesh; topology = StructuredGrid (R8P).

Arguments

Type IntentOptional AttributesName
integer(kind=I4P), intent(in) :: nx1
integer(kind=I4P), intent(in) :: nx2
integer(kind=I4P), intent(in) :: ny1
integer(kind=I4P), intent(in) :: ny2
integer(kind=I4P), intent(in) :: nz1
integer(kind=I4P), intent(in) :: nz2
integer(kind=I4P), intent(in) :: NN
real(kind=R8P), intent(in) :: X(1:NN)
real(kind=R8P), intent(in) :: Y(1:NN)
real(kind=R8P), intent(in) :: Z(1:NN)

Return Value integer(kind=I4P)


Calls

proc~~vtk_geo_xml_strg_r8~~CallsGraph proc~vtk_geo_xml_strg_r8 VTK_GEO_XML_STRG_R8 interface~str str proc~vtk_geo_xml_strg_r8->interface~str proc~str_i8p str_I8P interface~str->proc~str_i8p proc~str_r8p str_R8P interface~str->proc~str_r8p proc~str_i1p str_I1P interface~str->proc~str_i1p proc~str_i4p str_I4P interface~str->proc~str_i4p proc~str_r16p str_R16P interface~str->proc~str_r16p proc~str_i2p str_I2P interface~str->proc~str_i2p proc~str_r4p str_R4P interface~str->proc~str_r4p

Called by

proc~~vtk_geo_xml_strg_r8~~CalledByGraph proc~vtk_geo_xml_strg_r8 VTK_GEO_XML_STRG_R8 interface~vtk_geo_xml VTK_GEO_XML interface~vtk_geo_xml->proc~vtk_geo_xml_strg_r8

Contents

Source Code


Source Code

  function VTK_GEO_XML_STRG_R8(nx1,nx2,ny1,ny2,nz1,nz2,NN,X,Y,Z) result(E_IO)
  !---------------------------------------------------------------------------------------------------------------------------------
  !! Function for saving mesh; topology = StructuredGrid (R8P).
  !---------------------------------------------------------------------------------------------------------------------------------

  !---------------------------------------------------------------------------------------------------------------------------------
  implicit none
  integer(I4P), intent(IN):: nx1,nx2  ! initial and final nodes of x axis
  integer(I4P), intent(IN):: ny1,ny2  ! initial and final nodes of y axis
  integer(I4P), intent(IN):: nz1,nz2  ! initial and final nodes of z axis
  integer(I4P), intent(IN):: NN       ! number of all nodes
  real(R8P),    intent(IN):: X(1:NN)  ! x coordinates
  real(R8P),    intent(IN):: Y(1:NN)  ! y coordinates
  real(R8P),    intent(IN):: Z(1:NN)  ! z coordinates
  integer(I4P)::             E_IO     ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
  character(len=maxlen)::    s_buffer ! buffer string
  integer(I4P)::             n1       ! counter
  !---------------------------------------------------------------------------------------------------------------------------------

  !---------------------------------------------------------------------------------------------------------------------------------
  select case(f_out)
  case(f_out_ascii)
    write(unit=Unit_VTK,fmt='(A,6'//FI4P//',A)',iostat=E_IO)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
    indent = indent + 2
    write(unit=Unit_VTK,fmt='(A)',iostat=E_IO)repeat(' ',indent)//'<Points>'
    indent = indent + 2
    write(unit=Unit_VTK,fmt='(A)',iostat=E_IO)repeat(' ',indent)// &
                                              '<DataArray type="Float64" NumberOfComponents="3" Name="Point" format="ascii">'
    write(unit=Unit_VTK,fmt='(3'//FR8P//')',iostat=E_IO)(X(n1),Y(n1),Z(n1),n1=1,NN)
    write(unit=Unit_VTK,fmt='(A)',iostat=E_IO)repeat(' ',indent)//'</DataArray>'
    indent = indent - 2
    write(unit=Unit_VTK,fmt='(A)',iostat=E_IO)repeat(' ',indent)//'</Points>'
  case(f_out_binary)
    write(s_buffer,fmt='(A,6'//FI4P//',A)',iostat=E_IO)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
    indent = indent + 2
    write(unit=Unit_VTK,iostat=E_IO)trim(s_buffer)//end_rec
    write(unit=Unit_VTK,iostat=E_IO)repeat(' ',indent)//'<Points>'//end_rec
    indent = indent + 2
    write(unit=Unit_VTK,iostat=E_IO)repeat(' ',indent)//                                                                        &
                                    '<DataArray type="Float64" NumberOfComponents="3" Name="Point" format="appended" offset="', &
                                    trim(str(.true.,ioffset)),                                                                  &
                                    '">'//                                                                                      &
                                    end_rec
    N_Byte  = 3*NN*sizeof(Tipo_R8)
    ioffset = ioffset + sizeof(Tipo_I4) + N_Byte
    write(unit=Unit_VTK_Append,iostat=E_IO)N_Byte,'R8',3*NN
    write(unit=Unit_VTK_Append,iostat=E_IO)(X(n1),Y(n1),Z(n1),n1=1,NN)
    write(unit=Unit_VTK,iostat=E_IO)repeat(' ',indent)//'</DataArray>'//end_rec
    indent = indent - 2
    write(unit=Unit_VTK,iostat=E_IO)repeat(' ',indent)//'</Points>'//end_rec
  endselect
  return
  !---------------------------------------------------------------------------------------------------------------------------------
  endfunction VTK_GEO_XML_STRG_R8