source: LMDZ6/branches/DYNAMICO-conv/libf/dyn3dpar/bands.F90 @ 3024

Last change on this file since 3024 was 2351, checked in by Ehouarn Millour, 9 years ago

More on physics/dynamics separation and cleanup:

  • Set things up so that all physics-related initializations are done via iniphysiq.
  • Created a "geometry_mod.F90" module in phy_common to store information on the loacl grid (i.e. replaces comgeomphy) and moreover give these variables more obvious names (e.g.: rlond => longitude, rlatd => latitude, airephy => cell_area).
  • removed obsolete comgeomphy.h and comgeomphy.F90

EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 11.5 KB
Line 
1!
2! $Id: bands.F90 2351 2015-08-25 15:14:59Z fairhead $
3!
4  module Bands
5 
6    integer, parameter :: bands_caldyn=1
7    integer, parameter :: bands_vanleer=2
8    integer, parameter :: bands_dissip=3
9   
10    INTEGER,dimension(:),allocatable :: jj_Nb_Caldyn
11    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer
12    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer2
13    INTEGER,dimension(:),allocatable :: jj_Nb_dissip
14    INTEGER,dimension(:),allocatable :: jj_Nb_physic
15    INTEGER,dimension(:),allocatable :: jj_Nb_physic_bis
16    INTEGER,dimension(:),allocatable :: distrib_phys
17 
18  contains
19 
20  subroutine AllocateBands
21    USE parallel_lmdz
22    implicit none
23   
24    allocate(jj_Nb_Caldyn(0:MPI_Size-1))
25    allocate(jj_Nb_vanleer(0:MPI_Size-1))
26    allocate(jj_Nb_vanleer2(0:MPI_Size-1))
27    allocate(jj_Nb_dissip(0:MPI_Size-1))
28    allocate(jj_Nb_physic(0:MPI_Size-1))
29    allocate(jj_Nb_physic_bis(0:MPI_Size-1))
30    allocate(distrib_phys(0:MPI_Size-1))
31 
32  end subroutine AllocateBands
33 
34  subroutine Read_distrib
35    USE parallel_lmdz
36    implicit none
37
38    include "dimensions.h"
39      integer :: i,j
40      character (len=4) :: siim,sjjm,sllm,sproc
41      character (len=255) :: filename
42      integer :: unit_number=10
43      integer :: ierr
44   
45      call AllocateBands
46      write(siim,'(i3)') iim
47      write(sjjm,'(i3)') jjm
48      write(sllm,'(i3)') llm
49      write(sproc,'(i3)') mpi_size
50      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
51                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
52       
53      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='old',FORM='formatted',IOSTAT=ierr)
54     
55      if (ierr==0) then
56     
57         do i=0,mpi_size-1
58          read (unit_number,*) j,jj_nb_caldyn(i)
59        enddo
60     
61        do i=0,mpi_size-1
62          read (unit_number,*) j,jj_nb_vanleer(i)
63        enddo
64     
65        do i=0,mpi_size-1
66          read (unit_number,*) j,jj_nb_dissip(i)
67        enddo
68     
69        do i=0,mpi_size-1
70          read (unit_number,*) j,distrib_phys(i)
71        enddo
72       
73        CLOSE(unit_number) 
74 
75      else
76        do i=0,mpi_size-1
77          jj_nb_caldyn(i)=(jjm+1)/mpi_size
78          if (i<MOD(jjm+1,mpi_size)) jj_nb_caldyn(i)=jj_nb_caldyn(i)+1
79        enddo
80     
81        jj_nb_vanleer(:)=jj_nb_caldyn(:)
82        jj_nb_dissip(:)=jj_nb_caldyn(:)
83       
84        do i=0,mpi_size-1
85          distrib_phys(i)=(iim*(jjm-1)+2)/mpi_size
86          IF (i<MOD(iim*(jjm-1)+2,mpi_size)) distrib_phys(i)=distrib_phys(i)+1
87        enddo
88      endif
89 
90   end subroutine Read_distrib
91   
92   
93   SUBROUTINE  Set_Bands
94     USE parallel_lmdz
95     IMPLICIT NONE
96     INCLUDE 'dimensions.h'   
97     INTEGER :: i, ij
98     INTEGER :: jj_para_begin(0:mpi_size-1)
99     INTEGER :: jj_para_end(0:mpi_size-1)
100       
101      do i=0,mpi_size-1
102         jj_nb_vanleer2(i)=(jjm+1)/mpi_size
103         if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1
104      enddo
105         
106      jj_para_begin(0)=1
107      ij=distrib_phys(0)+iim-1
108      jj_para_end(0)=((ij-1)/iim)+1
109     
110      DO i=1,mpi_Size-1
111        ij=ij+1
112        jj_para_begin(i)=((ij-1)/iim)+1
113        ij=ij+distrib_phys(i)-1
114        jj_para_end(i)=((ij-1)/iim)+1
115      ENDDO
116 
117      do i=0,MPI_Size-1
118        jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1
119        if (i/=0) then
120          if (jj_para_begin(i)==jj_para_end(i-1)) then
121            jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1
122          endif
123        endif
124      enddo
125     
126      do i=0,MPI_Size-1
127        jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1
128        if (i/=0) then
129          if (jj_para_begin(i)==jj_para_end(i-1)) then
130            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
131          else
132            jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1
133            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
134          endif
135        endif
136      enddo
137     
138    end subroutine Set_Bands
139
140
141    subroutine AdjustBands_caldyn
142      use times
143      USE parallel_lmdz
144      implicit none
145
146      real :: minvalue,maxvalue
147      integer :: min_proc,max_proc
148      integer :: i,j
149      real,allocatable,dimension(:) :: value
150      integer,allocatable,dimension(:) :: index
151      real :: tmpvalue
152      integer :: tmpindex
153     
154      allocate(value(0:mpi_size-1))
155      allocate(index(0:mpi_size-1))
156       
157 
158      call allgather_timer_average
159
160      do i=0,mpi_size-1
161        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
162        index(i)=i
163      enddo
164     
165      do i=0,mpi_size-2
166        do j=i+1,mpi_size-1
167          if (value(i)>value(j)) then
168            tmpvalue=value(i)
169            value(i)=value(j)
170            value(j)=tmpvalue
171           
172            tmpindex=index(i)
173            index(i)=index(j)
174            index(j)=tmpindex
175           endif
176         enddo
177      enddo
178     
179      maxvalue=value(mpi_size-1)
180      max_proc=index(mpi_size-1)           
181           
182      do i=0,mpi_size-2
183        minvalue=value(i)
184        min_proc=index(i)
185        if (jj_nb_caldyn(max_proc)>3) then
186          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
187             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
188             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
189             exit
190           else
191             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
192                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
193               jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
194               jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
195               exit
196             endif
197           endif
198         endif
199      enddo
200     
201      deallocate(value)
202      deallocate(index)
203         
204    end subroutine AdjustBands_caldyn
205   
206    subroutine AdjustBands_vanleer
207      use times
208      USE parallel_lmdz
209      implicit none
210
211      real :: minvalue,maxvalue
212      integer :: min_proc,max_proc
213      integer :: i,j
214      real,allocatable,dimension(:) :: value
215      integer,allocatable,dimension(:) :: index
216      real :: tmpvalue
217      integer :: tmpindex
218     
219      allocate(value(0:mpi_size-1))
220      allocate(index(0:mpi_size-1))
221       
222 
223      call allgather_timer_average
224
225      do i=0,mpi_size-1
226        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
227        index(i)=i
228      enddo
229     
230      do i=0,mpi_size-2
231        do j=i+1,mpi_size-1
232          if (value(i)>value(j)) then
233            tmpvalue=value(i)
234            value(i)=value(j)
235            value(j)=tmpvalue
236           
237            tmpindex=index(i)
238            index(i)=index(j)
239            index(j)=tmpindex
240           endif
241         enddo
242      enddo
243     
244      maxvalue=value(mpi_size-1)
245      max_proc=index(mpi_size-1)           
246           
247      do i=0,mpi_size-2
248        minvalue=value(i)
249        min_proc=index(i)
250
251        if (jj_nb_vanleer(max_proc)>3) then
252          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
253             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
254             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
255             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
256             exit
257           else
258             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
259               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
260               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
261               exit
262             endif
263           endif
264         endif
265      enddo
266     
267      deallocate(value)
268      deallocate(index)
269         
270    end subroutine AdjustBands_vanleer
271
272    subroutine AdjustBands_dissip
273      use times
274      USE parallel_lmdz
275      implicit none
276
277      real :: minvalue,maxvalue
278      integer :: min_proc,max_proc
279      integer :: i,j
280      real,allocatable,dimension(:) :: value
281      integer,allocatable,dimension(:) :: index
282      real :: tmpvalue
283      integer :: tmpindex
284     
285      allocate(value(0:mpi_size-1))
286      allocate(index(0:mpi_size-1))
287       
288 
289      call allgather_timer_average
290
291      do i=0,mpi_size-1
292        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
293        index(i)=i
294      enddo
295     
296      do i=0,mpi_size-2
297        do j=i+1,mpi_size-1
298          if (value(i)>value(j)) then
299            tmpvalue=value(i)
300            value(i)=value(j)
301            value(j)=tmpvalue
302           
303            tmpindex=index(i)
304            index(i)=index(j)
305            index(j)=tmpindex
306           endif
307         enddo
308      enddo
309     
310      maxvalue=value(mpi_size-1)
311      max_proc=index(mpi_size-1)           
312           
313      do i=0,mpi_size-2
314        minvalue=value(i)
315        min_proc=index(i)
316
317        if (jj_nb_dissip(max_proc)>3) then
318          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
319             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
320             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
321             exit
322           else
323             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
324                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
325               jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
326               jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
327               exit
328             endif
329           endif
330         endif
331      enddo
332     
333      deallocate(value)
334      deallocate(index)
335         
336    end subroutine AdjustBands_dissip
337
338    subroutine AdjustBands_physic
339      use times
340#ifdef CPP_PHYS
341! Ehouarn: what follows is only related to // physics
342      USE mod_phys_lmdz_para, only : klon_mpi_para_nb
343#endif
344      USE parallel_lmdz
345      implicit none
346
347      integer :: i,Index
348      real,allocatable,dimension(:) :: value
349      integer,allocatable,dimension(:) :: Inc
350      real :: medium
351      integer :: NbTot,sgn
352     
353      allocate(value(0:mpi_size-1))
354      allocate(Inc(0:mpi_size-1))
355       
356 
357      call allgather_timer_average
358     
359      medium=0
360      do i=0,mpi_size-1
361        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
362        medium=medium+value(i)
363      enddo   
364     
365      medium=medium/mpi_size     
366      NbTot=0
367#ifdef CPP_PHYS
368      do i=0,mpi_size-1
369        Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i))
370        NbTot=NbTot+Inc(i) 
371      enddo
372     
373      if (NbTot>=0) then
374        Sgn=1
375      else
376        Sgn=-1
377        NbTot=-NbTot
378      endif
379     
380      Index=0
381      do i=1,NbTot
382        Inc(Index)=Inc(Index)-Sgn
383        Index=Index+1
384        if (Index>mpi_size-1) Index=0
385      enddo
386     
387      do i=0,mpi_size-1
388        distrib_phys(i)=klon_mpi_para_nb(i)+inc(i)
389      enddo
390#endif     
391    end subroutine AdjustBands_physic
392
393    subroutine WriteBands
394    USE parallel_lmdz
395    implicit none
396    include "dimensions.h"
397
398      integer :: i,j
399      character (len=4) :: siim,sjjm,sllm,sproc
400      character (len=255) :: filename
401      integer :: unit_number=10
402      integer :: ierr
403 
404      write(siim,'(i3)') iim
405      write(sjjm,'(i3)') jjm
406      write(sllm,'(i3)') llm
407      write(sproc,'(i3)') mpi_size
408
409      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
410                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
411     
412      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
413     
414      if (ierr==0) then
415       
416!       write (unit_number,*) '*** Bandes caldyn ***'
417        do i=0,mpi_size-1
418          write (unit_number,*) i,jj_nb_caldyn(i)
419        enddo
420       
421!       write (unit_number,*) '*** Bandes vanleer ***'
422        do i=0,mpi_size-1
423          write (unit_number,*) i,jj_nb_vanleer(i)
424        enddo
425       
426!        write (unit_number,*) '*** Bandes dissip ***'
427        do i=0,mpi_size-1
428          write (unit_number,*) i,jj_nb_dissip(i)
429        enddo
430       
431        do i=0,mpi_size-1
432          write (unit_number,*) i,distrib_phys(i)
433        enddo
434       
435        CLOSE(unit_number)   
436      else
437        print *,'probleme lors de l ecriture des bandes'
438      endif
439       
440    end subroutine WriteBands
441 
442  end module Bands
443 
444 
Note: See TracBrowser for help on using the repository browser.