Writes particle data to file in parallel using SILO.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(particle_set), | intent(inout) | :: | this |
Point particles to write |
||
| 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 |
impure subroutine particle_set_WriteSilo(this,iter,time,list) !> Writes particle data to file in parallel using SILO. implicit none class(particle_set), intent(inout) :: this !! Point particles to write 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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 (particle_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 particle_set_WriteSilo