write_surfnode Subroutine

public subroutine write_surfnode(files, nodes, control, bc, dims)

Extract and write the wall surface node points in a file shared by all the MPI processes

Arguments

Type IntentOptional AttributesName
type(filetype), intent(in) :: files
type(nodetype), intent(in), dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3):: nodes
type(controltype), intent(in) :: control
type(boundarytype), intent(in) :: bc
type(extent), intent(in) :: dims

Calls

proc~~write_surfnode~~CallsGraph proc~write_surfnode write_surfnode mpi_file_write_shared mpi_file_write_shared proc~write_surfnode->mpi_file_write_shared mpi_bcast mpi_bcast proc~write_surfnode->mpi_bcast mpi_type_commit mpi_type_commit proc~write_surfnode->mpi_type_commit mpi_barrier mpi_barrier proc~write_surfnode->mpi_barrier mpi_gather mpi_gather proc~write_surfnode->mpi_gather mpi_type_contiguous mpi_type_contiguous proc~write_surfnode->mpi_type_contiguous mpi_file_close mpi_file_close proc~write_surfnode->mpi_file_close proc~setup_surface setup_surface proc~write_surfnode->proc~setup_surface proc~surface_points surface_points proc~write_surfnode->proc~surface_points proc~setup_surface->mpi_barrier proc~allocate_memory~3 allocate_memory proc~setup_surface->proc~allocate_memory~3 debugcall debugcall proc~setup_surface->debugcall proc~find_wall find_wall proc~setup_surface->proc~find_wall mpi_file_open mpi_file_open proc~setup_surface->mpi_file_open proc~link_aliases link_aliases proc~setup_surface->proc~link_aliases proc~surface_points->debugcall proc~allocate_memory~3->debugcall interface~alloc alloc proc~allocate_memory~3->interface~alloc proc~link_aliases->debugcall proc~alloc_rank2_real alloc_rank2_real interface~alloc->proc~alloc_rank2_real proc~alloc_rank4_real alloc_rank4_real interface~alloc->proc~alloc_rank4_real proc~alloc_rank3_real alloc_rank3_real interface~alloc->proc~alloc_rank3_real proc~alloc_rank2_integer alloc_rank2_integer interface~alloc->proc~alloc_rank2_integer proc~alloc_rank1_integer alloc_rank1_integer interface~alloc->proc~alloc_rank1_integer proc~alloc_rank6_real alloc_rank6_real interface~alloc->proc~alloc_rank6_real proc~alloc_rank5_real alloc_rank5_real interface~alloc->proc~alloc_rank5_real proc~alloc_rank3_integer alloc_rank3_integer interface~alloc->proc~alloc_rank3_integer proc~alloc_rank1_real alloc_rank1_real interface~alloc->proc~alloc_rank1_real

Called by

proc~~write_surfnode~~CalledByGraph proc~write_surfnode write_surfnode proc~setup_solver setup_solver proc~setup_solver->proc~write_surfnode proc~start_run start_run proc~start_run->proc~setup_solver program~main main program~main->proc~start_run

Contents

Source Code


Source Code

    subroutine write_surfnode(files, nodes, control, bc, dims)
      !< Extract and write the wall surface node points
      !< in a file shared by all the MPI processes
      implicit none
      type(filetype), intent(in) :: files
      type(controltype), intent(in) :: control
      type(extent), intent(in) :: dims
      type(boundarytype), intent(in) :: bc
      type(nodetype), dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3), intent(in) :: nodes
      integer :: count
      integer :: i

      imx = dims%imx
      jmx = dims%jmx
      kmx = dims%kmx

      call setup_surface(files, control, bc)
      call surface_points(nodes)
      call MPI_GATHER(n_wall, 1, MPI_Integer, n_wall_buf, 1, &
                      MPI_integer,0, MPI_COMM_WORLD, ierr)
      total_n_wall = sum(n_wall_buf(:))
      call MPI_Bcast(total_n_wall,1, MPI_Integer, 0, &
                       MPI_COMM_WORLD, ierr)
      call MPI_Bcast(n_wall_buf, control%total_process, MPI_Integer, 0, &
                       MPI_COMM_WORLD, ierr)

      write_flag=0
      count=0
      do i=1,control%total_process
        if(n_wall_buf(i)>0) then
          write_flag(i)=count
          count = count+1
        end if
      end do
      call MPI_TYPE_CONTIGUOUS(maxlen,MPI_Character, new_type, ierr)
      call MPI_TYPE_COMMIT(new_type, ierr)
      if(process_id==0)then
      write(line, '(I0)') total_n_wall
      line(len(line):len(line))=lf
      call MPI_FILE_WRITE_shared(thisfile, line, 1, &
                              new_type, &
                              MPI_STATUS_IGNORE, ierr)
      end if
      call mpi_barrier(MPI_COMM_WORLD,ierr)
      if(n_wall>0)then
      do i=1,n_wall
        write(line, '(3(ES18.10E3,4x))') wall_x(i), wall_y(i), wall_z(i)
        line(len(line):len(line))=lf
        call MPI_FILE_WRITE_shared(thisfile, line, 1, &
                                new_type, &
                                MPI_STATUS_IGNORE, ierr)
      end do
      end if
      call mpi_barrier(MPI_COMM_WORLD,ierr)
      call MPI_FILE_CLOSE(thisfile, ierr)

    end subroutine write_surfnode