11#include "./base/base_uses.f90"
30 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_fb_buffer_types'
44 INTEGER :: ref_count = -1
46 INTEGER,
DIMENSION(:),
POINTER :: disps => null()
47 INTEGER,
DIMENSION(:),
POINTER :: data_1d => null()
48 END TYPE fb_buffer_i_data
56 TYPE(fb_buffer_i_data),
POINTER,
PRIVATE :: obj => null()
57 END TYPE fb_buffer_i_obj
71 INTEGER :: ref_count = -1
73 INTEGER,
DIMENSION(:),
POINTER :: disps => null()
74 REAL(KIND=
dp),
DIMENSION(:),
POINTER :: data_1d => null()
75 END TYPE fb_buffer_d_data
83 TYPE(fb_buffer_d_data),
POINTER,
PRIVATE :: obj => null()
88 MODULE PROCEDURE fb_buffer_i_add
89 MODULE PROCEDURE fb_buffer_d_add
92 INTERFACE fb_buffer_associate
93 MODULE PROCEDURE fb_buffer_i_associate
94 MODULE PROCEDURE fb_buffer_d_associate
95 END INTERFACE fb_buffer_associate
98 MODULE PROCEDURE fb_buffer_i_create
99 MODULE PROCEDURE fb_buffer_d_create
102 INTERFACE fb_buffer_calc_disps
103 MODULE PROCEDURE fb_buffer_i_calc_disps
104 MODULE PROCEDURE fb_buffer_d_calc_disps
105 END INTERFACE fb_buffer_calc_disps
107 INTERFACE fb_buffer_calc_sizes
108 MODULE PROCEDURE fb_buffer_i_calc_sizes
109 MODULE PROCEDURE fb_buffer_d_calc_sizes
110 END INTERFACE fb_buffer_calc_sizes
113 MODULE PROCEDURE fb_buffer_i_get
114 MODULE PROCEDURE fb_buffer_d_get
118 MODULE PROCEDURE fb_buffer_i_has_data
119 MODULE PROCEDURE fb_buffer_d_has_data
123 MODULE PROCEDURE fb_buffer_i_release
124 MODULE PROCEDURE fb_buffer_d_release
127 INTERFACE fb_buffer_retain
128 MODULE PROCEDURE fb_buffer_i_retain
129 MODULE PROCEDURE fb_buffer_d_retain
130 END INTERFACE fb_buffer_retain
133 MODULE PROCEDURE fb_buffer_i_nullify
134 MODULE PROCEDURE fb_buffer_d_nullify
138 MODULE PROCEDURE fb_buffer_i_replace
139 MODULE PROCEDURE fb_buffer_d_replace
151 SUBROUTINE fb_buffer_i_retain(buffer)
152 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
154 cpassert(
ASSOCIATED(buffer%obj))
155 buffer%obj%ref_count = buffer%obj%ref_count + 1
156 END SUBROUTINE fb_buffer_i_retain
163 SUBROUTINE fb_buffer_i_release(buffer)
164 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
166 IF (
ASSOCIATED(buffer%obj))
THEN
167 cpassert(buffer%obj%ref_count > 0)
168 buffer%obj%ref_count = buffer%obj%ref_count - 1
169 IF (buffer%obj%ref_count == 0)
THEN
170 buffer%obj%ref_count = 1
171 IF (
ASSOCIATED(buffer%obj%data_1d))
THEN
172 DEALLOCATE (buffer%obj%data_1d)
174 IF (
ASSOCIATED(buffer%obj%disps))
THEN
175 DEALLOCATE (buffer%obj%disps)
177 buffer%obj%ref_count = 0
178 DEALLOCATE (buffer%obj)
183 END SUBROUTINE fb_buffer_i_release
190 SUBROUTINE fb_buffer_i_nullify(buffer)
191 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
194 END SUBROUTINE fb_buffer_i_nullify
202 SUBROUTINE fb_buffer_i_associate(a, b)
203 TYPE(fb_buffer_i_obj),
INTENT(OUT) :: a
204 TYPE(fb_buffer_i_obj),
INTENT(IN) :: b
207 CALL fb_buffer_retain(a)
208 END SUBROUTINE fb_buffer_i_associate
216 PURE FUNCTION fb_buffer_i_has_data(buffer)
RESULT(res)
217 TYPE(fb_buffer_i_obj),
INTENT(IN) :: buffer
220 res =
ASSOCIATED(buffer%obj)
221 END FUNCTION fb_buffer_i_has_data
232 SUBROUTINE fb_buffer_i_create(buffer, &
237 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
238 INTEGER,
INTENT(IN),
OPTIONAL :: max_size, nslices
239 INTEGER,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: data_1d, sizes
241 INTEGER :: my_max_size, my_ndata, my_nslices
246 IF (
PRESENT(data_1d))
THEN
247 cpassert(
PRESENT(sizes))
250 cpassert(.NOT.
ASSOCIATED(buffer%obj))
251 ALLOCATE (buffer%obj)
256 NULLIFY (buffer%obj%data_1d, &
259 IF (
PRESENT(max_size)) my_max_size = max_size
260 IF (
PRESENT(nslices)) my_nslices = nslices
261 IF (
PRESENT(sizes))
THEN
262 my_nslices = min(my_nslices,
SIZE(sizes))
263 my_ndata = sum(sizes(1:my_nslices))
264 my_max_size = max(my_max_size, my_ndata)
267 ALLOCATE (buffer%obj%data_1d(my_max_size))
268 ALLOCATE (buffer%obj%disps(my_nslices))
269 buffer%obj%data_1d = 0
272 buffer%obj%n = my_nslices
274 IF (
PRESENT(sizes))
THEN
275 CALL fb_buffer_calc_disps(buffer, sizes)
278 IF (
PRESENT(data_1d))
THEN
279 check_ok =
SIZE(data_1d) .GE. my_max_size .AND. &
282 buffer%obj%data_1d(1:my_ndata) = data_1d(1:my_ndata)
285 buffer%obj%ref_count = 1
286 END SUBROUTINE fb_buffer_i_create
294 SUBROUTINE fb_buffer_i_add(buffer, data_1d)
295 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
296 INTEGER,
DIMENSION(:),
INTENT(IN) :: data_1d
298 INTEGER :: new_data_size, new_n, this_size
299 INTEGER,
DIMENSION(:),
POINTER :: new_data, new_disps
301 NULLIFY (new_disps, new_data)
303 this_size =
SIZE(data_1d)
304 new_n = buffer%obj%n + 1
305 new_data_size = buffer%obj%disps(new_n) + this_size
307 IF (
SIZE(buffer%obj%disps) .LT. new_n + 1)
THEN
308 ALLOCATE (new_disps(new_n*2))
310 new_disps(1:buffer%obj%n + 1) = buffer%obj%disps(1:buffer%obj%n + 1)
311 DEALLOCATE (buffer%obj%disps)
312 buffer%obj%disps => new_disps
314 IF (
SIZE(buffer%obj%data_1d) .LT. new_data_size)
THEN
315 ALLOCATE (new_data(new_data_size*2))
317 new_data(1:buffer%obj%disps(new_n)) = &
318 buffer%obj%data_1d(1:buffer%obj%disps(new_n))
319 DEALLOCATE (buffer%obj%data_1d)
320 buffer%obj%data_1d => new_data
323 buffer%obj%disps(new_n + 1) = new_data_size
324 buffer%obj%data_1d(buffer%obj%disps(new_n) + 1:new_data_size) = &
327 END SUBROUTINE fb_buffer_i_add
336 SUBROUTINE fb_buffer_i_calc_disps(buffer, sizes)
337 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
338 INTEGER,
DIMENSION(:),
INTENT(IN) :: sizes
342 cpassert(
SIZE(sizes) .GE. buffer%obj%n)
343 buffer%obj%disps(1) = 0
344 DO ii = 2, buffer%obj%n + 1
345 buffer%obj%disps(ii) = buffer%obj%disps(ii - 1) + sizes(ii - 1)
347 END SUBROUTINE fb_buffer_i_calc_disps
355 SUBROUTINE fb_buffer_i_calc_sizes(buffer, sizes)
356 TYPE(fb_buffer_i_obj),
INTENT(IN) :: buffer
357 INTEGER,
DIMENSION(:),
INTENT(OUT) :: sizes
361 cpassert(
SIZE(sizes) .GE. buffer%obj%n)
362 DO ii = 1, buffer%obj%n
363 sizes(ii) = buffer%obj%disps(ii + 1) - buffer%obj%disps(ii)
365 END SUBROUTINE fb_buffer_i_calc_sizes
385 SUBROUTINE fb_buffer_i_get(buffer, &
394 TYPE(fb_buffer_i_obj),
INTENT(IN) :: buffer
395 INTEGER,
INTENT(IN),
OPTIONAL :: i_slice
396 INTEGER,
INTENT(OUT),
OPTIONAL :: n, data_size
397 INTEGER,
DIMENSION(:),
INTENT(OUT),
OPTIONAL :: sizes, disps
398 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: data_1d
399 INTEGER,
DIMENSION(:, :),
OPTIONAL,
POINTER :: data_2d
400 INTEGER,
INTENT(IN),
OPTIONAL :: data_2d_ld
402 INTEGER :: ncols, slice_size
404 IF (
PRESENT(n)) n = buffer%obj%n
405 IF (
PRESENT(data_size)) data_size = buffer%obj%disps(buffer%obj%n + 1)
406 IF (
PRESENT(sizes))
THEN
407 CALL fb_buffer_calc_sizes(buffer, sizes)
409 IF (
PRESENT(disps))
THEN
410 cpassert(
SIZE(disps) .GE. buffer%obj%n)
411 disps(1:buffer%obj%n) = buffer%obj%disps(1:buffer%obj%n)
413 IF (
PRESENT(data_1d))
THEN
414 IF (
PRESENT(i_slice))
THEN
415 cpassert(i_slice .LE. buffer%obj%n)
416 data_1d => buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
417 buffer%obj%disps(i_slice + 1))
419 data_1d => buffer%obj%data_1d(1:buffer%obj%disps(buffer%obj%n + 1))
422 IF (
PRESENT(data_2d))
THEN
423 cpassert(
PRESENT(data_2d_ld))
424 cpassert(
PRESENT(i_slice))
430 slice_size = buffer%obj%disps(i_slice + 1) - buffer%obj%disps(i_slice)
431 ncols = slice_size/data_2d_ld
432 cpassert(slice_size == data_2d_ld*ncols)
433 data_2d(1:data_2d_ld, 1:ncols) => &
434 buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
435 buffer%obj%disps(i_slice + 1))
437 END SUBROUTINE fb_buffer_i_get
447 SUBROUTINE fb_buffer_i_replace(buffer, i_slice, data_1d)
448 TYPE(fb_buffer_i_obj),
INTENT(INOUT) :: buffer
449 INTEGER,
INTENT(IN) :: i_slice
450 INTEGER,
DIMENSION(:),
INTENT(IN) :: data_1d
452 INTEGER :: slice_size
454 slice_size = buffer%obj%disps(i_slice + 1) - buffer%obj%disps(i_slice)
455 cpassert(
SIZE(data_1d) == slice_size)
456 buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
457 buffer%obj%disps(i_slice + 1)) = data_1d
458 END SUBROUTINE fb_buffer_i_replace
467 SUBROUTINE fb_buffer_d_retain(buffer)
470 cpassert(
ASSOCIATED(buffer%obj))
471 buffer%obj%ref_count = buffer%obj%ref_count + 1
472 END SUBROUTINE fb_buffer_d_retain
479 SUBROUTINE fb_buffer_d_release(buffer)
482 IF (
ASSOCIATED(buffer%obj))
THEN
483 cpassert(buffer%obj%ref_count > 0)
484 buffer%obj%ref_count = buffer%obj%ref_count - 1
485 IF (buffer%obj%ref_count == 0)
THEN
486 buffer%obj%ref_count = 1
487 IF (
ASSOCIATED(buffer%obj%data_1d))
THEN
488 DEALLOCATE (buffer%obj%data_1d)
490 IF (
ASSOCIATED(buffer%obj%disps))
THEN
491 DEALLOCATE (buffer%obj%disps)
493 buffer%obj%ref_count = 0
494 DEALLOCATE (buffer%obj)
499 END SUBROUTINE fb_buffer_d_release
506 SUBROUTINE fb_buffer_d_nullify(buffer)
510 END SUBROUTINE fb_buffer_d_nullify
518 SUBROUTINE fb_buffer_d_associate(a, b)
523 CALL fb_buffer_retain(a)
524 END SUBROUTINE fb_buffer_d_associate
532 PURE FUNCTION fb_buffer_d_has_data(buffer)
RESULT(res)
536 res =
ASSOCIATED(buffer%obj)
537 END FUNCTION fb_buffer_d_has_data
548 SUBROUTINE fb_buffer_d_create(buffer, &
554 INTEGER,
INTENT(IN),
OPTIONAL :: max_size, nslices
555 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN),
OPTIONAL :: data_1d
556 INTEGER,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: sizes
558 INTEGER :: my_max_size, my_ndata, my_nslices
563 IF (
PRESENT(data_1d))
THEN
564 cpassert(
PRESENT(sizes))
567 cpassert(.NOT.
ASSOCIATED(buffer%obj))
568 ALLOCATE (buffer%obj)
573 NULLIFY (buffer%obj%data_1d, &
576 IF (
PRESENT(max_size)) my_max_size = max_size
577 IF (
PRESENT(nslices)) my_nslices = nslices
578 IF (
PRESENT(sizes))
THEN
579 my_nslices = min(my_nslices,
SIZE(sizes))
580 my_ndata = sum(sizes(1:my_nslices))
581 my_max_size = max(my_max_size, my_ndata)
584 ALLOCATE (buffer%obj%data_1d(my_max_size))
585 ALLOCATE (buffer%obj%disps(my_nslices + 1))
586 buffer%obj%data_1d = 0
589 buffer%obj%n = my_nslices
591 IF (
PRESENT(sizes))
THEN
592 CALL fb_buffer_calc_disps(buffer, sizes)
595 IF (
PRESENT(data_1d))
THEN
596 check_ok =
SIZE(data_1d) .GE. my_max_size .AND. &
599 buffer%obj%data_1d(1:my_ndata) = data_1d(1:my_ndata)
602 buffer%obj%ref_count = 1
603 END SUBROUTINE fb_buffer_d_create
611 SUBROUTINE fb_buffer_d_add(buffer, data_1d)
613 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: data_1d
615 INTEGER :: new_data_size, new_n, this_size
616 INTEGER,
DIMENSION(:),
POINTER :: new_disps
617 REAL(KIND=
dp),
DIMENSION(:),
POINTER :: new_data
619 NULLIFY (new_disps, new_data)
621 this_size =
SIZE(data_1d)
622 new_n = buffer%obj%n + 1
623 new_data_size = buffer%obj%disps(new_n) + this_size
625 IF (
SIZE(buffer%obj%disps) .LT. new_n + 1)
THEN
626 ALLOCATE (new_disps(new_n*2))
628 new_disps(1:buffer%obj%n + 1) = buffer%obj%disps(1:buffer%obj%n + 1)
629 DEALLOCATE (buffer%obj%disps)
630 buffer%obj%disps => new_disps
632 IF (
SIZE(buffer%obj%data_1d) .LT. new_data_size)
THEN
633 ALLOCATE (new_data(new_data_size*2))
635 new_data(1:buffer%obj%disps(new_n)) = &
636 buffer%obj%data_1d(1:buffer%obj%disps(new_n))
637 DEALLOCATE (buffer%obj%data_1d)
638 buffer%obj%data_1d => new_data
641 buffer%obj%disps(new_n + 1) = new_data_size
642 buffer%obj%data_1d(buffer%obj%disps(new_n) + 1:new_data_size) = &
645 END SUBROUTINE fb_buffer_d_add
654 SUBROUTINE fb_buffer_d_calc_disps(buffer, sizes)
656 INTEGER,
DIMENSION(:),
INTENT(IN) :: sizes
660 cpassert(
SIZE(sizes) .GE. buffer%obj%n)
661 buffer%obj%disps(1) = 0
662 DO ii = 2, buffer%obj%n + 1
663 buffer%obj%disps(ii) = buffer%obj%disps(ii - 1) + sizes(ii - 1)
665 END SUBROUTINE fb_buffer_d_calc_disps
673 SUBROUTINE fb_buffer_d_calc_sizes(buffer, sizes)
675 INTEGER,
DIMENSION(:),
INTENT(OUT) :: sizes
679 cpassert(
SIZE(sizes) .GE. buffer%obj%n)
680 DO ii = 1, buffer%obj%n
681 sizes(ii) = buffer%obj%disps(ii + 1) - buffer%obj%disps(ii)
683 END SUBROUTINE fb_buffer_d_calc_sizes
703 SUBROUTINE fb_buffer_d_get(buffer, &
713 INTEGER,
INTENT(IN),
OPTIONAL :: i_slice
714 INTEGER,
INTENT(OUT),
OPTIONAL :: n, data_size
715 INTEGER,
DIMENSION(:),
INTENT(OUT),
OPTIONAL :: sizes, disps
716 REAL(KIND=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: data_1d
717 REAL(KIND=
dp),
DIMENSION(:, :),
OPTIONAL,
POINTER :: data_2d
718 INTEGER,
INTENT(IN),
OPTIONAL :: data_2d_ld
720 INTEGER :: ncols, slice_size
722 IF (
PRESENT(n)) n = buffer%obj%n
723 IF (
PRESENT(data_size)) data_size = buffer%obj%disps(buffer%obj%n + 1)
724 IF (
PRESENT(sizes))
THEN
725 CALL fb_buffer_calc_sizes(buffer, sizes)
727 IF (
PRESENT(disps))
THEN
728 cpassert(
SIZE(disps) .GE. buffer%obj%n)
729 disps(1:buffer%obj%n) = buffer%obj%disps(1:buffer%obj%n)
731 IF (
PRESENT(data_1d))
THEN
732 IF (
PRESENT(i_slice))
THEN
733 cpassert(i_slice .LE. buffer%obj%n)
734 data_1d => buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
735 buffer%obj%disps(i_slice + 1))
737 data_1d => buffer%obj%data_1d(1:buffer%obj%disps(buffer%obj%n + 1))
740 IF (
PRESENT(data_2d))
THEN
741 cpassert(
PRESENT(data_2d_ld))
742 cpassert(
PRESENT(i_slice))
748 slice_size = buffer%obj%disps(i_slice + 1) - buffer%obj%disps(i_slice)
749 ncols = slice_size/data_2d_ld
750 cpassert(slice_size == data_2d_ld*ncols)
751 data_2d(1:data_2d_ld, 1:ncols) => &
752 buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
753 buffer%obj%disps(i_slice + 1))
755 END SUBROUTINE fb_buffer_d_get
765 SUBROUTINE fb_buffer_d_replace(buffer, i_slice, data_1d)
767 INTEGER,
INTENT(IN) :: i_slice
768 REAL(KIND=
dp),
DIMENSION(:),
INTENT(IN) :: data_1d
770 INTEGER :: slice_size
772 slice_size = buffer%obj%disps(i_slice + 1) - buffer%obj%disps(i_slice)
773 cpassert(
SIZE(data_1d) == slice_size)
774 buffer%obj%data_1d(buffer%obj%disps(i_slice) + 1: &
775 buffer%obj%disps(i_slice + 1)) = data_1d
776 END SUBROUTINE fb_buffer_d_replace
Defines the basic variable types.
integer, parameter, public dp
object/pointer wrapper for fb_buffer object