source: LMDZ6/branches/Optimisation_LMDZ/libf/phylmd/cv3a_compress.f90 @ 3707

Last change on this file since 3707 was 3707, checked in by adurocher, 4 years ago

New generic implementation for compression

  • Property svn:executable set to *
File size: 11.8 KB
Line 
1module cv3a_compress_mod
2  implicit none
3  private
4 
5  integer, parameter :: max_array_count = 100 ! Max number of arrays (of each type) added to a array_list
6  integer, parameter, public :: COMPRESS_MODE_COMPRESS=1, COMPRESS_MODE_COPY=2 ! Constants for different compress modes
7  integer, save, public :: compress_mode = COMPRESS_MODE_COMPRESS ! Current compress mode. Can be changed by user, but be sure to use only one type of compression.
8  !$omp threadprivate(compress_mode)
9 
10  ! Containers for arrays (Internal structure)
11  type array_r_t
12    real, pointer, contiguous :: buf(:,:,:) => null() ! max dim for real is 3
13  end type
14  type array_i_t
15    integer, pointer :: buf(:) => null() ! max dim for integer is 1
16  end type
17 
18  ! List of arrays (Black box structure)
19  type array_list
20    integer :: array_count_r = 0 ! number of real arrays
21    type(array_r_t) :: arrays_r(2,max_array_count) ! array of real arrays
22    logical :: arrays_r_init(max_array_count) = .false. ! Should arrays (2,:) be initialized before uncompress?
23    integer :: array_count_i = 0 ! number of int arrays
24    type(array_i_t) :: arrays_i(2,max_array_count) ! array of int arrays
25    logical :: arrays_i_init(max_array_count) = .false. ! Should arrays (2,:) be initialized before uncompress?
26  end type
27 
28  ! Information initialized by compress to transmit instructions for uncompress
29  type compress_data_t
30    integer :: ncum ! Number of cells in compressed arrays (first dimension)
31    integer, allocatable :: idcum(:) ! Correspondance between compressed and uncompressed arrays uncomp(idcum(i)) == comp(i)
32    logical, allocatable :: mask(:) ! Mask used for compression
33  end type
34 
35  public :: array_list, add_array_i1, add_array_r1, add_array_r2, add_array_r3, &
36            cv3a_compress, cv3a_uncompress, compress_data_t
37 
38  contains
39 
40  ! Add a pair of array to an array list (integer(:))
41  subroutine add_array_i1(arrays, array_in, array_out, init)
42    type(array_list) :: arrays ! List of arrays to add arrays to
43    integer, target :: array_in(:) ! Input array (should not be modified)
44    integer, target :: array_out(:) ! Output array (will be modified)
45    logical, optional, intent(in) :: init ! For uncompress : zero-initialize output array before uncompress
46    integer :: c   
47    arrays%array_count_i = arrays%array_count_i + 1
48    c = arrays%array_count_i
49    if( c > max_array_count ) call abort_physic('cv3a_compress', "Too many arrays : increase cv3a_compress_mod%max_array_count",1)
50    if( present(init) ) arrays%arrays_i_init(c)=init
51    arrays%arrays_i(1,c)%buf => array_in
52    arrays%arrays_i(2,c)%buf => array_out
53  end subroutine
54
55  ! Add a pair of array to an array list (real(:))
56  subroutine add_array_r1(arrays, array_in, array_out, init)
57    type(array_list) :: arrays
58    real, target :: array_in(:), array_out(:)
59    logical, optional, intent(in) :: init
60    integer :: c
61    arrays%array_count_r = arrays%array_count_r + 1
62    c = arrays%array_count_r
63    if( c > max_array_count ) call abort_physic('cv3a_compress', "Too many arrays : increase cv3a_compress_mod%max_array_count",1)
64    if( present(init) ) arrays%arrays_r_init(c)=init
65    arrays%arrays_r(1,c)%buf(1:size(array_in),1:1,1:1) => array_in
66    arrays%arrays_r(2,c)%buf(1:size(array_out),1:1,1:1) => array_out
67  end subroutine
68 
69  ! Add a pair of array to an array list (real(:,:))
70  subroutine add_array_r2(arrays, array_in, array_out, init)
71    type(array_list) :: arrays
72    real, target, contiguous :: array_in(:,:), array_out(:,:)
73    logical, optional, intent(in) :: init
74    integer :: c
75    arrays%array_count_r = arrays%array_count_r + 1
76    c = arrays%array_count_r
77    if( c > max_array_count ) call abort_physic('cv3a_compress', "Too many arrays : increase cv3a_compress_mod%max_array_count",1)
78    if( present(init) ) arrays%arrays_r_init(c)=init   
79    arrays%arrays_r(1,c)%buf(1:size(array_in,1),1:size(array_in,2),1:1) => array_in
80    arrays%arrays_r(2,c)%buf(1:size(array_out,1),1:size(array_out,2),1:1) => array_out
81  end subroutine
82 
83  ! Add a pair of array to an array list (real(:,:,:))
84  subroutine add_array_r3(arrays, array_in, array_out, init)
85    type(array_list) :: arrays
86    real, target :: array_in(:,:,:), array_out(:,:,:)
87    logical, optional, intent(in) :: init
88    integer :: c
89    arrays%array_count_r = arrays%array_count_r + 1
90    c = arrays%array_count_r
91    if( c > max_array_count ) call abort_physic('cv3a_compress', "Too many arrays : increase cv3a_compress_mod%max_array_count",1)
92    if( present(init) ) arrays%arrays_r_init(c)=init   
93    arrays%arrays_r(1,c)%buf => array_in
94    arrays%arrays_r(2,c)%buf => array_out
95  end subroutine
96   
97  ! Compress arrays from 'arrays' according to 'mask'.
98  ! Input arrays in 'arrays' should have 'len' as first dimension
99  ! 2D arrays (real only) are compressed on every vertical layer
100  ! Note : compressing 3d arrays is not supported yet
101  subroutine cv3a_compress(len, mask, arrays, d)
102    integer, INTENT (IN) :: len ! lenght of the contiguous dimension of arrays to compress
103    logical, INTENT (IN) :: mask(len) ! Mask of convective cells
104    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
105    type(compress_data_t), intent(out) :: d ! Used to store data about compression to be able to uncompress
106   
107    !real, parameter :: comp_threshold = 1
108    integer :: i, ncum, ncum_max
109   
110    ! Output arrays are supposed to be able to contain compressed data
111    ncum_max = max( size(arrays%arrays_i(2,1)%buf), len ) 
112    allocate(d%idcum(ncum_max))
113    allocate(d%mask(len))
114    d%mask(:) = mask(:)
115   
116    ! Generate permutation idcum() and number of active elements 'ncum' for compression
117    ncum = 0
118    do i = 1, len
119      if (mask(i)) then
120        ncum = ncum + 1
121        d%idcum(ncum) = i
122      end if
123    end do
124    if( ncum > ncum_max ) call abort_physic("cv3a_compress", "internal error : idcum too small", 1)
125   
126    !if (ncum < len * comp_threshold) compress_mode=COMPRESS_MODE_COMPRESS   
127    select case (compress_mode)
128      case (COMPRESS_MODE_COMPRESS)
129        call cv3a_compress_compress(arrays, ncum, d%idcum)
130      case (COMPRESS_MODE_COPY)
131        call cv3a_compress_copy(len, arrays, ncum, d%idcum)
132      case default
133        call abort_physic("cv3a_compress", "Unknown compress mode", 1)
134    end select
135   
136    !if( ncum /= get_compress_size(len, mask) ) call abort_physic("cv3a_compress", "internal error : get_compress_size mismatch", 1)
137   
138    d%ncum = ncum
139
140  end subroutine
141
142  ! cv3a_compress in mode COMPRESS_MODE_COMPRESS
143  ! Output arrays have ncum elements that correspond to the idcum(i)
144  ! elements of the input arrays : a_out(i) = a_in(idcum(i))
145  ! compressing 3D arrays is not supported
146  subroutine cv3a_compress_compress(arrays, ncum, idcum)
147    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
148    integer, intent(in) :: ncum ! Number of cells in compressed arrays (first dimension)
149    integer, intent(in) :: idcum(ncum) ! Correspondance between compressed and uncompressed arrays uncomp(idcum(i)) == comp(i)
150   
151    integer :: i,j,k,a
152     
153    ! Compress int arrays
154    do a = 1, arrays%array_count_i
155      do i = 1, ncum
156        arrays%arrays_i(2,a)%buf(i) = arrays%arrays_i(1,a)%buf(idcum(i))
157      end do
158    end do
159    ! Compress real arrays
160    do a = 1, arrays%array_count_r
161      if( size(arrays%arrays_r(2,a)%buf, 3) /= 1 ) call abort_physic("cv3a_compress", "cannot compress 3D array", 1)
162      do k = 1, size(arrays%arrays_r(2,a)%buf, 3)
163        do j = 1, size(arrays%arrays_r(2,a)%buf, 2)
164         do i = 1, ncum
165           arrays%arrays_r(2,a)%buf(i,j,k) = arrays%arrays_r(1,a)%buf(idcum(i),j,k)
166         end do
167        end do
168      end do           
169    end do 
170  end subroutine
171 
172  ! cv3a_compress in mode COMPRESS_MODE_COPY
173  subroutine cv3a_compress_copy(len, arrays, ncum, idcum)
174    integer, INTENT (IN) :: len ! lenght of the contiguous dimension of arrays to compress
175    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
176    integer, intent(out) :: ncum ! Number of cells in compressed arrays (first dimension)
177    integer, intent(out) :: idcum(len) ! Correspondance between compressed and uncompressed arrays uncomp(idcum(i)) == comp(i)
178    integer :: i 
179
180    ! Ignore and ovewrite computed permutation
181    ncum = len
182    DO i = 1,len
183      idcum(i) = i
184    ENDDO 
185    ! Copy int arrays
186    do i = 1, arrays%array_count_i
187      arrays%arrays_i(2,i)%buf(:) = arrays%arrays_i(1,i)%buf(:)
188    end do
189    ! Copy real arrays
190    do i = 1, arrays%array_count_r
191      arrays%arrays_r(2,i)%buf(:,:,:) = arrays%arrays_r(1,i)%buf(:,:,:)
192    end do
193  end subroutine
194 
195  ! Arrays in 'arrays' are uncompressed,
196  ! non-convective cells are zero-initialized if init=.true. was used when adding the array
197  ! 2D arrays are uncompressed up to level 'nl' (cv3param.h)
198  ! 3D arrays are uncompressed up to level (:,nl,nl) (cv3param.h)
199  ! after those levels cells are zero-initialized when init=.true.
200  subroutine cv3a_uncompress(len, d, arrays)
201    integer, intent(in) :: len ! Number of cells in uncompressed arrays (first dimension)
202    type(compress_data_t), intent(in) :: d ! Data to uncompress arrays
203    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
204   
205    select case (compress_mode)
206      case (COMPRESS_MODE_COMPRESS)
207        call cv3a_uncompress_uncompress(d%ncum, d%idcum, arrays)
208      case (COMPRESS_MODE_COPY)
209        call cv3a_uncompress_copy(d%ncum, d%mask, arrays)
210      case default
211        call abort_physic("cv3a_uncompress", "Unknown uncompress mode", 1)
212    end select 
213  end subroutine
214 
215  ! cv3a_uncompress in mode COMPRESS_MODE_COMPRESS
216  ! a_out(idcum(i)) = a_in(i)
217  subroutine cv3a_uncompress_uncompress(ncum, idcum, arrays)
218    integer, intent(in) :: ncum ! Number of cells in compressed arrays (first dimension)
219    integer, intent(in) :: idcum(ncum) ! Correspondance between compressed and uncompressed arrays uncomp(idcum(i)) == comp(i)
220    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
221   
222    include "cv3param.h"
223   
224    integer :: a, k, j, i 
225         
226    ! uncompress int arrays
227    do a = 1, arrays%array_count_i
228        if(arrays%arrays_i_init(a)) arrays%arrays_i(2,a)%buf(:) = 0
229        do i = 1, ncum
230          arrays%arrays_i(2,a)%buf(idcum(i)) = arrays%arrays_i(1,a)%buf(i)
231        end do
232     end do
233    ! Compress real arrays
234    do a = 1, arrays%array_count_r
235      if(arrays%arrays_r_init(a))  arrays%arrays_r(2,a)%buf(:,:,:) = 0
236      do k = 1, min(size(arrays%arrays_r(2,a)%buf, 3), nl)
237        do j = 1, min(size(arrays%arrays_r(2,a)%buf, 2), nl)
238          do i = 1, ncum
239            arrays%arrays_r(2,a)%buf(idcum(i),j,k) = arrays%arrays_r(1,a)%buf(i,j,k)
240          end do
241        end do
242      end do
243    end do               
244  end subroutine
245
246  subroutine cv3a_uncompress_copy(len, mask, arrays)
247    integer, intent(in) :: len
248    logical, intent(in) :: mask(:)
249    type(array_list), intent(inout) :: arrays ! List of arrays to compress (first compressed into second)
250   
251    include "cv3param.h"
252   
253    integer :: a, i, j, k
254   
255    ! Copy int arrays
256    do a = 1, arrays%array_count_i
257      if(arrays%arrays_i_init(a)) arrays%arrays_i(2,a)%buf(:) = 0
258      do i = 1, len
259        if(mask(i)) arrays%arrays_i(2,a)%buf(i) = arrays%arrays_i(1,a)%buf(i)
260      end do
261    end do
262    ! Copy real arrays
263    do a = 1, arrays%array_count_r
264      if(arrays%arrays_i_init(a)) arrays%arrays_r(2,a)%buf(:,:,:) = 0       
265      do k = 1, min(size(arrays%arrays_r(2,a)%buf, 3), nl)
266        do j = 1, min(size(arrays%arrays_r(2,a)%buf, 2), nl)
267          do i = 1, len
268            if(mask(i)) arrays%arrays_r(2,a)%buf(:,j,k) = merge( arrays%arrays_r(1,a)%buf(:,j,k), 0., mask)
269          end do
270        end do
271      end do
272    end do
273  end subroutine
274 
275end module
Note: See TracBrowser for help on using the repository browser.