Forráskód Böngészése

- add support for cl_args
- add a native Fortran version of the former f90_example code

Olivier Aumage 9 éve
szülő
commit
15c23fb8fe

+ 49 - 0
examples/native_fortran/Makefile.nf_example

@@ -0,0 +1,49 @@
+# StarPU --- Runtime system for heterogeneous multicore architectures.
+#
+# Copyright (C) 2015  ONERA
+# Copyright (C) 2015-2016  Inria
+#
+# StarPU is free software; you can redistribute it and/or modify
+# it under the terms of the GNU Lesser General Public License as published by
+# the Free Software Foundation; either version 2.1 of the License, or (at
+# your option) any later version.
+#
+# StarPU is distributed in the hope that it will be useful, but
+# WITHOUT ANY WARRANTY; without even the implied warranty of
+# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
+#
+# See the GNU Lesser General Public License in COPYING.LGPL for more details.
+
+PROG = nf_example
+
+FSTARPU_MOD = $(shell pkg-config --cflags-only-I starpu-1.3|sed -e 's/^\([^ ]*starpu\/1.3\).*$$/\1/;s/^.* //;s/^-I//')/fstarpu_mod.f90
+
+SRCSF = mod_types.f90		\
+	mod_compute.f90		\
+	nf_example.f90
+
+FC = gfortran
+
+FCFLAGS = -fdefault-real-8 -J. -g
+LDLIBS =  $(shell pkg-config --libs starpu-1.3)
+
+OBJS = fstarpu_mod.o $(SRCSF:%.f90=%.o)
+
+.phony: all clean
+all: $(PROG)
+
+$(PROG): $(OBJS)
+	$(FC) $(LDFLAGS) -o $@ $^ $(LDLIBS)
+
+fstarpu_mod.o: $(FSTARPU_MOD)
+	$(FC) $(FCFLAGS) -c -o $@ $<
+
+%.o: %.f90
+	$(FC) $(FCFLAGS) -c -o $@ $<
+
+clean:
+	rm -fv *.o *.mod $(PROG)
+
+# modfiles generation dependences
+mod_compute.o: mod_compute.f90 mod_types.o fstarpu_mod.o
+nf_example.o: nf_example.f90 mod_types.o mod_compute.o fstarpu_mod.o

+ 51 - 0
examples/native_fortran/Makefile.nf_matrix

@@ -0,0 +1,51 @@
+# StarPU --- Runtime system for heterogeneous multicore architectures.
+#
+# Copyright (C) 2016  Inria
+#
+# StarPU is free software; you can redistribute it and/or modify
+# it under the terms of the GNU Lesser General Public License as published by
+# the Free Software Foundation; either version 2.1 of the License, or (at
+# your option) any later version.
+#
+# StarPU is distributed in the hope that it will be useful, but
+# WITHOUT ANY WARRANTY; without even the implied warranty of
+# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
+#
+# See the GNU Lesser General Public License in COPYING.LGPL for more details.
+
+PROG = nf_matrix
+
+FSTARPU_MOD = $(shell pkg-config --cflags-only-I starpu-1.3|sed -e 's/^\([^ ]*starpu\/1.3\).*$$/\1/;s/^.* //;s/^-I//')/fstarpu_mod.f90
+
+SRCSF = nf_matrix.f90		\
+	nf_codelets.f90
+
+FC = gfortran
+CC = gcc
+
+CFLAGS = -g $(shell pkg-config --cflags starpu-1.3)
+FCFLAGS = -fdefault-real-8 -J. -g
+LDLIBS =  $(shell pkg-config --libs starpu-1.3)
+
+OBJS = $(SRCSC:%.c=%.o) fstarpu_mod.o $(SRCSF:%.f90=%.o)
+
+.phony: all clean
+all: $(PROG)
+
+$(PROG): $(OBJS)
+	$(FC) $(LDFLAGS) -o $@ $^ $(LDLIBS)
+
+%.o: %.c
+	$(CC) $(CFLAGS) -c -o $@ $<
+
+fstarpu_mod.o: $(FSTARPU_MOD)
+	$(FC) $(FCFLAGS) -c -o $@ $<
+
+%.o: %.f90
+	$(FC) $(FCFLAGS) -c -o $@ $<
+
+clean:
+	rm -fv *.o *.mod $(PROG)
+
+nf_matrix.o: nf_matrix.f90 nf_codelets.o fstarpu_mod.o
+nf_codelets.o: fstarpu_mod.o

examples/native_fortran/Makefile → examples/native_fortran/Makefile.nf_vector


+ 132 - 0
examples/native_fortran/mod_compute.f90

@@ -0,0 +1,132 @@
+! StarPU --- Runtime system for heterogeneous multicore architectures.
+!
+! Copyright (C) 2015  ONERA
+! Copyright (C) 2015-2016  Inria
+!
+! StarPU is free software; you can redistribute it and/or modify
+! it under the terms of the GNU Lesser General Public License as published by
+! the Free Software Foundation; either version 2.1 of the License, or (at
+! your option) any later version.
+!
+! StarPU is distributed in the hope that it will be useful, but
+! WITHOUT ANY WARRANTY; without even the implied warranty of
+! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
+!
+! See the GNU Lesser General Public License in COPYING.LGPL for more details.
+
+! Computation kernels for the simulation
+
+MODULE mod_compute
+
+  USE mod_types
+  USE fstarpu_mod
+  USE iso_c_binding
+
+  IMPLICIT NONE
+
+CONTAINS
+
+  !--------------------------------------------------------------!
+  SUBROUTINE init_element(ro,dro,basis,Neq_max,Np,Ng,i)
+    INTEGER(KIND=C_INT),INTENT(IN)                           :: Neq_max,Np,Ng,i
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER,INTENT(INOUT) :: ro,basis,dro
+    !Local variables
+    INTEGER(KIND=C_INT)                                      :: n,nb,neq
+
+    DO nb=1,Np
+       DO neq= 1,Neq_max
+          ro(neq,nb)  = 0.01*(nb+neq)*i
+       END DO
+    END DO
+
+    DO nb=1,Np
+       DO neq= 1,Neq_max
+          dro(neq,nb) = 0.05*(nb-neq)*i
+       END DO
+    END DO
+
+    DO n=1,Ng
+       DO nb=1,Np
+          basis(nb,n) = 0.05*(n+nb)*i
+       END DO
+    END DO
+
+  END SUBROUTINE init_element
+
+  !--------------------------------------------------------------!
+  RECURSIVE SUBROUTINE loop_element_cpu_fortran(buffers, cl_args) BIND(C)
+    TYPE(C_PTR), VALUE, INTENT(IN)              :: buffers, cl_args
+
+    INTEGER(KIND=C_INT)                         :: Neq_max,Np,Ng
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER  :: ro,dro,basis
+    REAL(KIND=C_DOUBLE),TARGET                  :: coeff
+
+    Neq_max = fstarpu_matrix_get_ny(buffers, 0)
+    Np = fstarpu_matrix_get_ny(buffers, 2)
+    Ng = fstarpu_matrix_get_nx(buffers, 2)
+
+    CALL fstarpu_unpack_arg(cl_args,(/ c_loc(coeff) /))
+
+    CALL c_f_pointer(fstarpu_matrix_get_ptr(buffers, 0), ro, shape=[Neq_max,Np])
+    CALL c_f_pointer(fstarpu_matrix_get_ptr(buffers, 1), dro, shape=[Neq_max,Np])
+    CALL c_f_pointer(fstarpu_matrix_get_ptr(buffers, 2), basis, shape=[Np,Ng])
+
+    CALL loop_element_cpu(ro,dro,basis,coeff,Neq_max,Ng,Np)
+  END SUBROUTINE loop_element_cpu_fortran
+
+  !--------------------------------------------------------------!
+  RECURSIVE SUBROUTINE loop_element_cpu(ro,dro,basis,coeff,Neq_max,Ng,Np)
+    REAL(KIND=C_DOUBLE),INTENT(IN)                           :: coeff
+    INTEGER(KIND=C_INT),INTENT(IN)                           :: Neq_max,Ng,Np
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER,INTENT(IN)    :: ro,basis
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER,INTENT(INOUT) :: dro
+    !Local variables
+    REAL(KIND=C_DOUBLE)                                      :: coeff2,r
+    INTEGER(KIND=C_INT)                                      :: n,nb,neq
+
+    DO n=1,Ng
+       r = 0.
+       DO nb=1,Np
+          DO neq= 1,Neq_max
+             r = r + basis(nb,n) * ro(neq,nb)
+          ENDDO
+       ENDDO
+
+       coeff2 = r + coeff
+
+       DO nb=1,Np
+          DO neq = 1,Neq_max
+             dro(neq,nb) = coeff2 + dro(neq,nb)
+          ENDDO
+       ENDDO
+    ENDDO
+
+  END SUBROUTINE loop_element_cpu
+
+  !--------------------------------------------------------------!
+  RECURSIVE SUBROUTINE copy_element_cpu_fortran(buffers, cl_args) BIND(C)
+    TYPE(C_PTR), VALUE, INTENT(IN)              :: buffers, cl_args
+
+    INTEGER(KIND=C_INT)                         :: Neq_max,Np
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER  :: ro,dro
+
+    Neq_max = fstarpu_matrix_get_ny(buffers, 0)
+    Np = fstarpu_matrix_get_nx(buffers, 0)
+
+    CALL c_f_pointer(fstarpu_matrix_get_ptr(buffers, 0), ro, shape=[Neq_max,Np])
+    CALL c_f_pointer(fstarpu_matrix_get_ptr(buffers, 1), dro, shape=[Neq_max,Np])
+
+    CALL copy_element_cpu(ro,dro)
+
+  END SUBROUTINE copy_element_cpu_fortran
+
+  !--------------------------------------------------------------!
+  RECURSIVE SUBROUTINE copy_element_cpu(ro,dro)
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER,INTENT(INOUT) :: ro
+    REAL(KIND=C_DOUBLE),DIMENSION(:,:),POINTER,INTENT(IN)    :: dro
+
+    ro = ro + dro
+
+  END SUBROUTINE copy_element_cpu
+
+END MODULE mod_compute

+ 37 - 0
examples/native_fortran/mod_types.f90

@@ -0,0 +1,37 @@
+! StarPU --- Runtime system for heterogeneous multicore architectures.
+!
+! Copyright (C) 2015  ONERA
+! Copyright (C) 2015-2016  Inria
+!
+! StarPU is free software; you can redistribute it and/or modify
+! it under the terms of the GNU Lesser General Public License as published by
+! the Free Software Foundation; either version 2.1 of the License, or (at
+! your option) any later version.
+!
+! StarPU is distributed in the hope that it will be useful, but
+! WITHOUT ANY WARRANTY; without even the implied warranty of
+! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
+!
+! See the GNU Lesser General Public License in COPYING.LGPL for more details.
+
+MODULE mod_types
+
+  USE iso_c_binding
+
+  TYPE type_numpar
+     REAL(KIND=C_DOUBLE)                        :: coeff
+     INTEGER(KIND=C_INT)                        :: Neq_max
+  END TYPE type_numpar
+
+  TYPE type_mesh_elt
+     INTEGER(KIND=C_INT)                        :: Ng, Np
+     REAL(KIND=C_DOUBLE),POINTER,DIMENSION(:,:) :: ro, dro
+     REAL(KIND=C_DOUBLE),POINTER,DIMENSION(:,:) :: basis
+     TYPE(C_PTR)                                :: ro_h, dro_h, basis_h
+  END TYPE type_mesh_elt
+
+  TYPE type_mesh
+     TYPE(type_mesh_elt), POINTER, DIMENSION(:) :: elt
+  END TYPE type_mesh
+
+END MODULE mod_types

+ 176 - 0
examples/native_fortran/nf_example.f90

@@ -0,0 +1,176 @@
+! StarPU --- Runtime system for heterogeneous multicore architectures.
+!
+! Copyright (C) 2015  ONERA
+! Copyright (C) 2015-2016  Inria
+!
+! StarPU is free software; you can redistribute it and/or modify
+! it under the terms of the GNU Lesser General Public License as published by
+! the Free Software Foundation; either version 2.1 of the License, or (at
+! your option) any later version.
+!
+! StarPU is distributed in the hope that it will be useful, but
+! WITHOUT ANY WARRANTY; without even the implied warranty of
+! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
+!
+! See the GNU Lesser General Public License in COPYING.LGPL for more details.
+
+! This is an example of Fortran90 program making use of StarPU.
+! It registers a few matrices for each element of a domain, performs
+! update computations on them, and checks the result.
+
+PROGRAM f90_example
+
+  USE mod_types
+  USE fstarpu_mod
+  USE mod_compute
+  USE iso_c_binding
+
+  IMPLICIT NONE
+
+  TYPE(type_mesh)                :: mesh
+  TYPE(type_numpar),TARGET       :: numpar
+  TYPE(type_mesh_elt),POINTER    :: elt   => NULL()
+  INTEGER(KIND=C_INT)            :: i,Nelt,res,cpus
+  INTEGER(KIND=C_INT)            :: starpu_maj,starpu_min,starpu_rev
+  INTEGER(KIND=C_INT)            :: neq,ng,nb,it,it_tot
+  REAL(KIND=C_DOUBLE)            :: r, coeff2
+
+  TYPE(C_PTR) :: cl_loop_element = C_NULL_PTR ! loop codelet
+  TYPE(C_PTR) :: cl_copy_element = C_NULL_PTR ! copy codelet
+
+  !Initialization with arbitrary data
+  Nelt           = 2
+  it_tot = 2
+  numpar%Neq_max = 5
+  numpar%coeff   = 1.0
+  ALLOCATE(mesh%elt(Nelt))
+  DO i = 1,Nelt
+     elt => mesh%elt(i)
+     elt%Ng  = 4
+     elt%Np  = 2
+     ALLOCATE(elt%ro(numpar%Neq_max,elt%Np))
+     ALLOCATE(elt%dro(numpar%Neq_max,elt%Np))
+     ALLOCATE(elt%basis(elt%Np,elt%Ng))
+     CALL init_element(elt%ro,elt%dro,elt%basis,numpar%Neq_max,elt%Np,elt%Ng,i)
+  ENDDO
+
+  !Initialization of StarPU
+  res = fstarpu_init(C_NULL_PTR)
+  IF (res == -19) THEN
+     STOP 77
+  END IF
+  CALL fstarpu_get_version(starpu_maj,starpu_min,starpu_rev)
+  WRITE(6,'(a,i4,a,i4,a,i4)')      "StarPU version: ", starpu_maj , "." , starpu_min , "." , starpu_rev
+  cpus = fstarpu_cpu_worker_get_count()
+  IF (cpus == 0) THEN
+     CALL fstarpu_shutdown()
+     STOP 77
+  END IF
+
+  cl_loop_element = fstarpu_codelet_allocate()
+  CALL fstarpu_codelet_add_cpu_func(cl_loop_element, C_FUNLOC(loop_element_cpu_fortran))
+  CALL fstarpu_codelet_add_buffer(cl_loop_element, FSTARPU_R)
+  CALL fstarpu_codelet_add_buffer(cl_loop_element, FSTARPU_RW)
+  CALL fstarpu_codelet_add_buffer(cl_loop_element, FSTARPU_R)
+  ! TODO: add name "LOOP_ELEMENT"
+
+  cl_copy_element = fstarpu_codelet_allocate()
+  CALL fstarpu_codelet_add_cpu_func(cl_copy_element, C_FUNLOC(copy_element_cpu_fortran))
+  CALL fstarpu_codelet_add_buffer(cl_copy_element, FSTARPU_RW)
+  CALL fstarpu_codelet_add_buffer(cl_copy_element, FSTARPU_R)
+  ! TODO: add name "COPY_ELEMENT"
+
+  !Registration of elements
+  DO i = 1,Nelt
+     elt => mesh%elt(i)
+     elt%ro_h = fstarpu_matrix_data_register(c_loc(elt%ro), numpar%Neq_max, numpar%Neq_max, elt%Np, c_sizeof(elt%ro(1,1)), 0)
+     elt%dro_h = fstarpu_matrix_data_register(c_loc(elt%dro), numpar%Neq_max, numpar%Neq_max, elt%Np, c_sizeof(elt%dro(1,1)), 0)
+     elt%basis_h = fstarpu_matrix_data_register(c_loc(elt%basis), elt%Np, elt%Np, elt%Ng, c_sizeof(elt%basis(1,1)), 0)
+  ENDDO
+  !Compute
+  DO it = 1,it_tot
+
+     ! compute new dro for each element
+     DO i = 1,Nelt
+        elt => mesh%elt(i)
+        CALL fstarpu_insert_task((/ cl_loop_element,    &
+                FSTARPU_VALUE, c_loc(numpar%coeff), FSTARPU_SZ_REAL8, &
+                FSTARPU_DATA, elt%ro_h,                 &
+                FSTARPU_DATA, elt%dro_h,                &
+                FSTARPU_DATA, elt%basis_h,              &
+                C_NULL_PTR /))
+     ENDDO
+     ! sync (if needed by the algorithm)
+     CALL fstarpu_task_wait_for_all()
+
+     ! - - - - -
+
+     ! copy dro to ro for each element
+     DO i = 1,Nelt
+        elt => mesh%elt(i)
+        CALL fstarpu_insert_task((/ cl_copy_element,    &
+                FSTARPU_DATA, elt%ro_h,                 &
+                FSTARPU_DATA, elt%dro_h,                &
+                C_NULL_PTR /))
+     ENDDO
+     ! sync (if needed by the algorithm)
+     CALL fstarpu_task_wait_for_all()
+
+  ENDDO
+  !Unregistration of elements
+  DO i = 1,Nelt
+     elt => mesh%elt(i)
+     CALL fstarpu_data_unregister(elt%ro_h)
+     CALL fstarpu_data_unregister(elt%dro_h)
+     CALL fstarpu_data_unregister(elt%basis_h)
+  ENDDO
+
+  !Terminate StarPU, no task can be submitted after
+  CALL fstarpu_shutdown()
+
+  !Check data with StarPU
+  WRITE(6,'(a)') " "
+  WRITE(6,'(a)') " %%%% RESULTS STARPU %%%% "
+  WRITE(6,'(a)') " "
+  DO i = 1,Nelt
+     WRITE(6,'(a,i4,a)')      " elt ", i , "  ;  elt%ro = "
+     WRITE(6,'(10(1x,F11.2))') mesh%elt(i)%ro
+     WRITE(6,'(a)')           " ------------------------ "
+  ENDDO
+
+  !Same compute without StarPU
+  DO i = 1,Nelt
+     elt => mesh%elt(i)
+     CALL init_element(elt%ro,elt%dro,elt%basis,numpar%Neq_max,elt%Np,elt%Ng,i)
+  ENDDO
+
+  DO it = 1, it_tot
+     DO i = 1,Nelt
+        elt => mesh%elt(i)
+        CALL loop_element_cpu(elt%ro,elt%dro,elt%basis,numpar%coeff,numpar%Neq_max,elt%Ng,elt%Np)
+        elt%ro = elt%ro + elt%dro
+     ENDDO
+  ENDDO
+
+  WRITE(6,'(a)') " "
+  WRITE(6,'(a)') " %%%% RESULTS VERIFICATION %%%% "
+  WRITE(6,'(a)') " "
+
+  DO i = 1,Nelt
+     WRITE(6,'(a,i4,a)')      " elt ", i , "  ;  elt%ro = "
+     WRITE(6,'(10(1x,F11.2))') mesh%elt(i)%ro
+     WRITE(6,'(a)')           " ------------------------ "
+  ENDDO
+
+  WRITE(6,'(a)') " "
+
+  !Deallocation
+  DO i = 1,Nelt
+     elt => mesh%elt(i)
+     DEALLOCATE(elt%ro)
+     DEALLOCATE(elt%dro)
+     DEALLOCATE(elt%basis)
+  ENDDO
+  DEALLOCATE(mesh%elt)
+
+END PROGRAM f90_example

+ 1 - 1
examples/native_fortran/nf_matrix.f90

@@ -1,7 +1,7 @@
 program nf_matrix
         use iso_c_binding       ! C interfacing module
         use fstarpu_mod         ! StarPU interfacing module
-        use codelets
+        use nf_codelets
         implicit none
 
         real(8), dimension(:,:), allocatable, target :: ma

+ 37 - 0
include/fstarpu_mod.f90

@@ -23,7 +23,13 @@ module fstarpu_mod
         integer(c_int), bind(C) :: FSTARPU_REDUX
 
         type(c_ptr), bind(C) :: FSTARPU_DATA
+        type(c_ptr), bind(C) :: FSTARPU_VALUE
 
+        type(c_ptr), bind(C) :: FSTARPU_SZ_INT4
+        type(c_ptr), bind(C) :: FSTARPU_SZ_INT8
+
+        type(c_ptr), bind(C) :: FSTARPU_SZ_REAL4
+        type(c_ptr), bind(C) :: FSTARPU_SZ_REAL8
         interface
                 ! == starpu.h ==
 
@@ -215,6 +221,12 @@ module fstarpu_mod
                         type(c_ptr), dimension(:), intent(in) :: arglist
                 end subroutine fstarpu_insert_task
 
+                subroutine fstarpu_unpack_arg(cl_arg,bufferlist) bind(C)
+                        use iso_c_binding, only: c_ptr
+                        type(c_ptr), value, intent(in) :: cl_arg
+                        type(c_ptr), dimension(:), intent(in) :: bufferlist
+                end subroutine fstarpu_unpack_arg
+
         end interface
 
         contains
@@ -224,6 +236,11 @@ module fstarpu_mod
                         integer(c_int) :: fstarpu_init
                         type(c_ptr), value, intent(in) :: conf
 
+                        integer(4) :: FSTARPU_SZ_INT4_dummy
+                        integer(8) :: FSTARPU_SZ_INT8_dummy
+                        real(4) :: FSTARPU_SZ_REAL4_dummy
+                        real(8) :: FSTARPU_SZ_REAL8_dummy
+
                         ! Note: Referencing global C constants from Fortran has
                         ! been found unreliable on some architectures, notably
                         ! on Darwin. The get_integer/get_pointer_constant
@@ -258,6 +275,12 @@ module fstarpu_mod
                         FSTARPU_REDUX = fstarpu_get_integer_constant(C_CHAR_"FSTARPU_REDUX"//C_NULL_CHAR)
                         ! Initialize Fortran 'pointer' constants from C peers
                         FSTARPU_DATA = fstarpu_get_pointer_constant(C_CHAR_"FSTARPU_DATA"//C_NULL_CHAR)
+                        FSTARPU_VALUE = fstarpu_get_pointer_constant(C_CHAR_"FSTARPU_VALUE"//C_NULL_CHAR)
+                        ! Initialize size constants as 'c_ptr'
+                        FSTARPU_SZ_INT4 = transfer(int(c_sizeof(FSTARPU_SZ_INT4_dummy),kind=c_intptr_t),C_NULL_PTR)
+                        FSTARPU_SZ_INT8 = transfer(int(c_sizeof(FSTARPU_SZ_INT8_dummy),kind=c_intptr_t),C_NULL_PTR)
+                        FSTARPU_SZ_REAL4 = transfer(int(c_sizeof(FSTARPU_SZ_REAL4_dummy),kind=c_intptr_t),C_NULL_PTR)
+                        FSTARPU_SZ_REAL8 = transfer(int(c_sizeof(FSTARPU_SZ_REAL8_dummy),kind=c_intptr_t),C_NULL_PTR)
                         ! Initialize StarPU
                         if (c_associated(conf)) then 
                                 fstarpu_init = fstarpu_init_internal(conf)
@@ -265,4 +288,18 @@ module fstarpu_mod
                                 fstarpu_init = fstarpu_init_internal(C_NULL_PTR)
                         end if
                 end function fstarpu_init
+
+                function fstarpu_csizet_to_cptr(i) bind(C)
+                        use iso_c_binding
+                        type(c_ptr) :: fstarpu_csizet_to_cptr
+                        integer(c_size_t) :: i
+                        fstarpu_csizet_to_cptr = transfer(int(i,kind=c_intptr_t),C_NULL_PTR)
+                end function fstarpu_csizet_to_cptr
+
+                function fstarpu_int_to_cptr(i) bind(C)
+                        use iso_c_binding
+                        type(c_ptr) :: fstarpu_int_to_cptr
+                        integer :: i
+                        fstarpu_int_to_cptr = transfer(int(i,kind=c_intptr_t),C_NULL_PTR)
+                end function fstarpu_int_to_cptr
 end module fstarpu_mod

+ 59 - 1
src/util/fstarpu.c

@@ -30,6 +30,11 @@ static const int fstarpu_redux	= STARPU_REDUX;
 static const int _fstarpu_data = STARPU_R | STARPU_W | STARPU_SCRATCH | STARPU_REDUX;
 static const void * const fstarpu_data = &_fstarpu_data;
 
+static const int _fstarpu_value = STARPU_VALUE;
+static const void * const fstarpu_value = &_fstarpu_value;
+
+extern void _starpu_pack_arguments(size_t *current_offset, size_t *arg_buffer_size_, char **arg_buffer_, void *ptr, size_t ptr_size);
+
 int fstarpu_get_integer_constant(char *s)
 {
 	if	(!strcmp(s, "FSTARPU_R"))	{ return fstarpu_r; }
@@ -43,6 +48,7 @@ int fstarpu_get_integer_constant(char *s)
 const void *fstarpu_get_pointer_constant(char *s)
 {
 	if (!strcmp(s, "FSTARPU_DATA")) { return fstarpu_data; }
+	if (!strcmp(s, "FSTARPU_VALUE")) { return fstarpu_value; }
 	else { _FSTARPU_ERROR("unknown pointer constant"); }
 }
 
@@ -168,11 +174,42 @@ int fstarpu_matrix_get_nx(void *buffers[], int i)
 	return STARPU_MATRIX_GET_NY(buffers[i]);
 }
 
+void fstarpu_unpack_arg(char *cl_arg, void ***_buffer_list)
+{
+	void **buffer_list = *_buffer_list;
+	size_t current_arg_offset = 0;
+	int nargs, arg;
+
+	/* We fill the different pointers with the appropriate arguments */
+	memcpy(&nargs, cl_arg, sizeof(nargs));
+	current_arg_offset += sizeof(nargs);
+
+	for (arg = 0; arg < nargs; arg++)
+	{
+		void *argptr = buffer_list[arg];
+
+		/* If not reading all cl_args */
+		if(argptr == NULL)
+			break;
+
+		size_t arg_size;
+		memcpy(&arg_size, cl_arg+current_arg_offset, sizeof(arg_size));
+		current_arg_offset += sizeof(arg_size);
+
+		memcpy(argptr, cl_arg+current_arg_offset, arg_size);
+		current_arg_offset += arg_size;
+	}
+}
+
 void fstarpu_insert_task(void ***_arglist)
 {
 	void **arglist = *_arglist;
 	int i = 0;
+	char *arg_buffer_ = NULL;
+	size_t arg_buffer_size_ = 0;
+	size_t current_offset = sizeof(int);
 	int current_buffer = 0;
+	int nargs = 0;
 	struct starpu_task *task = NULL;
 	struct starpu_codelet *cl = arglist[i++];
 	if (cl == NULL)
@@ -182,7 +219,6 @@ void fstarpu_insert_task(void ***_arglist)
 	task = starpu_task_create();
 	task->cl = cl;
 	task->name = NULL;
-	task->cl_arg_free = 0;
 	while (arglist[i] != NULL)
 	{
 		if (arglist[i] == fstarpu_data)
@@ -200,12 +236,34 @@ void fstarpu_insert_task(void ***_arglist)
 			}
 			current_buffer++;
 		}
+		else if (arglist[i] == fstarpu_value)
+		{
+			i++;
+			void *ptr = arglist[i];
+			i++;
+			size_t ptr_size = (size_t)(intptr_t)arglist[i];
+			nargs++;
+			_starpu_pack_arguments(&current_offset, &arg_buffer_size_, &arg_buffer_, ptr, ptr_size);
+		}
 		else
 		{
 			_FSTARPU_ERROR("unknown/unsupported argument type");
 		}
 		i++;
 	}
+
+	if (nargs)
+	{
+		memcpy(arg_buffer_, (int *)&nargs, sizeof(nargs));
+		task->cl_arg = arg_buffer_;
+		task->cl_arg_size = arg_buffer_size_;
+	}
+	else
+	{
+		free(arg_buffer_);
+		arg_buffer_ = NULL;
+	}
+
 	int ret = starpu_task_submit(task);
 	if (ret != 0)
 	{

+ 0 - 1
src/util/starpu_task_insert_utils.c

@@ -21,7 +21,6 @@
 #include <common/utils.h>
 #include <core/task.h>
 
-static
 void _starpu_pack_arguments(size_t *current_offset, size_t *arg_buffer_size_, char **arg_buffer_, void *ptr, size_t ptr_size)
 {
 	if (*current_offset + sizeof(ptr_size) + ptr_size > *arg_buffer_size_)