PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
plasma_f90.f90
Go to the documentation of this file.
1 !
2 ! Copyright © 2011 The Numerical Algorithms Group Ltd. All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions are
6 ! met:
7 ! - Redistributions of source code must retain the above copyright notice,
8 ! this list of conditions, and the following disclaimer.
9 ! - Redistributions in binary form must reproduce the above copyright
10 ! notice, this list of conditions and the following disclaimer listed in
11 ! this license in the documentation and/or other materials provided with
12 ! the distribution.
13 ! - Neither the name of the copyright holders nor the names of its
14 ! contributors may be used to endorse or promote products derived from
15 ! this software without specific prior written permission.
16 !
17 ! This software is provided by the copyright holders and contributors "as
18 ! is" and any express or implied warranties, including, but not limited
19 ! to, the implied warranties of merchantability and fitness for a
20 ! particular purpose are disclaimed. in no event shall the copyright owner
21 ! or contributors be liable for any direct, indirect, incidental, special,
22 ! exemplary, or consequential damages (including, but not limited to,
23 ! procurement of substitute goods or services; loss of use, data, or
24 ! profits; or business interruption) however caused and on any theory of
25 ! liability, whether in contract, strict liability, or tort (including
26 ! negligence or otherwise) arising in any way out of the use of this
27 ! software, even if advised of the possibility of such damage.
28 !
29 !
30 ! @file plasma_f90.f90
31 !
32 ! PLASMA fortran 90 interface
33 ! PLASMA is a software package provided by Univ. of Tennessee,
34 ! Univ. of California Berkeley and Univ. of Colorado Denver
35 !
36 ! @version 2.4.2
37 ! @author Numerical Algorithm Group
38 ! @date 2011-09-15
39 ! @precisions normal z -> c d s
40 !
41 module plasma
42 
43  use plasma_s
44  use plasma_d
45  use plasma_c
46  use plasma_z
47  include 'plasmaf.h'
48 
49  logical :: plasma_initialized = .false.
50  integer, parameter :: sp = kind(0.0)
51  integer, parameter :: dp = kind(0.0d0)
52 
53 ! private PLASMA_Init_c
54  interface
55  function plasma_init_c(cores) &
56  & bind(c, name='PLASMA_Init')
57  use iso_c_binding
58  implicit none
59  integer(kind=c_int) :: plasma_init_c
60  integer(kind=c_int), value :: cores
61  end function plasma_init_c
62  end interface
63 
64 ! private PLASMA_Finalize_c
65  interface
66  function plasma_finalize_c() &
67  & bind(c, name='PLASMA_Finalize')
68  use iso_c_binding
69  implicit none
70  integer(kind=c_int) :: plasma_finalize_c
71  end function plasma_finalize_c
72  end interface
73 
74 ! private PLASMA_Set_c
75  interface
76  function plasma_set_c(param, pval) &
77  & bind(c, name='PLASMA_Set')
78  use iso_c_binding
79  implicit none
80  integer(kind=c_int) :: plasma_set_c
81  integer(kind=c_int), value :: param
82  integer(kind=c_int), value :: pval
83  end function plasma_set_c
84  end interface
85 
86 ! private PLASMA_Get_c
87  interface
88  function plasma_get_c(param, pval) &
89  & bind(c, name='PLASMA_Get')
90  use iso_c_binding
91  implicit none
92  integer(kind=c_int) :: plasma_get_c
93  integer(kind=c_int), value :: param
94  type(c_ptr), value :: pval
95  end function plasma_get_c
96  end interface
97 
98 ! private PLASMA_Enable_c
99  interface
100  function plasma_enable_c(param) &
101  & bind(c, name='PLASMA_Enable')
102  use iso_c_binding
103  implicit none
104  integer(kind=c_int) :: plasma_enable_c
105  integer(kind=c_int), value :: param
106  end function plasma_enable_c
107  end interface
108 
109 ! private PLASMA_Disable_c
110  interface
111  function plasma_disable_c(param) &
112  & bind(c, name='PLASMA_Disable')
113  use iso_c_binding
114  implicit none
115  integer(kind=c_int) :: plasma_disable_c
116  integer(kind=c_int), value :: param
117  end function plasma_disable_c
118  end interface
119 
120 ! private PLASMA_Lapack_to_Tile_c
121  interface
122  function plasma_lapack_to_tile_c(a_lpk,lda,a_pma) &
123  & bind(c, name='PLASMA_Lapack_to_Tile')
124  use iso_c_binding
125  integer(kind=c_int) :: plasma_lapack_to_tile_c
126  type(c_ptr), value :: a_lpk, a_pma
127  integer(kind=c_int), value :: lda
128  end function plasma_lapack_to_tile_c
129  end interface
130 
131 ! private PLASMA_Tile_to_Lapack_c
132  interface
133  function plasma_tile_to_lapack_c(a_pma,a_lpk,lda) &
134  & bind(c, name='PLASMA_Tile_to_Lapack')
135  use iso_c_binding
136  integer(kind=c_int) :: plasma_tile_to_lapack_c
137  type(c_ptr), value :: a_lpk, a_pma
138  integer(kind=c_int), value :: lda
139  end function plasma_tile_to_lapack_c
140  end interface
141 
142 ! private PLASMA_Desc_Create_c
143  interface
144  function plasma_desc_create_c(desc, mat, dtyp, mb, nb, bsiz, lm, ln, i, j, m, n) &
145  & bind(c, name='PLASMA_Desc_Create')
146  use iso_c_binding
147  integer(kind=c_int) :: plasma_desc_create_c
148  type(c_ptr) :: desc
149  type(c_ptr), value :: mat
150  integer(kind=c_int), value :: dtyp
151  integer(kind=c_int), value :: mb, nb, bsiz, lm, ln, i, j, m, n
152  end function plasma_desc_create_c
153  end interface
154 
155 ! private PLASMA_Desc_Destroy_c
156  interface
157  function plasma_desc_destroy_c(desc) &
158  & bind(c, name='PLASMA_Desc_Destroy')
159  use iso_c_binding
160  integer(kind=c_int) :: plasma_desc_destroy_c
161  type(c_ptr) :: desc
162  end function plasma_desc_destroy_c
163  end interface
164 
165 ! private free_c
166  interface
167  subroutine free_c(ptr) bind(c, name='free')
168  use iso_c_binding
169  implicit none
170  type(c_ptr), value :: ptr
171  end subroutine free_c
172  end interface
173 
175  module procedure plasma_lapack_to_tile_s
176  module procedure plasma_lapack_to_tile_d
177  module procedure plasma_lapack_to_tile_cpx
178  module procedure plasma_lapack_to_tile_z
179  end interface plasma_lapack_to_tile
180 
182  module procedure plasma_tile_to_lapack_s
183  module procedure plasma_tile_to_lapack_d
184  module procedure plasma_tile_to_lapack_cpx
185  module procedure plasma_tile_to_lapack_z
186  end interface plasma_tile_to_lapack
187 
189  module procedure plasma_desc_create_s
190  module procedure plasma_desc_create_d
191  module procedure plasma_desc_create_cpx
192  module procedure plasma_desc_create_z
193  end interface plasma_desc_create
194 
195  contains
196 
197  subroutine plasma_init(ncores,info)
198  use iso_c_binding
199  implicit none
200  integer(kind=c_int), intent(in) :: ncores
201  integer(kind=c_int), intent(out) :: info
202  info = plasma_init_c(ncores)
203  plasma_initialized = .true.
204  end subroutine plasma_init
205 
206  subroutine plasma_finalize(info)
207  use iso_c_binding
208  implicit none
209  integer(kind=c_int), intent(out) :: info
210  info = plasma_finalize_c()
211  plasma_initialized = .false.
212  end subroutine plasma_finalize
213 
214  subroutine plasma_set(param,pval,info)
215  use iso_c_binding
216  implicit none
217  integer(kind=c_int), intent(in) :: param
218  integer(kind=c_int), intent(in) :: pval
219  integer(kind=c_int), intent(out) :: info
220  info = plasma_set_c(param,pval)
221  end subroutine plasma_set
222 
223  subroutine plasma_get(param,pval,info)
224  use iso_c_binding
225  implicit none
226  integer(kind=c_int), intent(in) :: param
227  integer(kind=c_int), intent(out), target :: pval
228  integer(kind=c_int), intent(out) :: info
229  info = plasma_get_c(param,c_loc(pval))
230  end subroutine plasma_get
231 
232  subroutine plasma_enable(param,info)
233  use iso_c_binding
234  implicit none
235  integer(kind=c_int), intent(in) :: param
236  integer(kind=c_int), intent(out) :: info
237  info = plasma_enable_c(param)
238  end subroutine plasma_enable
239 
240  subroutine plasma_disable(param,info)
241  use iso_c_binding
242  implicit none
243  integer(kind=c_int), intent(in) :: param
244  integer(kind=c_int), intent(out) :: info
245  info = plasma_disable_c(param)
246  end subroutine plasma_disable
247 
248 ! overloaded: single precision
249  subroutine plasma_lapack_to_tile_s(a_lpk,lda,a_pma,info)
250  use iso_c_binding
251  implicit none
252  integer(kind=c_int), intent(in) :: lda
253  real(kind=sp), intent(out), target :: a_lpk(lda,*)
254  type(c_ptr), intent(out) :: a_pma
255  integer(kind=c_int), intent(out) :: info
256  info = plasma_lapack_to_tile_c(c_loc(a_lpk),lda,a_pma)
257  end subroutine plasma_lapack_to_tile_s
258 ! overloaded: double precision
259  subroutine plasma_lapack_to_tile_d(a_lpk,lda,a_pma,info)
260  use iso_c_binding
261  implicit none
262  integer(kind=c_int), intent(in) :: lda
263  real(kind=dp), intent(out), target :: a_lpk(lda,*)
264  type(c_ptr), intent(out) :: a_pma
265  integer(kind=c_int), intent(out) :: info
266  info = plasma_lapack_to_tile_c(c_loc(a_lpk),lda,a_pma)
267  end subroutine plasma_lapack_to_tile_d
268 ! overloaded: single precision complex
269  subroutine plasma_lapack_to_tile_cpx(a_lpk,lda,a_pma,info)
270  use iso_c_binding
271  implicit none
272  integer(kind=c_int), intent(in) :: lda
273  complex(kind=sp), intent(out), target :: a_lpk(lda,*)
274  type(c_ptr), intent(out) :: a_pma
275  integer(kind=c_int), intent(out) :: info
276  info = plasma_lapack_to_tile_c(c_loc(a_lpk),lda,a_pma)
277  end subroutine plasma_lapack_to_tile_cpx
278 ! overloaded: double precision complex
279  subroutine plasma_lapack_to_tile_z(a_lpk,lda,a_pma,info)
280  use iso_c_binding
281  implicit none
282  integer(kind=c_int), intent(in) :: lda
283  complex(kind=dp), intent(out), target :: a_lpk(lda,*)
284  type(c_ptr), intent(out) :: a_pma
285  integer(kind=c_int), intent(out) :: info
286  info = plasma_lapack_to_tile_c(c_loc(a_lpk),lda,a_pma)
287  end subroutine plasma_lapack_to_tile_z
288 
289 ! overloaded: single precision
290  subroutine plasma_tile_to_lapack_s(a_pma,a_lpk,lda,info)
291  use iso_c_binding
292  implicit none
293  integer(kind=c_int), intent(in) :: lda
294  real(kind=sp), intent(out), target :: a_lpk(lda,*)
295  type(c_ptr), intent(in) :: a_pma
296  integer(kind=c_int), intent(out) :: info
297  info = plasma_tile_to_lapack_c(a_pma,c_loc(a_lpk),lda)
298  end subroutine plasma_tile_to_lapack_s
299 ! overloaded: double precision
300  subroutine plasma_tile_to_lapack_d(a_pma,a_lpk,lda,info)
301  use iso_c_binding
302  implicit none
303  integer(kind=c_int), intent(in) :: lda
304  real(kind=dp), intent(out), target :: a_lpk(lda,*)
305  type(c_ptr), intent(in) :: a_pma
306  integer(kind=c_int), intent(out) :: info
307  info = plasma_tile_to_lapack_c(a_pma,c_loc(a_lpk),lda)
308  end subroutine plasma_tile_to_lapack_d
309 ! overloaded: single precision complex
310  subroutine plasma_tile_to_lapack_cpx(a_pma,a_lpk,lda,info)
311  use iso_c_binding
312  implicit none
313  integer(kind=c_int), intent(in) :: lda
314  complex(kind=sp), intent(out), target :: a_lpk(lda,*)
315  type(c_ptr), intent(in) :: a_pma
316  integer(kind=c_int), intent(out) :: info
317  info = plasma_tile_to_lapack_c(a_pma,c_loc(a_lpk),lda)
318  end subroutine plasma_tile_to_lapack_cpx
319 ! overloaded: double precision complex
320  subroutine plasma_tile_to_lapack_z(a_pma,a_lpk,lda,info)
321  use iso_c_binding
322  implicit none
323  integer(kind=c_int), intent(in) :: lda
324  complex(kind=dp), intent(out), target :: a_lpk(lda,*)
325  type(c_ptr), intent(in) :: a_pma
326  integer(kind=c_int), intent(out) :: info
327  info = plasma_tile_to_lapack_c(a_pma,c_loc(a_lpk),lda)
328  end subroutine plasma_tile_to_lapack_z
329 
330 ! overloaded: single precision
331  subroutine plasma_desc_create_s(desc,mat,dtyp,mb,nb,bsiz,lm,ln,i,j,m,n,info)
332  use iso_c_binding
333  implicit none
334  type(c_ptr), intent(out) :: desc
335  integer(kind=c_int), intent(in) :: mb, nb, bsiz, lm, ln, i, j, m, n
336  real(kind=sp), intent(in), target :: mat(lm,*)
337  integer(kind=c_int), intent(in) :: dtyp
338  integer(kind=c_int), intent(out) :: info
339  info = plasma_desc_create_c(desc,c_loc(mat),dtyp,mb,nb,bsiz,lm,ln,i,j,m,n)
340  end subroutine plasma_desc_create_s
341 ! overloaded: double precision
342  subroutine plasma_desc_create_d(desc,mat,dtyp,mb,nb,bsiz,lm,ln,i,j,m,n,info)
343  use iso_c_binding
344  implicit none
345  type(c_ptr), intent(out) :: desc
346  integer(kind=c_int), intent(in) :: mb, nb, bsiz, lm, ln, i, j, m, n
347  real(kind=dp), intent(in), target :: mat(lm,*)
348  integer(kind=c_int), intent(in) :: dtyp
349  integer(kind=c_int), intent(out) :: info
350  info = plasma_desc_create_c(desc,c_loc(mat),dtyp,mb,nb,bsiz,lm,ln,i,j,m,n)
351  end subroutine plasma_desc_create_d
352 ! overloaded: single precision complex
353  subroutine plasma_desc_create_cpx(desc,mat,dtyp,mb,nb,bsiz,lm,ln,i,j,m,n,info)
354  use iso_c_binding
355  implicit none
356  type(c_ptr), intent(out) :: desc
357  integer(kind=c_int), intent(in) :: mb, nb, bsiz, lm, ln, i, j, m, n
358  complex(kind=sp), intent(in), target :: mat(lm,*)
359  integer(kind=c_int), intent(in) :: dtyp
360  integer(kind=c_int), intent(out) :: info
361  info = plasma_desc_create_c(desc,c_loc(mat),dtyp,mb,nb,bsiz,lm,ln,i,j,m,n)
362  end subroutine plasma_desc_create_cpx
363 ! overloaded: double precision complex
364  subroutine plasma_desc_create_z(desc,mat,dtyp,mb,nb,bsiz,lm,ln,i,j,m,n,info)
365  use iso_c_binding
366  implicit none
367  type(c_ptr), intent(out) :: desc
368  integer(kind=c_int), intent(in) :: mb, nb, bsiz, lm, ln, i, j, m, n
369  complex(kind=dp), intent(in), target :: mat(lm,*)
370  integer(kind=c_int), intent(in) :: dtyp
371  integer(kind=c_int), intent(out) :: info
372  info = plasma_desc_create_c(desc,c_loc(mat),dtyp,mb,nb,bsiz,lm,ln,i,j,m,n)
373  end subroutine plasma_desc_create_z
374 
375  subroutine plasma_desc_destroy(desc,info)
376  use iso_c_binding
377  implicit none
378  type(c_ptr), intent(inout) :: desc
379  real(kind=sp), intent(out) :: info
380  info = plasma_desc_destroy_c(desc)
381  end subroutine plasma_desc_destroy
382 
383  subroutine plasma_free(ptr)
384  use iso_c_binding
385  implicit none
386  type(c_ptr), intent(in) :: ptr
387  call free_c(ptr)
388  end subroutine plasma_free
389 
390 end module plasma