source: LMDZ5/branches/LMDZ5-DOFOCO/libf/dyn3dmem/write_field_loc.F90 @ 4019

Last change on this file since 4019 was 1632, checked in by Laurent Fairhead, 12 years ago

Import initial du répertoire dyn3dmem

Attention! ceci n'est qu'une version préliminaire du code "basse mémoire":
le code contenu dans ce répertoire est basé sur la r1320 et a donc besoin
d'être mis à jour par rapport à la dynamique parallèle d'aujourd'hui.
Ce code est toutefois mis à disposition pour circonvenir à des problèmes
de mémoire que certaines configurations du modèle pourraient rencontrer.
Dans l'état, il compile et tourne sur vargas et au CCRT


Initial import of dyn3dmem

Warning! this is just a preliminary version of the memory light code:
it is based on r1320 of the code and thus needs to be updated before
it can replace the present dyn3dpar code. It is nevertheless put at your
disposal to circumvent some memory problems some LMDZ configurations may
encounter. In its present state, it will compile and run on vargas and CCRT

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