source: LMDZ4/branches/V3_test/libf/phylmd/mod_phys_openmp.F90 @ 711

Last change on this file since 711 was 711, checked in by Laurent Fairhead, 18 years ago

Import initial YM
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 2.9 KB
Line 
1module mod_phys_openmp
2
3  integer,save :: omp_size
4  integer,save :: omp_rank
5 
6  INTEGER,SAVE,dimension(:),allocatable :: klon_omp_nb
7  INTEGER,SAVE,dimension(:),allocatable :: klon_omp_begin
8  INTEGER,SAVE,dimension(:),allocatable :: klon_omp_end   
9  INTEGER,SAVE :: klon_omp
10!$OMP  THREADPRIVATE(omp_rank,klon_omp)
11
12  REAL,save,allocatable,dimension(:) :: zmasq
13!$OMP THREADPRIVATE(zmasq)   
14
15contains
16 
17  subroutine init_phys_openmp
18  USE mod_phys_mpi, only : klon_mpi,kfdia,kidia,kdlon
19    implicit none
20    integer :: i
21#ifdef _OPENMP   
22    integer :: OMP_GET_NUM_THREADS
23    external OMP_GET_NUM_THREADS
24    integer :: OMP_GET_THREAD_NUM
25    external OMP_GET_THREAD_NUM
26#endif 
27
28#ifdef _OPENMP
29!$OMP MASTER
30        omp_size=OMP_GET_NUM_THREADS()
31!$OMP END MASTER
32        omp_rank=OMP_GET_THREAD_NUM()   
33#else   
34    omp_size=1
35    omp_rank=0
36#endif
37!$OMP MASTER
38     print *,'MASTER :omp_rank',omp_rank
39    allocate(klon_omp_nb(0:omp_size-1))
40    allocate(klon_omp_begin(0:omp_size-1))
41    allocate(klon_omp_end(0:omp_size-1))
42   
43    do i=0,omp_size-1
44      klon_omp_nb(i)=klon_mpi/omp_size
45      if (i<MOD(klon_mpi,omp_size)) klon_omp_nb(i)=klon_omp_nb(i)+1
46    enddo
47   
48    klon_omp_begin(0)=1
49    klon_omp_end(0)=klon_omp_nb(0)
50   
51    do i=1,omp_size-1
52      klon_omp_begin(i)=klon_omp_end(i-1)+1
53      klon_omp_end(i)=klon_omp_begin(i)+klon_omp_nb(i)-1
54    enddo
55!$OMP END MASTER
56!$OMP BARRIER
57   
58    klon_omp=klon_omp_nb(omp_rank)
59    allocate(zmasq(klon_omp)) 
60 
61    kidia=1
62    kfdia=klon_omp
63    kdlon=klon_omp
64    print *,omp_rank,' : klon_omp_nb',klon_omp_nb
65  end subroutine init_phys_openmp
66 
67
68  subroutine ScatterField_omp(Fields,Fieldr,ll)
69  USE mod_phys_mpi, only : klon_mpi
70  implicit none
71    INTEGER :: ll
72    REAL, dimension(klon_mpi,ll) :: Fields
73    REAL, dimension(klon_omp,ll) :: Fieldr     
74    REAL, dimension(:,:),ALLOCATABLE,SAVE :: Field_tmp     
75   
76    INTEGER :: i,l,offset
77   
78!$OMP BARRIER
79!$OMP MASTER
80    ALLOCATE(Field_tmp(klon_mpi,ll))
81
82    DO l=1,ll
83      DO i=1,klon_mpi
84        Field_tmp(i,l)=fields(i,l)
85      ENDDO
86    ENDDO
87   
88!$OMP END MASTER
89
90!$OMP BARRIER   
91    offset=klon_omp_begin(omp_rank)-1
92    do l=1,ll
93      do i=1,klon_omp
94        Fieldr(i,l)=Field_tmp(offset+i,l)
95      enddo
96    enddo
97!$OMP BARRIER
98
99!$OMP MASTER
100    DEALLOCATE(Field_tmp)
101!$OMP END MASTER
102
103  end subroutine ScatterField_omp
104
105  subroutine GatherField_omp(Fields,Fieldr,ll)
106  USE mod_phys_mpi, only : klon_mpi
107  implicit none
108    INTEGER :: ll
109    REAL, dimension(klon_omp,ll) :: Fields
110    REAL, dimension(klon_mpi,ll) :: Fieldr     
111    REAL, dimension(klon_omp,ll) :: Field_tmp
112   
113    INTEGER :: i,l,offset
114
115    Field_tmp(:,:)=Fields(:,:)
116!$OMP BARRIER     
117    offset=klon_omp_begin(omp_rank)-1
118    do l=1,ll
119      do i=1,klon_omp
120        Fieldr(offset+i,l)=Field_tmp(i,l)
121      enddo
122    enddo
123!$OMP BARRIER 
124  end subroutine GatherField_omp
125   
126end module mod_phys_openmp
Note: See TracBrowser for help on using the repository browser.