source: LMDZ6/trunk/libf/dyn3dmem/write_field_loc.f90 @ 5456

Last change on this file since 5456 was 5285, checked in by abarral, 2 months ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

  • 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.6 KB
RevLine 
[1632]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)
[1823]36    USE parallel_lmdz
[1632]37    USE write_field
38    USE mod_hallo
[5271]39    USE dimensions_mod, ONLY: iim, jjm, llm, ndm
[5285]40USE paramet_mod_h
[5271]41implicit none
42
[5272]43
[1632]44     
45    character(len=*)   :: name
46    real, dimension(ijb_u:ije_u,ll) :: Field
47    real, allocatable,SAVE :: New_Field(:,:,:)
48    integer,dimension(0:mpi_size-1) :: jj_nb_master
[1848]49    type(Request),SAVE :: Request_write
50!$OMP THREADPRIVATE(Request_write)
[1632]51    integer :: ll,i
52   
53   
54    jj_nb_master(:)=0
55    jj_nb_master(0)=jjp1
56!$OMP BARRIER
57!$OMP MASTER
58    allocate(New_Field(iip1,jjp1,ll))
59!$OMP END MASTER
60!$OMP BARRIER
61
62!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
63    DO i=1,ll   
64      New_Field(:,jj_begin:jj_end,i)=reshape(Field(ij_begin:ij_end,i),(/iip1,jj_nb/))
65    ENDDO
[1848]66!$OMP BARRIER   
[1632]67    call Register_SwapField(new_field,new_field,ip1jmp1,ll,jj_Nb_master,Request_write)
68    call SendRequest(Request_write)
69!$OMP BARRIER
70    call WaitRequest(Request_write)     
71!$OMP BARRIER
72
73!$OMP MASTER
74    if (MPI_Rank==0) call WriteField(name,New_Field)
75    DEALLOCATE(New_Field)
76!$OMP END MASTER       
77!$OMP BARRIER
78    END SUBROUTINE write_field_u_gen
79
80
81  subroutine write_field1D_v(name,Field)
82    character(len=*)   :: name
83    real, dimension(:) :: Field
84
85    CALL write_field_v_gen(name,Field,1)
86
87  end subroutine write_field1D_v
88
89  subroutine write_field2D_v(name,Field)
90    implicit none
91     
92    character(len=*)   :: name
93    real, dimension(:,:) :: Field
94    integer :: ll
95   
96    ll=size(field,2)   
97    CALL write_field_v_gen(name,Field,ll)
98   
99    end subroutine write_field2D_v
100
101
102   SUBROUTINE write_field_v_gen(name,Field,ll)
[1823]103    USE parallel_lmdz
[1632]104    USE write_field
105    USE mod_hallo
[5271]106    USE dimensions_mod, ONLY: iim, jjm, llm, ndm
[5285]107USE paramet_mod_h
[5271]108implicit none
109
[5272]110
[1632]111     
112    character(len=*)   :: name
113    real, dimension(ijb_v:ije_v,ll) :: Field
114    real, allocatable,SAVE :: New_Field(:,:,:)
115    integer,dimension(0:mpi_size-1) :: jj_nb_master
[1848]116    type(Request),SAVE :: Request_write
117!$OMP THREADPRIVATE(Request_write)   
[1632]118    integer :: ll,i,jje,ije,jjn
119   
120   
121    jj_nb_master(:)=0
122    jj_nb_master(0)=jjp1
123
124!$OMP BARRIER
125!$OMP MASTER
126    allocate(New_Field(iip1,jjm,ll))
127!$OMP END MASTER
128!$OMP BARRIER
129
130   IF (pole_sud) THEN
131     jje=jj_end-1
132     ije=ij_end-iip1
133     jjn=jj_nb-1
134   ELSE
135     jje=jj_end
136     ije=ij_end
137     jjn=jj_nb
138   ENDIF
139   
140!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
141    DO i=1,ll   
142      New_Field(:,jj_begin:jje,i)=reshape(Field(ij_begin:ije,i),(/iip1,jjn/))
143    ENDDO
[1848]144!$OMP BARRIER   
[1632]145    call Register_SwapField(new_field,new_field,ip1jm,ll,jj_Nb_master,Request_write)
146    call SendRequest(Request_write)
147!$OMP BARRIER
148    call WaitRequest(Request_write)     
149!$OMP BARRIER
150
151!$OMP MASTER
152    if (MPI_Rank==0) call WriteField(name,New_Field)
153    DEALLOCATE(New_Field)
154!$OMP END MASTER       
155!$OMP BARRIER
156    END SUBROUTINE write_field_v_gen
157   
158end module write_field_loc
159 
Note: See TracBrowser for help on using the repository browser.