~fluidity-core/fluidity/embedded_models

« back to all changes in this revision

Viewing changes to libadaptivity/adapt3d/src/AdaptProgress.F90

  • Committer: Timothy Bond
  • Date: 2011-04-14 15:40:24 UTC
  • Revision ID: timothy.bond@imperial.ac.uk-20110414154024-116ci9gq6mwigmaw
Following the move from svn to bzr we change the nature of inclusion of these
four software libraries. Previously, they were included as svn externals and
pulled at latest version for trunk, pinned to specific versions for release
and stable trunk. Since bzr is less elegant at dealing with externals we have
made the decision to include the packages directly into the trunk instead.

At this import the versions are:

libadaptivity: r163
libvtkfortran: r67
libspud: r545
libmba2d: r28

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
! Copyright (C) 2009 Imperial College London and others.
 
2
 
3
! Please see the AUTHORS file in the main source directory for a full list
 
4
! of copyright holders.
 
5
 
6
! Gerard Gorman
 
7
! Applied Modelling and Computation Group
 
8
! Department of Earth Science and Engineering
 
9
! Imperial College London
 
10
 
11
! g.gorman@imperial.ac.uk
 
12
 
13
! This library is free software; you can redistribute it and/or
 
14
! modify it under the terms of the GNU Lesser General Public
 
15
! License as published by the Free Software Foundation; either
 
16
! version 2.1 of the License.
 
17
 
18
! This library is distributed in the hope that it will be useful,
 
19
! but WITHOUT ANY WARRANTY; without even the implied warranty of
 
20
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
 
21
! Lesser General Public License for more details.
 
22
 
23
! You should have received a copy of the GNU Lesser General Public
 
24
! License along with this library; if not, write to the Free Software
 
25
! Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307
 
26
! USA
 
27
 
 
28
#include "confdefs.h"
 
29
 
 
30
module AdaptProgress
 
31
  implicit none
 
32
  
 
33
  private
 
34
  
 
35
  public::initialise, finalize, should_exit
 
36
#ifdef HAVE_MPI
 
37
  include 'mpif.h'
 
38
#endif
 
39
  logical::initialised=.false.
 
40
  integer, dimension(:), allocatable::load, rrequest
 
41
  integer::myrank, nprocs
 
42
  real::imbalance_tol=0.5
 
43
  integer win
 
44
 
 
45
contains
 
46
  subroutine initialise(count, tol)
 
47
    integer, intent(in)::count
 
48
    real, intent(in)::tol
 
49
#ifdef HAVE_MPI
 
50
    integer have_mpi_init, i, ierr
 
51
    if(.not.initialised) then
 
52
       call MPI_Initialized(have_mpi_init, ierr)
 
53
       if(have_mpi_init.eq.0) then
 
54
          nprocs=1
 
55
       else
 
56
          call MPI_Comm_size(MPI_COMM_WORLD, nprocs, ierr)
 
57
       endif
 
58
 
 
59
       if(nprocs.gt.1) then
 
60
          call MPI_Comm_rank(MPI_COMM_WORLD, myrank, ierr)
 
61
          allocate(load(0:nprocs-1))
 
62
          call MPI_Allgather(count, 1, MPI_INTEGER, load, 1, &
 
63
               MPI_INTEGER, MPI_COMM_WORLD, ierr) 
 
64
          
 
65
          allocate(rrequest(0:nprocs-1))
 
66
          
 
67
          do i=0, nprocs-1
 
68
             if(i.ne.myrank) then
 
69
                call MPI_Irecv(load(i), 1, MPI_INTEGER, i, 1, &
 
70
                     MPI_COMM_WORLD, rrequest(i), ierr)
 
71
             else
 
72
                rrequest(i) = MPI_REQUEST_NULL
 
73
             end if
 
74
          end do
 
75
 
 
76
          imbalance_tol = tol
 
77
       end if
 
78
 
 
79
       initialised = .true.
 
80
    end if
 
81
#endif
 
82
  end subroutine initialise
 
83
 
 
84
  subroutine finalize(count)
 
85
    integer, intent(in)::count
 
86
#ifdef HAVE_MPI
 
87
    integer i, ierr
 
88
    integer, allocatable, dimension(:)::request
 
89
    integer, allocatable, dimension(:, :)::status
 
90
 
 
91
    if(nprocs.gt.1) then
 
92
       allocate(request(0:nprocs-1))
 
93
       allocate(status(MPI_STATUS_SIZE, 0:nprocs-1))
 
94
       
 
95
       do i=0, nprocs-1
 
96
          if(i.ne.myrank) then
 
97
             call MPI_Isend(count, 1, MPI_INTEGER, i, 1, &
 
98
                  MPI_COMM_WORLD, request(i), ierr)
 
99
          else
 
100
             request(i) = MPI_REQUEST_NULL
 
101
          end if
 
102
       end do
 
103
       
 
104
       call MPI_Waitall(nprocs, rrequest, status, ierr)
 
105
       call MPI_Waitall(nprocs, request, status, ierr)
 
106
       
 
107
       deallocate(request, rrequest, status, load)
 
108
    end if
 
109
    initialised = .false.
 
110
#endif
 
111
  end subroutine finalize
 
112
  
 
113
  logical function should_exit(count)
 
114
    integer, intent(in)::count
 
115
#ifdef HAVE_MPI
 
116
    real imbalance
 
117
    integer ierr
 
118
    
 
119
    if(nprocs.gt.1) then
 
120
       load(myrank) = count
 
121
       
 
122
       imbalance = 1.0 - real(sum(load))/(nprocs*maxval(load))
 
123
       
 
124
       should_exit = (imbalance>imbalance_tol)
 
125
    else
 
126
       should_exit = .false.
 
127
    end if
 
128
#else
 
129
    should_exit = .false.
 
130
#endif
 
131
  end function should_exit
 
132
  
 
133
end module AdaptProgress