ResPart_set_WriteSILO Subroutine

private impure subroutine ResPart_set_WriteSILO(this, iter, time, list)

Writes data to disk in Silo format.

Type Bound

ResPart_set

Arguments

Type IntentOptional Attributes Name
class(ResPart_set), intent(inout) :: this

Set of resolved particles

integer, intent(in) :: iter

Iteration at write

real(kind=wp), intent(in) :: time

Time at write

character(len=str8), intent(in), optional :: list(:)

Names of components to write


Calls

proc~~respart_set_writesilo~~CallsGraph proc~respart_set_writesilo ResPart_set%ResPart_set_WriteSILO particles particles proc~respart_set_writesilo->particles proc~silo_obj_final silo_obj%silo_obj_Final proc~respart_set_writesilo->proc~silo_obj_final proc~silo_obj_init silo_obj%silo_obj_Init proc~respart_set_writesilo->proc~silo_obj_init proc~silo_obj_newtimestep silo_obj%silo_obj_NewTimeStep proc~respart_set_writesilo->proc~silo_obj_newtimestep proc~silo_obj_writelagrangianmesh silo_obj%silo_obj_WriteLagrangianMesh proc~respart_set_writesilo->proc~silo_obj_writelagrangianmesh v v proc~respart_set_writesilo->v w w proc~respart_set_writesilo->w dbclose dbclose proc~silo_obj_final->dbclose proc~parallel_obj_rankisroot parallel_obj%parallel_obj_RankIsRoot proc~silo_obj_final->proc~parallel_obj_rankisroot dbsetemptyok dbsetemptyok proc~silo_obj_init->dbsetemptyok proc~silo_obj_init->proc~parallel_obj_rankisroot proc~silo_obj_creategroups silo_obj%silo_obj_CreateGroups proc~silo_obj_init->proc~silo_obj_creategroups proc~stringtool_obj_removeextension stringtool_obj%stringtool_obj_RemoveExtension proc~silo_obj_init->proc~stringtool_obj_removeextension proc~sysutils_obj_createdirectory sysutils_obj%sysutils_obj_CreateDirectory proc~silo_obj_init->proc~sysutils_obj_createdirectory dbcreate dbcreate proc~silo_obj_newtimestep->dbcreate dbset2dstrlen dbset2dstrlen proc~silo_obj_newtimestep->dbset2dstrlen mpi_barrier mpi_barrier proc~silo_obj_newtimestep->mpi_barrier proc~silo_obj_newtimestep->proc~parallel_obj_rankisroot proc~silo_obj_setupgroupfiles silo_obj%silo_obj_SetupGroupFiles proc~silo_obj_newtimestep->proc~silo_obj_setupgroupfiles proc~silo_obj_newtimestep->proc~sysutils_obj_createdirectory dbaddiopt dbaddiopt proc~silo_obj_writelagrangianmesh->dbaddiopt proc~silo_obj_writelagrangianmesh->dbclose dbfreeoptlist dbfreeoptlist proc~silo_obj_writelagrangianmesh->dbfreeoptlist dbmkoptlist dbmkoptlist proc~silo_obj_writelagrangianmesh->dbmkoptlist dbopen dbopen proc~silo_obj_writelagrangianmesh->dbopen dbputmmesh dbputmmesh proc~silo_obj_writelagrangianmesh->dbputmmesh dbputpm dbputpm proc~silo_obj_writelagrangianmesh->dbputpm dbsetdir dbsetdir proc~silo_obj_writelagrangianmesh->dbsetdir proc~silo_obj_writelagrangianmesh->mpi_barrier proc~silo_obj_writelagrangianmesh->proc~parallel_obj_rankisroot mpi_comm_rank mpi_comm_rank proc~silo_obj_creategroups->mpi_comm_rank mpi_comm_split mpi_comm_split proc~silo_obj_creategroups->mpi_comm_split proc~silo_obj_setupgroupfiles->dbclose proc~silo_obj_setupgroupfiles->dbcreate proc~silo_obj_setupgroupfiles->mpi_barrier dbmkdir dbmkdir proc~silo_obj_setupgroupfiles->dbmkdir mpi_comm_size mpi_comm_size proc~silo_obj_setupgroupfiles->mpi_comm_size

Called by

proc~~respart_set_writesilo~~CalledByGraph proc~respart_set_writesilo ResPart_set%ResPart_set_WriteSILO proc~cdifs_obj_writeoutputdata cdifs_obj_WriteOutputData proc~cdifs_obj_writeoutputdata->proc~respart_set_writesilo proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->proc~respart_set_writesilo interface~cdifs_obj_writeoutputdata cdifs_obj%cdifs_obj_WriteOutputData interface~cdifs_obj_writeoutputdata->proc~cdifs_obj_writeoutputdata interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData interface~grans_obj_writeoutputdata->proc~grans_obj_writeoutputdata proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->interface~cdifs_obj_writeoutputdata proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver

Source Code

    impure subroutine ResPart_set_WriteSILO(this,iter,time,list)
      !> Writes data to disk in Silo format.
      implicit none
      class(ResPart_set), intent(inout) :: this                                !! Set of resolved particles
      integer,            intent(in)    :: iter                                !! Iteration at write
      real(wp),           intent(in)    :: time                                !! Time at write
      character(str8),    intent(in),    &
                               optional :: list(:)                             !! Names of components to write
      ! Work variables
      type(silo_obj)       :: silo
      real(wp), allocatable:: x1(:)
      real(wp), allocatable:: x2(:)
      real(wp), allocatable:: x3(:)
      real(wp), allocatable:: var(:)
      integer              :: n

      call silo%Initialize(trim(adjustl(this%write_file)),"W",this%parallel)

      ! Add new timestep
      call silo%NewTimeStep(time)

      ! Get Lagrangian positions
      allocate(x1(this%count_))
      allocate(x2(this%count_))
      allocate(x3(this%count_))
      if (this%count_.gt.0) then
        x1 = this%p(1:this%count_)%p(1)
        x2 = this%p(1:this%count_)%p(2)
        x3 = this%p(1:this%count_)%p(3)
      end if

      ! Write Lagrangian Mesh
      call silo%WriteLagrangianMesh(x1,x2,x3,iter,time)

      deallocate(x1,x2,x3)

      ! Allocate buffer
      allocate(var(this%count_))

      if (.not.present(list)) then
        ! Write only IDs
        if (this%count_.gt.0) then
          select type (particles =>this%p)
          type is (ResPart_obj)
            var = particles(1:this%count_)%id
          end select
        end if

        call silo%Write('id',var)
      else
        do  n=1,size(list)
          select case (trim(adjustl(list(n))))
          case ('id' )
            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%id
              end select
            end if

            call silo%Write('id', var)
          case ('d'  )
            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%d
              end select
            end if

            call silo%Write('d',  var)
          case ('rho')
            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%rho
              end select
            end if

            call silo%Write('rho',var)
          case ('v'  )
            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%v(1)
              end select
            end if

            call silo%Write('v1', var)

            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%v(2)
              end select
            end if

            call silo%Write('v2', var)

            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%v(3)
              end select
            end if

            call silo%Write('v3' ,var)
          case ('w'  )
            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%w(1)
              end select
            end if

            call silo%Write('w1', var)

            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%w(2)
              end select
            end if

            call silo%Write('w2', var)

            if (this%count_.gt.0) then
              select type (particles =>this%p)
              type is (ResPart_obj)
                var = particles(1:this%count_)%w(3)
              end select
            end if

            call silo%Write('w3' ,var)
          end select

        end do
      end if

      ! Deallocate buffer
      deallocate(var)

      call silo%Finalize()

      return
    end subroutine ResPart_set_WriteSILO