ObserverUpdateTdelay Subroutine

public subroutine ObserverUpdateTdelay(self, totalNumberOfFaces)

Uses

    • MPI_Process_Info

Type Bound

ObserverClass

Arguments

Type IntentOptional Attributes Name
class(ObserverClass) :: self
integer, intent(in) :: totalNumberOfFaces

Source Code

   Subroutine ObserverUpdateTdelay(self, totalNumberOfFaces)

!     *******************************************************************
!        This subroutine updates the observer time delay. For static surfaces it 
!        doesn't need to be updated at every iteration.
!        Minimum time is used, for interpolation procedures
!     *******************************************************************
!
      use MPI_Process_Info
      implicit none
      class   (ObserverClass)                              :: self
      integer, intent(in)                                  :: totalNumberOfFaces
      
      ! local variables
      integer                                              :: i, ierr
      real(kind=RP)                                        :: t, tmax
      real(kind=RP), dimension(:), allocatable             :: alltDelay
      real(kind=RP), dimension(self % numberOfFaces)       :: tDelayArray
      integer, dimension(MPI_Process % nProcs)             :: no_of_faces_p, displs

      do i =1, self % numberOfFaces
        tDelayArray(i) = self % sourcePair(i) % tDelay
      end do

      if ( (MPI_Process % doMPIAction) ) then
#ifdef _HAS_MPI_
          call mpi_gather(self % numberOfFaces,1,MPI_INT,no_of_faces_p,1,MPI_INT,0,MPI_COMM_WORLD,ierr)

          if (MPI_Process % isRoot) then
              displs=0
              do i = 2, MPI_Process % nProcs 
                  displs(i) = displs(i-1) + no_of_faces_p(i-1)
              end do
          end if

          allocate(alltDelay(totalNumberOfFaces))
          call mpi_gatherv(tDelayArray, self % numberOfFaces, MPI_DOUBLE, &
                                     alltDelay, no_of_faces_p, displs, MPI_DOUBLE, 0, MPI_COMM_WORLD, ierr)

          if (MPI_Process % isRoot) then
              t = minval(alltDelay)
              tmax = maxval(alltDelay)
          end if

          call mpi_Bcast(t, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD, ierr)
          call mpi_Bcast(tmax, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD, ierr)
#endif
      else
          t = minval(tDelayArray)
          tmax = maxval(tDelayArray)
      end if

      self % tDelay = t
      self % tDelayMax = tmax

   End Subroutine ObserverUpdateTdelay