MPI_IO_Pack Subroutine

public subroutine MPI_IO_Pack(nout)

Uses

  • proc~~mpi_io_pack~~UsesGraph proc~mpi_io_pack MPI_IO_Pack module~basicmod basicmod proc~mpi_io_pack->module~basicmod module~mpimod mpimod proc~mpi_io_pack->module~mpimod module~config config module~basicmod->module~config module~mpimod->module~config mpi mpi module~mpimod->mpi

Output the grid and variables. The grid data contain the information of cell center and cell edge. The variable data contains that of cell center. Include the information of the cell edge! Information of meta data. the final grid point is not necessary but outputed.
the final grid is necessary.

the final grid point is not necessary but outputed. the final grid is necessary.

the final grid point is not necessary but outputed. the final grid is necessary.

The cell center value

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nout

Calls

proc~~mpi_io_pack~~CallsGraph proc~mpi_io_pack MPI_IO_Pack makedirs makedirs proc~mpi_io_pack->makedirs

Called by

proc~~mpi_io_pack~~CalledByGraph proc~mpi_io_pack MPI_IO_Pack proc~output Output proc~output->proc~mpi_io_pack

Source Code

  subroutine MPI_IO_Pack(nout)
    !! Output the grid and variables.
    !! The grid data contain the information of cell center and cell edge.
    !! The variable data contains that of  cell center.
    use basicmod
    use mpimod
    implicit none
    integer::i,j,k,n
    integer,intent(in):: nout
    integer::iee,jee,kee
    character(20),parameter::dirname=datadir
    character(40)::filename
    integer,parameter::unitout=17
    integer,parameter:: gs=1
     
    logical, save:: is_inited
    data is_inited /.false./

    iee = ie
    jee = je
    kee = ke
  
    !> Include the information of the cell edge!
    if(coords(1) .eq. ntiles(1)-1) iee = ie+1
    if(coords(2) .eq. ntiles(2)-1) jee = je+1
    if(coords(3) .eq. ntiles(3)-1) kee = ke+1
    
    if (.not. is_inited) then
       npart(1) = ngrid1
       npart(2) = ngrid2
       npart(3) = ngrid3
       
       ntotal(1) = ngrid1*ntiles(1)
       ntotal(2) = ngrid2*ntiles(2)
       ntotal(3) = ngrid3*ntiles(3)
     
       allocate(gridX(1:iee-is+1,nvarg))
       allocate(gridY(1:jee-js+1,nvarg))
       allocate(gridZ(1:kee-ks+1,nvarg))

       allocate(data3D(ngrid1,ngrid2,ngrid3,nvars))
         
       call makedirs(dirname)
       is_inited =.true.
      endif
      
      !> Information of meta data.
      if(myid_w == 0)then
         write(filename,'(a3,i5.5,a4)')"unf",nout,".dat"
         filename = trim(dirname)//filename
         
         open(unitout,file=filename,status='replace',form='formatted')
         write(unitout,'(a2,2(1x,E12.3))') "# ",time,dt
         write(unitout,'(a2,1x,i5)') "# ",ngrid1*ntiles(1)
         write(unitout,'(a2,1x,i5)') "# ",ngrid2*ntiles(2)
         write(unitout,'(a2,1x,i5)') "# ",ngrid3*ntiles(3)
         close(unitout)
      endif
      
      gridX(1:iee-is+1,1) = x1b(is:iee) !! the final grid point is not necessary but outputed.  
      gridX(1:iee-is+1,2) = x1a(is:iee) !! the final grid is necessary. 
      
      gridY(1:jee-js+1,1) = x2b(js:jee) !! the final grid point is not necessary but outputed.
      gridY(1:jee-js+1,2) = x2a(js:jee) !! the final grid is necessary.
      
      gridZ(1:kee-ks+1,1) = x3b(ks:kee) !! the final grid point is not necessary but outputed.
      gridZ(1:kee-ks+1,2) = x3a(ks:kee) !! the final grid is necessary.

      !> The cell center value 
      data3D(1:npart(1),1:npart(2),1:npart(3), 1) =  d(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 2) = v1(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 3) = v2(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 4) = v3(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 5) = b1(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 6) = b2(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 7) = b3(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 8) = bp(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3), 9) =  p(is:ie,js:je,ks:ke)
      data3D(1:npart(1),1:npart(2),1:npart(3),10) = gp(is:ie,js:je,ks:ke)
      do n=1,ncomp
         data3D(1:npart(1),1:npart(2),1:npart(3),nvars-ncomp+n) = Xcomp(n,is:ie,js:je,ks:ke)
      enddo
      return
      !         write(6,*) "bpf2",nflux2(mbps,i,j,k)
  end subroutine MPI_IO_Pack