Line data Source code
1 : module buffer
2 :
3 : !-----------------------------------------------------------------------
4 : !
5 : ! Purpose:
6 : ! LOW level handler for f90 arrays.
7 : !
8 : ! Author: J. Edwards
9 : !
10 : ! This file is used with genf90.pl to generate buffer.F90
11 : !
12 : !-----------------------------------------------------------------------
13 :
14 : use shr_kind_mod, only: r8 => shr_kind_r8, r4=> shr_kind_r4, i4=> shr_kind_i4
15 : use cam_logfile, only: iulog
16 : use cam_abortutils, only: endrun
17 :
18 : implicit none
19 : private
20 : ! The maximum number of dims in a fortran array
21 : #define MAXDIMS 7
22 :
23 :
24 : type buffer_field_default_type
25 : private
26 : real(r8), pointer :: data(:,:,:,:,:,:,:) => null()
27 : end type buffer_field_default_type
28 :
29 : ! TYPE int,double,real
30 : type buffer_field_{TYPE}
31 : private
32 : {VTYPE}, pointer :: data(:,:,:,:,:,:,:) => null()
33 : end type buffer_field_{TYPE}
34 :
35 : integer(i4), parameter,public :: dtype_i4=1
36 : real(r8), parameter,public :: dtype_r8=1_r8
37 : real(r4), parameter,public :: dtype_r4=1_r4
38 :
39 : interface buffer_field_deallocate
40 : ! TYPE int,double,real
41 : module procedure buffer_field_deallocate_{TYPE}
42 : end interface
43 :
44 : interface buffer_field_allocate
45 : ! TYPE int,double,real
46 : module procedure buffer_field_allocate_{TYPE}
47 : end interface
48 :
49 : interface buffer_set_field
50 : ! TYPE int,double,real
51 : module procedure buffer_set_field_const_{TYPE}
52 : ! DIMS 1,2,3,4,5,6,7
53 : ! TYPE int,double,real
54 : module procedure buffer_set_field_{DIMS}d_{TYPE}
55 : end interface
56 :
57 : interface buffer_get_field_ptr
58 : ! DIMS 1,2,3,4,5,6,7
59 : ! TYPE int,double,real
60 : module procedure buffer_get_field_ptr_{DIMS}d_{TYPE}
61 : end interface
62 :
63 : public :: buffer_field_deallocate, buffer_field_allocate, buffer_set_field, buffer_get_field_ptr, buffer_field_default_type
64 : public :: buffer_field_is_alloc
65 :
66 :
67 : CONTAINS
68 :
69 : ! TYPE int,double,real
70 1448928 : subroutine buffer_field_deallocate_{TYPE}(bfg, dtype)
71 :
72 : type(buffer_field_default_type),intent(inout) :: bfg
73 : {VTYPE}, intent(in) :: dtype
74 :
75 : type(buffer_field_{TYPE}) :: b1
76 :
77 1448928 : b1 = transfer(bfg, b1)
78 :
79 1448928 : if(.not.associated(b1%data)) then
80 0 : call endrun('Attempt to deallocate unassociated array ptr')
81 : end if
82 :
83 1448928 : deallocate(b1%data)
84 :
85 1448928 : nullify(bfg%data)
86 :
87 1448928 : end subroutine buffer_field_deallocate_{TYPE}
88 :
89 1204168465 : logical function buffer_field_is_alloc(bfg)
90 : type(buffer_field_default_type),intent(in) :: bfg
91 :
92 1204168465 : buffer_field_is_alloc = associated(bfg%data)
93 :
94 1204168465 : end function buffer_field_is_alloc
95 :
96 :
97 : ! TYPE int,double,real
98 1448928 : subroutine buffer_field_allocate_{TYPE} (bfg, dimsizes, dtype)
99 :
100 : type(buffer_field_default_type),intent(inout) :: bfg
101 : integer, intent(in) :: dimsizes(:)
102 : integer :: alldimsizes( MAXDIMS )
103 : {VTYPE}, intent(in) :: dtype
104 : integer :: ierr
105 :
106 : type(buffer_field_{TYPE}) :: b1
107 :
108 11591424 : alldimsizes(:) = 1
109 11591424 : alldimsizes(1:size(dimsizes)) = dimsizes
110 :
111 1448928 : if(associated(bfg%data)) then
112 0 : call endrun('Attempt to allocate array to associated ptr')
113 : end if
114 :
115 : allocate(b1%data(alldimsizes(1),alldimsizes(2),alldimsizes(3),alldimsizes(4),&
116 13040352 : alldimsizes(5),alldimsizes(6),alldimsizes(7)),stat=ierr)
117 :
118 1448928 : if(ierr/=0) then
119 0 : call endrun("allocate failed")
120 : end if
121 :
122 1448928 : bfg = transfer(b1,bfg)
123 :
124 1448928 : end subroutine buffer_field_allocate_{TYPE}
125 :
126 : ! TYPE int,double,real
127 297216 : subroutine buffer_set_field_const_{TYPE}(bfg, const, start,kount)
128 : type(buffer_field_default_type) :: bfg
129 : {VTYPE}, intent(in) :: const
130 : integer, intent(in), optional :: start(:),kount(:)
131 :
132 : type(buffer_field_{TYPE}) :: ptr
133 :
134 : integer :: i, ns, strt(7), fin(7), cnt(7)
135 :
136 297216 : ptr = transfer(bfg, ptr)
137 :
138 297216 : if(present(start).and.present(kount)) then
139 0 : strt(:) = 1
140 0 : cnt = shape(ptr%data)
141 :
142 0 : ns = size(start)
143 0 : strt(1:ns) = start
144 :
145 0 : fin = strt+cnt-1
146 :
147 0 : do i=1,ns
148 0 : fin(i) = strt(i)+kount(i)-1
149 0 : if(strt(i)<1 .or. fin(i)>cnt(i)) then
150 0 : call endrun('Start plus kount exceeds dimension bounds')
151 : endif
152 : enddo
153 :
154 :
155 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),&
156 0 : strt(5):fin(5),strt(6):fin(6),strt(7):fin(7))=const
157 : else
158 440734176 : ptr%data = const
159 : endif
160 :
161 297216 : end subroutine buffer_set_field_const_{TYPE}
162 :
163 : !=========================================================================================
164 : !
165 : ! Given a physics_buffer chunk and an index return a pointer to a field chunk
166 : !
167 : !
168 :
169 : ! TYPE int,double,real
170 : ! DIMS 1,2,3,4,5,6,7
171 1193047057 : subroutine buffer_get_field_ptr_{DIMS}d_{TYPE}(bfg, field, start,kount)
172 : type(buffer_field_default_type), intent(in) :: bfg
173 : {VTYPE}, pointer :: field{DIMSTR}
174 : integer, intent(in), optional :: start(:), kount(:)
175 : type(buffer_field_{TYPE}), target :: ptr
176 :
177 : integer :: ns, strt(7), fin(7), cnt(7)
178 :
179 1193047057 : ptr = transfer(bfg, ptr)
180 :
181 :
182 9544376456 : strt(:) = 1
183 9544376456 : cnt = shape(ptr%data)
184 :
185 1193047057 : if(present(start)) then
186 243153648 : ns = size(start)
187 1215313128 : strt(1:ns) = start
188 : end if
189 1193047057 : if(present(kount)) then
190 1215313128 : cnt(1:ns) = kount
191 : end if
192 9544376456 : fin = strt+cnt-1
193 :
194 : #if ({DIMS}==1)
195 188639280 : field => ptr%data(strt(1):fin(1),strt(2),strt(3),strt(4),strt(5),strt(6),strt(7))
196 : #elif ({DIMS}==2)
197 981017497 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3),strt(4),strt(5),strt(6),strt(7))
198 : #elif ({DIMS}==3)
199 22350024 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4),strt(5),strt(6),strt(7))
200 : #elif ({DIMS}==4)
201 0 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5),strt(6),strt(7))
202 : #elif ({DIMS}==5)
203 693504 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5):fin(5),strt(6),strt(7))
204 : #elif ({DIMS}==6)
205 0 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),&
206 346752 : strt(5):fin(5),strt(6):fin(6),strt(7))
207 : #else
208 0 : field => ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),&
209 0 : strt(5):fin(5),strt(6):fin(6),strt(7):fin(7))
210 : #endif
211 :
212 1193047057 : end subroutine buffer_get_field_ptr_{DIMS}d_{TYPE}
213 :
214 : ! TYPE int,double,real
215 : ! DIMS 1,2,3,4,5,6,7
216 10452096 : subroutine buffer_set_field_{DIMS}d_{TYPE}(bfg,field,start,kount)
217 : type(buffer_field_default_type) :: bfg
218 : {VTYPE},intent(in) :: field{DIMSTR}
219 : integer,intent(in),optional :: start(:),kount(:)
220 : type(buffer_field_{TYPE}) :: ptr
221 :
222 : integer :: i, nc, strt(7), fin(7), cnt(7)
223 :
224 10452096 : ptr = transfer(bfg,ptr)
225 10452096 : if(present(start).and.present(kount)) then
226 23876352 : strt(:) = 1
227 23876352 : cnt = shape(ptr%data)
228 :
229 2984544 : nc=size(start)
230 14922720 : strt(1:nc) = start
231 23876352 : fin = strt+cnt-1
232 :
233 8953632 : do i=1,nc
234 5969088 : fin(i) = strt(i)+kount(i)-1
235 8953632 : if(strt(i)<1 .or. fin(i)>cnt(i)) then
236 0 : call endrun('Start plus kount exceeds dimension bounds')
237 : endif
238 : enddo
239 :
240 :
241 : #if ({DIMS}==1)
242 25315992 : ptr%data(strt(1):fin(1),strt(2),strt(3),strt(4),strt(5),strt(6),strt(7))=field
243 : #elif ({DIMS}==2)
244 2323627992 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3),strt(4),strt(5),strt(6),strt(7))=field
245 : #elif ({DIMS}==3)
246 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4),strt(5),strt(6),strt(7))=field
247 : #elif ({DIMS}==4)
248 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5),strt(6),strt(7))=field
249 : #elif ({DIMS}==5)
250 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),strt(5):fin(5),strt(6),strt(7))=field
251 : #elif ({DIMS}==6)
252 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),&
253 0 : strt(5):fin(5),strt(6):fin(6),strt(7))=field
254 : #else
255 0 : ptr%data(strt(1):fin(1),strt(2):fin(2),strt(3):fin(3),strt(4):fin(4),&
256 0 : strt(5):fin(5),strt(6):fin(6),strt(7):fin(7))=field
257 : #endif
258 : else
259 : #if ({DIMS}==1)
260 76263768 : ptr%data(:,1,1,1,1,1,1) = field
261 : #elif({DIMS}==2)
262 4762966896 : ptr%data(:,:,1,1,1,1,1) = field
263 : #elif({DIMS}==3)
264 0 : ptr%data(:,:,:,1,1,1,1) = field
265 : #elif({DIMS}==4)
266 0 : ptr%data(:,:,:,:,1,1,1) = field
267 : #elif({DIMS}==5)
268 0 : ptr%data(:,:,:,:,:,1,1) = field
269 : #elif({DIMS}==6)
270 0 : ptr%data(:,:,:,:,:,:,1) = field
271 : #else
272 0 : ptr%data = field
273 : #endif
274 : end if
275 10452096 : end subroutine buffer_set_field_{DIMS}d_{TYPE}
276 :
277 :
278 :
279 0 : end module buffer
280 :
281 :
282 :
283 :
284 :
|