FV3 Bundle
m_distribution.f90
Go to the documentation of this file.
2  use fckit_module
3  use fckit_log_module, only: fckit_log
4 
5  implicit none
6 
7  private
8  public:: random_distribution
9 
11  type(fckit_mpi_comm),private :: comm !> parallel communicator
12  integer, private :: nobs !> the counter of the number of jobs
13  integer, private :: myproc !> size of parallel communicator
14  integer, private :: nproc !> my rank
15  integer, allocatable, public :: indx(:) !> distribution layout on each PE
16  contains
17  procedure :: received
18  procedure :: reset
19  procedure :: nobs_pe
20  procedure :: mype
21  final :: destructor
22  end type
23 
24  interface random_distribution
25  module procedure constructor
26  end interface
27 
28  contains
29 
30  type(random_distribution) function constructor(gnobs)
31  integer, intent(in) :: gnobs
32 
33  ! Local variables
34  integer :: info, i, ic
35  logical :: init
36 
37  call fckit_log%debug('random_distribution object created')
38  constructor%nobs=0
39 
40  ! Get the nproc and myproc
41  constructor%nproc=constructor%comm%size()
42  constructor%myproc=constructor%comm%rank()
43 
44  ! Count the obs. number on each PE
45  do i=1,gnobs
46  if (mod(i-1, constructor%nproc) .eq. constructor%myproc) then
47  constructor%nobs=constructor%nobs+1
48  end if
49  end do
50 
51  allocate(constructor%indx(constructor%nobs))
52 
53  ! Create a default layout
54  ic=0
55  do i=1,gnobs
56  if (mod(i-1, constructor%nproc) .eq. constructor%myproc) then
57  ic=ic+1
58  constructor%indx(ic)=i
59  end if
60  end do
61  end function
62 
63  subroutine destructor(this)
64  type(random_distribution), intent(inout):: this
65  if (allocated (this%indx)) then
66  deallocate(this%indx)
67  end if
68  call fckit_log%debug('random_distribution object finalized')
69  end subroutine
70 
71  logical function received(this, seqno)
72  class(random_distribution), intent(inout):: this
73  integer, intent(in) :: seqno
74 
75  received=.false.
76  !> seqno should start from 1
77  if (mod(seqno-1, this%nproc) .eq. this%myproc) then
78  this%nobs=this%nobs+1
79  received=.true.
80  end if
81  call fckit_log%debug('distribute Obs Data Functionality')
82  end function
83 
84  integer function nobs_pe(this)
85  class(random_distribution), intent(in):: this
86 
87  call fckit_log%debug('nobs_pe Data Functionality')
88  nobs_pe=this%nobs
89  end function
90 
91  integer function mype(this)
92  class(random_distribution), intent(in):: this
93 
94  mype=this%myproc
95  end function
96 
97  subroutine reset(this)
98  class(random_distribution), intent(inout):: this
99 
100  call fckit_log%debug('reset Data Functionality')
101  this%nobs=0
102  end subroutine
103 
104 end module type_distribution
integer function mype(this)
integer function nobs_pe(this)
subroutine destructor(this)
logical function received(this, seqno)
subroutine reset(this)
type(random_distribution) function constructor(gnobs)