source: LMDZ5/branches/testing/libf/dyn3dmem/write_field_loc.F90 @ 5423

Last change on this file since 5423 was 1910, checked in by Laurent Fairhead, 11 years ago

Merged trunk changes r1860:1909 into testing branch

  • 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
File size: 3.5 KB
Line 
1module write_field_loc
2implicit none
3 
4  interface WriteField_u
5    module procedure Write_field1d_u,Write_Field2d_u
6  end interface WriteField_u
7
8  interface WriteField_v
9    module procedure Write_field1d_v,Write_Field2d_v
10  end interface WriteField_v
11 
12  contains
13 
14  subroutine write_field1D_u(name,Field)
15    character(len=*)   :: name
16    real, dimension(:) :: Field
17
18    CALL write_field_u_gen(name,Field,1)
19
20  end subroutine write_field1D_u
21
22  subroutine write_field2D_u(name,Field)
23    implicit none
24     
25    character(len=*)   :: name
26    real, dimension(:,:) :: Field
27    integer :: ll
28   
29    ll=size(field,2)   
30    CALL write_field_u_gen(name,Field,ll)
31   
32    end subroutine write_field2D_u
33
34
35   SUBROUTINE write_field_u_gen(name,Field,ll)
36    USE parallel_lmdz
37    USE write_field
38    USE mod_hallo
39    implicit none
40    include 'dimensions.h'
41    include 'paramet.h'
42     
43    character(len=*)   :: name
44    real, dimension(ijb_u:ije_u,ll) :: Field
45    real, allocatable,SAVE :: New_Field(:,:,:)
46    integer,dimension(0:mpi_size-1) :: jj_nb_master
47    type(Request),SAVE :: Request_write
48!$OMP THREADPRIVATE(Request_write)
49    integer :: ll,i
50   
51   
52    jj_nb_master(:)=0
53    jj_nb_master(0)=jjp1
54!$OMP BARRIER
55!$OMP MASTER
56    allocate(New_Field(iip1,jjp1,ll))
57!$OMP END MASTER
58!$OMP BARRIER
59
60!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
61    DO i=1,ll   
62      New_Field(:,jj_begin:jj_end,i)=reshape(Field(ij_begin:ij_end,i),(/iip1,jj_nb/))
63    ENDDO
64!$OMP BARRIER   
65    call Register_SwapField(new_field,new_field,ip1jmp1,ll,jj_Nb_master,Request_write)
66    call SendRequest(Request_write)
67!$OMP BARRIER
68    call WaitRequest(Request_write)     
69!$OMP BARRIER
70
71!$OMP MASTER
72    if (MPI_Rank==0) call WriteField(name,New_Field)
73    DEALLOCATE(New_Field)
74!$OMP END MASTER       
75!$OMP BARRIER
76    END SUBROUTINE write_field_u_gen
77
78
79  subroutine write_field1D_v(name,Field)
80    character(len=*)   :: name
81    real, dimension(:) :: Field
82
83    CALL write_field_v_gen(name,Field,1)
84
85  end subroutine write_field1D_v
86
87  subroutine write_field2D_v(name,Field)
88    implicit none
89     
90    character(len=*)   :: name
91    real, dimension(:,:) :: Field
92    integer :: ll
93   
94    ll=size(field,2)   
95    CALL write_field_v_gen(name,Field,ll)
96   
97    end subroutine write_field2D_v
98
99
100   SUBROUTINE write_field_v_gen(name,Field,ll)
101    USE parallel_lmdz
102    USE write_field
103    USE mod_hallo
104    implicit none
105    include 'dimensions.h'
106    include 'paramet.h'
107     
108    character(len=*)   :: name
109    real, dimension(ijb_v:ije_v,ll) :: Field
110    real, allocatable,SAVE :: New_Field(:,:,:)
111    integer,dimension(0:mpi_size-1) :: jj_nb_master
112    type(Request),SAVE :: Request_write
113!$OMP THREADPRIVATE(Request_write)   
114    integer :: ll,i,jje,ije,jjn
115   
116   
117    jj_nb_master(:)=0
118    jj_nb_master(0)=jjp1
119
120!$OMP BARRIER
121!$OMP MASTER
122    allocate(New_Field(iip1,jjm,ll))
123!$OMP END MASTER
124!$OMP BARRIER
125
126   IF (pole_sud) THEN
127     jje=jj_end-1
128     ije=ij_end-iip1
129     jjn=jj_nb-1
130   ELSE
131     jje=jj_end
132     ije=ij_end
133     jjn=jj_nb
134   ENDIF
135   
136!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
137    DO i=1,ll   
138      New_Field(:,jj_begin:jje,i)=reshape(Field(ij_begin:ije,i),(/iip1,jjn/))
139    ENDDO
140!$OMP BARRIER   
141    call Register_SwapField(new_field,new_field,ip1jm,ll,jj_Nb_master,Request_write)
142    call SendRequest(Request_write)
143!$OMP BARRIER
144    call WaitRequest(Request_write)     
145!$OMP BARRIER
146
147!$OMP MASTER
148    if (MPI_Rank==0) call WriteField(name,New_Field)
149    DEALLOCATE(New_Field)
150!$OMP END MASTER       
151!$OMP BARRIER
152    END SUBROUTINE write_field_v_gen
153   
154end module write_field_loc
155 
Note: See TracBrowser for help on using the repository browser.