source: LMDZ4/branches/LMDZ4-dev-20091210/libf/dyn3dpar/bands.F90 @ 5066

Last change on this file since 5066 was 1222, checked in by Ehouarn Millour, 15 years ago

Changes and cleanups to enable compiling without physics
and without ioipsl.

IOIPSL related cleanups:

  • bibio/writehist.F encapsulate the routine (which needs IOIPSL to function)

with #ifdef IOIPSL flag.

  • dyn3d/abort_gcm.F, dyn3dpar/abort_gcm.F and dyn3dpar/getparam.F90: use ioipsl_getincom module when not compiling with IOIPSL library, in order to always be able to use getin() routine.
  • removed unused "use IOIPSL" in dyn3dpar/guide_p_mod.F90
  • calendar related issue: Initialize day_ref and annee_ref in iniacademic.F (i.e. when they are not read from start.nc file)

Earth-specific programs/routines/modules:
create_etat0.F, fluxstokenc.F, limit_netcdf.F, startvar.F
(versions in dyn3d and dyn3dpar)
These routines and modules, which by design and porpose are made to function with
Earth physics are encapsulated with #CPP_EARTH cpp flag.

Earth-specific instructions:

  • calls to qminimum (specific treatment of first 2 tracers, i.e. water) in dyn3d/caladvtrac.F, dyn3d/integrd.F, dyn3dpar/caladvtrac_p.F, dyn3dpar/integrd_p.F only if (planet_type == 'earth')

Interaction with parallel physics:

  • routine dyn3dpar/parallel.F90 uses "surface_data" module (which is in the physics ...) to know value of "type_ocean" . Encapsulated that with #ifdef CPP_EARTH and set to a default type_ocean="dummy" otherwise.
  • So far, only Earth physics are parallelized, so all the interaction between parallel dynamics and parallel physics are encapsulated with #ifdef CCP_EARTH (this way we can run parallel without any physics). The (dyn3dpar) routines which contains such interaction are: bands.F90, gr_dyn_fi_p.F, gr_fi_dyn_p.F, mod_interface_dyn_phys.F90 This should later (when improving dyn/phys interface) be encapsulated with a more general and appropriate #ifdef CPP_PHYS cpp flag.

I checked that these changes do not alter results (on the simple
32x24x11 bench) on Ciclad (seq & mpi), Brodie (seq, mpi & omp) and
Vargas (seq, mpi & omp).

EM

  • 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 1222 2009-08-07 11:48:33Z abarral $
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
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
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
95#ifdef CPP_EARTH
96! Ehouarn: what follows is only related to // physics; for now only for Earth
97     USE mod_phys_lmdz_para, ONLY : jj_para_begin,jj_para_end
98#endif
99     IMPLICIT NONE
100     INCLUDE 'dimensions.h'   
101     INTEGER :: i
102       
103      do i=0,mpi_size-1
104         jj_nb_vanleer2(i)=(jjm+1)/mpi_size
105         if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1
106      enddo
107         
108#ifdef CPP_EARTH
109! Ehouarn: what follows is only related to // physics; for now only for Earth         
110      do i=0,MPI_Size-1
111        jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1
112        if (i/=0) then
113          if (jj_para_begin(i)==jj_para_end(i-1)) then
114            jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1
115          endif
116        endif
117      enddo
118     
119      do i=0,MPI_Size-1
120        jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1
121        if (i/=0) then
122          if (jj_para_begin(i)==jj_para_end(i-1)) then
123            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
124          else
125            jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1
126            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
127          endif
128        endif
129      enddo
130#endif     
131     
132    end subroutine Set_Bands
133
134
135    subroutine AdjustBands_caldyn
136      use times
137      use parallel
138      implicit none
139
140      real :: minvalue,maxvalue
141      integer :: min_proc,max_proc
142      integer :: i,j
143      real,allocatable,dimension(:) :: value
144      integer,allocatable,dimension(:) :: index
145      real :: tmpvalue
146      integer :: tmpindex
147     
148      allocate(value(0:mpi_size-1))
149      allocate(index(0:mpi_size-1))
150       
151 
152      call allgather_timer_average
153
154      do i=0,mpi_size-1
155        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
156        index(i)=i
157      enddo
158     
159      do i=0,mpi_size-2
160        do j=i+1,mpi_size-1
161          if (value(i)>value(j)) then
162            tmpvalue=value(i)
163            value(i)=value(j)
164            value(j)=tmpvalue
165           
166            tmpindex=index(i)
167            index(i)=index(j)
168            index(j)=tmpindex
169           endif
170         enddo
171      enddo
172     
173      maxvalue=value(mpi_size-1)
174      max_proc=index(mpi_size-1)           
175           
176      do i=0,mpi_size-2
177        minvalue=value(i)
178        min_proc=index(i)
179        if (jj_nb_caldyn(max_proc)>3) then
180          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
181             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
182             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
183             exit
184           else
185             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
186                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) 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             endif
191           endif
192         endif
193      enddo
194     
195      deallocate(value)
196      deallocate(index)
197         
198    end subroutine AdjustBands_caldyn
199   
200    subroutine AdjustBands_vanleer
201      use times
202      use parallel
203      implicit none
204
205      real :: minvalue,maxvalue
206      integer :: min_proc,max_proc
207      integer :: i,j
208      real,allocatable,dimension(:) :: value
209      integer,allocatable,dimension(:) :: index
210      real :: tmpvalue
211      integer :: tmpindex
212     
213      allocate(value(0:mpi_size-1))
214      allocate(index(0:mpi_size-1))
215       
216 
217      call allgather_timer_average
218
219      do i=0,mpi_size-1
220        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
221        index(i)=i
222      enddo
223     
224      do i=0,mpi_size-2
225        do j=i+1,mpi_size-1
226          if (value(i)>value(j)) then
227            tmpvalue=value(i)
228            value(i)=value(j)
229            value(j)=tmpvalue
230           
231            tmpindex=index(i)
232            index(i)=index(j)
233            index(j)=tmpindex
234           endif
235         enddo
236      enddo
237     
238      maxvalue=value(mpi_size-1)
239      max_proc=index(mpi_size-1)           
240           
241      do i=0,mpi_size-2
242        minvalue=value(i)
243        min_proc=index(i)
244
245        if (jj_nb_vanleer(max_proc)>3) then
246          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
247             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
248             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
249             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
250             exit
251           else
252             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
253               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
254               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
255               exit
256             endif
257           endif
258         endif
259      enddo
260     
261      deallocate(value)
262      deallocate(index)
263         
264    end subroutine AdjustBands_vanleer
265
266    subroutine AdjustBands_dissip
267      use times
268      use parallel
269      implicit none
270
271      real :: minvalue,maxvalue
272      integer :: min_proc,max_proc
273      integer :: i,j
274      real,allocatable,dimension(:) :: value
275      integer,allocatable,dimension(:) :: index
276      real :: tmpvalue
277      integer :: tmpindex
278     
279      allocate(value(0:mpi_size-1))
280      allocate(index(0:mpi_size-1))
281       
282 
283      call allgather_timer_average
284
285      do i=0,mpi_size-1
286        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
287        index(i)=i
288      enddo
289     
290      do i=0,mpi_size-2
291        do j=i+1,mpi_size-1
292          if (value(i)>value(j)) then
293            tmpvalue=value(i)
294            value(i)=value(j)
295            value(j)=tmpvalue
296           
297            tmpindex=index(i)
298            index(i)=index(j)
299            index(j)=tmpindex
300           endif
301         enddo
302      enddo
303     
304      maxvalue=value(mpi_size-1)
305      max_proc=index(mpi_size-1)           
306           
307      do i=0,mpi_size-2
308        minvalue=value(i)
309        min_proc=index(i)
310
311        if (jj_nb_dissip(max_proc)>3) then
312          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
313             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
314             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
315             exit
316           else
317             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
318                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) 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             endif
323           endif
324         endif
325      enddo
326     
327      deallocate(value)
328      deallocate(index)
329         
330    end subroutine AdjustBands_dissip
331
332    subroutine AdjustBands_physic
333      use times
334#ifdef CPP_EARTH
335! Ehouarn: what follows is only related to // physics; for now only for Earth
336      USE mod_phys_lmdz_para, only : klon_mpi_para_nb
337#endif
338      USE parallel
339      implicit none
340
341      integer :: i,Index
342      real,allocatable,dimension(:) :: value
343      integer,allocatable,dimension(:) :: Inc
344      real :: medium
345      integer :: NbTot,sgn
346     
347      allocate(value(0:mpi_size-1))
348      allocate(Inc(0:mpi_size-1))
349       
350 
351      call allgather_timer_average
352     
353      medium=0
354      do i=0,mpi_size-1
355        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
356        medium=medium+value(i)
357      enddo   
358     
359      medium=medium/mpi_size     
360      NbTot=0
361#ifdef CPP_EARTH
362! Ehouarn: what follows is only related to // physics; for now only for Earth
363      do i=0,mpi_size-1
364        Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i))
365        NbTot=NbTot+Inc(i) 
366      enddo
367     
368      if (NbTot>=0) then
369        Sgn=1
370      else
371        Sgn=-1
372        NbTot=-NbTot
373      endif
374     
375      Index=0
376      do i=1,NbTot
377        Inc(Index)=Inc(Index)-Sgn
378        Index=Index+1
379        if (Index>mpi_size-1) Index=0
380      enddo
381     
382      do i=0,mpi_size-1
383        distrib_phys(i)=klon_mpi_para_nb(i)+inc(i)
384      enddo
385#endif     
386    end subroutine AdjustBands_physic
387
388    subroutine WriteBands
389    USE parallel
390    implicit none
391    include "dimensions.h"
392
393      integer :: i,j
394      character (len=4) :: siim,sjjm,sllm,sproc
395      character (len=255) :: filename
396      integer :: unit_number=10
397      integer :: ierr
398 
399      write(siim,'(i3)') iim
400      write(sjjm,'(i3)') jjm
401      write(sllm,'(i3)') llm
402      write(sproc,'(i3)') mpi_size
403
404      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
405                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
406     
407      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
408     
409      if (ierr==0) then
410       
411!       write (unit_number,*) '*** Bandes caldyn ***'
412        do i=0,mpi_size-1
413          write (unit_number,*) i,jj_nb_caldyn(i)
414        enddo
415       
416!       write (unit_number,*) '*** Bandes vanleer ***'
417        do i=0,mpi_size-1
418          write (unit_number,*) i,jj_nb_vanleer(i)
419        enddo
420       
421!        write (unit_number,*) '*** Bandes dissip ***'
422        do i=0,mpi_size-1
423          write (unit_number,*) i,jj_nb_dissip(i)
424        enddo
425       
426        do i=0,mpi_size-1
427          write (unit_number,*) i,distrib_phys(i)
428        enddo
429       
430        CLOSE(unit_number)   
431      else
432        print *,'probleme lors de l ecriture des bandes'
433      endif
434       
435    end subroutine WriteBands
436 
437  end module Bands
438 
439 
Note: See TracBrowser for help on using the repository browser.