ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_blas.F90
Go to the documentation of this file.
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!! !!!!
4!!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5!!!! !!!!
6!!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7!!!! !!!!
8!!!! This file is part of the ParaMonte library. !!!!
9!!!! !!!!
10!!!! LICENSE !!!!
11!!!! !!!!
12!!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13!!!! !!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16
32
33!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
34
35module pm_blas
36
37 use pm_kind, only: SK, IK, LK, RKS, RKD
38
39 implicit none
40
41 character(*, SK), parameter :: MODULE_NAME = "@pm_blas"
42
43!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
44
60 interface blasSWAP
61 pure subroutine zswap(n, x, incx, y, incy)
62 use pm_kind , only: IK, TKG => RKD
63 integer(IK) , intent(in) :: n, incx, incy
64 complex(TKG), intent(inout) :: x(*), y(*)
65 end subroutine zswap
66 pure subroutine cswap(n, x, incx, y, incy)
67 use pm_kind , only: IK, TKG => RKS
68 integer(IK) , intent(in) :: n, incx, incy
69 complex(TKG), intent(inout) :: x(*), y(*)
70 end subroutine cswap
71 pure subroutine dswap(n, x, incx, y, incy)
72 use pm_kind , only: IK, TKG => RKD
73 integer(IK) , intent(in) :: n, incx, incy
74 real(TKG) , intent(inout) :: x(*), y(*)
75 end subroutine dswap
76 pure subroutine sswap(n, x, incx, y, incy)
77 use pm_kind , only: IK, TKG => RKS
78 integer(IK) , intent(in) :: n, incx, incy
79 real(TKG) , intent(inout) :: x(*), y(*)
80 end subroutine sswap
81 end interface
82
83!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
84
102 interface blasGEMV
103 pure subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
104 use pm_kind , only: IK, TKG => RKD
105 implicit none
106 complex(TKG), intent(inout) :: y
107 complex(TKG), intent(in) :: alpha, beta, a, x
108 integer(IK) , intent(in) :: m, n, lda, incx, incy
109 character , intent(in) :: trans
110 end subroutine
111 pure subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
112 use pm_kind , only: IK, TKG => RKS
113 implicit none
114 complex(TKG), intent(inout) :: y
115 complex(TKG), intent(in) :: alpha, beta, a, x
116 integer(IK) , intent(in) :: m, n, lda, incx, incy
117 character , intent(in) :: trans
118 end subroutine
119 pure subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
120 use pm_kind , only: IK, TKG => RKD
121 implicit none
122 real(TKG) , intent(inout) :: y
123 real(TKG) , intent(in) :: alpha, beta, a, x
124 integer(IK) , intent(in) :: m, n, lda, incx, incy
125 character , intent(in) :: trans
126 end subroutine
127 pure subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
128 use pm_kind , only: IK, TKG => RKS
129 implicit none
130 real(TKG) , intent(inout) :: y
131 real(TKG) , intent(in) :: alpha, beta, a, x
132 integer(IK) , intent(in) :: m, n, lda, incx, incy
133 character , intent(in) :: trans
134 end subroutine
135 end interface
136
137!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
138
156 interface blasHPMV
157 pure subroutine zhpmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
158 use pm_kind , only: IK, TKG => RKD
159 implicit none
160 complex(TKG), intent(inout) :: y
161 complex(TKG), intent(in) :: alpha, beta, ap, x
162 integer(IK) , intent(in) :: n, incx, incy
163 character , intent(in) :: uplo
164 end subroutine
165 pure subroutine chpmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
166 use pm_kind , only: IK, TKG => RKS
167 implicit none
168 complex(TKG), intent(inout) :: y
169 complex(TKG), intent(in) :: alpha, beta, ap, x
170 integer(IK) , intent(in) :: n, incx, incy
171 character , intent(in) :: uplo
172 end subroutine
173 end interface
174
175!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
176
194 interface blasSPMV
195 pure subroutine dspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
196 use pm_kind , only: IK, TKG => RKD
197 implicit none
198 real(TKG) , intent(inout) :: y
199 real(TKG) , intent(in) :: alpha, beta, ap, x
200 integer(IK) , intent(in) :: n, incx, incy
201 character , intent(in) :: uplo
202 end subroutine
203 pure subroutine sspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
204 use pm_kind , only: IK, TKG => RKS
205 implicit none
206 real(TKG) , intent(inout) :: y
207 real(TKG) , intent(in) :: alpha, beta, ap, x
208 integer(IK) , intent(in) :: n, incx, incy
209 character , intent(in) :: uplo
210 end subroutine
211 end interface
212
213!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
214
232 interface blasHEMV
233 pure subroutine zhemv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)
234 use pm_kind , only: IK, TKG => RKD
235 implicit none
236 complex(TKG), intent(inout) :: y
237 complex(TKG), intent(in) :: alpha, beta, a, x
238 integer(IK) , intent(in) :: n, lda, incx, incy
239 character , intent(in) :: uplo
240 end subroutine
241 pure subroutine chemv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)
242 use pm_kind , only: IK, TKG => RKS
243 implicit none
244 complex(TKG), intent(inout) :: y
245 complex(TKG), intent(in) :: alpha, beta, a, x
246 integer(IK) , intent(in) :: n, lda, incx, incy
247 character , intent(in) :: uplo
248 end subroutine
249 end interface
250
251!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
252
270 interface blasSYMV
271 pure subroutine dsymv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)
272 use pm_kind , only: IK, TKG => RKD
273 implicit none
274 real(TKG) , intent(inout) :: y
275 real(TKG) , intent(in) :: alpha, beta, a, x
276 integer(IK) , intent(in) :: n, lda, incx, incy
277 character , intent(in) :: uplo
278 end subroutine
279 pure subroutine ssymv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)
280 use pm_kind , only: IK, TKG => RKS
281 implicit none
282 real(TKG) , intent(inout) :: y
283 real(TKG) , intent(in) :: alpha, beta, a, x
284 integer(IK) , intent(in) :: n, lda, incx, incy
285 character , intent(in) :: uplo
286 end subroutine
287 end interface
288
289!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
290
307 interface blasSYMM
308 pure subroutine zsymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
309 use pm_kind , only: IK, TKG => RKD
310 implicit none
311 complex(TKG), intent(inout) :: c
312 complex(TKG), intent(in) :: alpha, beta, a, b
313 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
314 character , intent(in) :: side, uplo
315 end subroutine
316 pure subroutine csymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
317 use pm_kind , only: IK, TKG => RKS
318 implicit none
319 complex(TKG), intent(inout) :: c
320 complex(TKG), intent(in) :: alpha, beta, a, b
321 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
322 character , intent(in) :: side, uplo
323 end subroutine
324 pure subroutine dsymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
325 use pm_kind , only: IK, TKG => RKD
326 implicit none
327 real(TKG) , intent(inout) :: c
328 real(TKG) , intent(in) :: alpha, beta, a, b
329 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
330 character , intent(in) :: side, uplo
331 end subroutine
332 pure subroutine ssymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
333 use pm_kind , only: IK, TKG => RKS
334 implicit none
335 real(TKG) , intent(inout) :: c
336 real(TKG) , intent(in) :: alpha, beta, a, b
337 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
338 character , intent(in) :: side, uplo
339 end subroutine
340 end interface
341
342!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
343
360 interface blasHEMM
361 pure subroutine zhemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
362 use pm_kind , only: IK, TKG => RKD
363 implicit none
364 complex(TKG), intent(inout) :: c
365 complex(TKG), intent(in) :: alpha, beta, a, b
366 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
367 character , intent(in) :: side, uplo
368 end subroutine zhemm
369 pure subroutine chemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
370 use pm_kind , only: IK, TKG => RKS
371 implicit none
372 complex(TKG), intent(inout) :: c
373 complex(TKG), intent(in) :: alpha, beta, a, b
374 integer(IK) , intent(in) :: m, n, lda, ldb, ldc
375 character , intent(in) :: side, uplo
376 end subroutine chemm
377 end interface
378
379!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
380
396 interface blasGEMM
397 pure subroutine zgemm(transa, transb, l, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
398 use pm_kind , only: IK, TKG => RKD
399 implicit none
400 complex(TKG), intent(inout) :: c
401 complex(TKG), intent(in) :: alpha, beta, a, b
402 integer(IK) , intent(in) :: l, m, n, lda, ldb, ldc
403 character , intent(in) :: transa, transb
404 end subroutine
405 pure subroutine cgemm(transa, transb, l, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
406 use pm_kind , only: IK, TKG => RKS
407 implicit none
408 complex(TKG), intent(inout) :: c
409 complex(TKG), intent(in) :: alpha, beta, a, b
410 integer(IK) , intent(in) :: l, m, n, lda, ldb, ldc
411 character , intent(in) :: transa, transb
412 end subroutine
413 pure subroutine dgemm(transa, transb, l, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
414 use pm_kind , only: IK, TKG => RKD
415 implicit none
416 real(TKG) , intent(inout) :: c
417 real(TKG) , intent(in) :: alpha, beta, a, b
418 integer(IK) , intent(in) :: l, m, n, lda, ldb, ldc
419 character , intent(in) :: transa, transb
420 end subroutine
421 pure subroutine sgemm(transa, transb, l, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
422 use pm_kind , only: IK, TKG => RKS
423 implicit none
424 real(TKG) , intent(inout) :: c
425 real(TKG) , intent(in) :: alpha, beta, a, b
426 integer(IK) , intent(in) :: l, m, n, lda, ldb, ldc
427 character , intent(in) :: transa, transb
428 end subroutine
429 end interface
430
431!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
432
453 interface blasTRMV
454 pure subroutine ztrmv(uplo, trans, diag, n, a, lda, x, incx)
455 use pm_kind , only: IK, TKG => RKD
456 implicit none
457 character , intent(in) :: uplo, trans, diag
458 integer(IK) , intent(in) :: n, lda, incx
459 complex(TKG), intent(inout) :: x
460 complex(TKG), intent(in) :: a
461 end subroutine
462 pure subroutine ctrmv(uplo, trans, diag, n, a, lda, x, incx)
463 use pm_kind , only: IK, TKG => RKS
464 implicit none
465 character , intent(in) :: uplo, trans, diag
466 integer(IK) , intent(in) :: n, lda, incx
467 complex(TKG), intent(inout) :: x
468 complex(TKG), intent(in) :: a
469 end subroutine
470 pure subroutine dtrmv(uplo, trans, diag, n, a, lda, x, incx)
471 use pm_kind , only: IK, TKG => RKD
472 implicit none
473 character , intent(in) :: uplo, trans, diag
474 integer(IK) , intent(in) :: n, lda, incx
475 real(TKG) , intent(inout) :: x
476 real(TKG) , intent(in) :: a
477 end subroutine
478 pure subroutine strmv(uplo, trans, diag, n, a, lda, x, incx)
479 use pm_kind , only: IK, TKG => RKS
480 implicit none
481 character , intent(in) :: uplo, trans, diag
482 integer(IK) , intent(in) :: n, lda, incx
483 real(TKG) , intent(inout) :: x
484 real(TKG) , intent(in) :: a
485 end subroutine
486 end interface
487
488!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
489
510 interface blasTRSV
511 pure subroutine ztrsv(uplo, trans, diag, n, a, lda, x, incx)
512 use pm_kind , only: IK, TKG => RKD
513 implicit none
514 character , intent(in) :: uplo, trans, diag
515 integer(IK) , intent(in) :: n, lda, incx
516 complex(TKG), intent(inout) :: x
517 complex(TKG), intent(in) :: a
518 end subroutine
519 pure subroutine ctrsv(uplo, trans, diag, n, a, lda, x, incx)
520 use pm_kind , only: IK, TKG => RKS
521 implicit none
522 character , intent(in) :: uplo, trans, diag
523 integer(IK) , intent(in) :: n, lda, incx
524 complex(TKG), intent(inout) :: x
525 complex(TKG), intent(in) :: a
526 end subroutine
527 pure subroutine dtrsv(uplo, trans, diag, n, a, lda, x, incx)
528 use pm_kind , only: IK, TKG => RKD
529 implicit none
530 character , intent(in) :: uplo, trans, diag
531 integer(IK) , intent(in) :: n, lda, incx
532 real(TKG) , intent(inout) :: x
533 real(TKG) , intent(in) :: a
534 end subroutine
535 pure subroutine strsv(uplo, trans, diag, n, a, lda, x, incx)
536 use pm_kind , only: IK, TKG => RKS
537 implicit none
538 character , intent(in) :: uplo, trans, diag
539 integer(IK) , intent(in) :: n, lda, incx
540 real(TKG) , intent(inout) :: x
541 real(TKG) , intent(in) :: a
542 end subroutine
543 end interface
544
545!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
546
563 interface blasTRMM
564 pure subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
565 use pm_kind , only: IK, TKG => RKD
566 implicit none
567 integer(IK) , intent(in) :: m, n, lda, ldb
568 complex(TKG), intent(inout) :: b
569 complex(TKG), intent(in) :: alpha, a
570 character , intent(in) :: side, uplo, transa, diag
571 end subroutine
572 pure subroutine ctrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
573 use pm_kind , only: IK, TKG => RKS
574 implicit none
575 integer(IK) , intent(in) :: m, n, lda, ldb
576 complex(TKG), intent(inout) :: b
577 complex(TKG), intent(in) :: alpha, a
578 character , intent(in) :: side, uplo, transa, diag
579 end subroutine
580 pure subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
581 use pm_kind , only: IK, TKG => RKD
582 implicit none
583 integer(IK) , intent(in) :: m, n, lda, ldb
584 real(TKG) , intent(inout) :: b
585 real(TKG) , intent(in) :: alpha, a
586 character , intent(in) :: side, uplo, transa, diag
587 end subroutine
588 pure subroutine strmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
589 use pm_kind , only: IK, TKG => RKS
590 implicit none
591 integer(IK) , intent(in) :: m, n, lda, ldb
592 real(TKG) , intent(inout) :: b
593 real(TKG) , intent(in) :: alpha, a
594 character , intent(in) :: side, uplo, transa, diag
595 end subroutine
596 end interface
597
598!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
599
616 interface blasTRSM
617 pure subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
618 use pm_kind , only: IK, TKG => RKD
619 implicit none
620 integer(IK) , intent(in) :: m, n, lda, ldb
621 complex(TKG), intent(inout) :: b
622 complex(TKG), intent(in) :: alpha, a
623 character , intent(in) :: side, uplo, transa, diag
624 end subroutine
625 pure subroutine ctrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
626 use pm_kind , only: IK, TKG => RKS
627 implicit none
628 integer(IK) , intent(in) :: m, n, lda, ldb
629 complex(TKG), intent(inout) :: b
630 complex(TKG), intent(in) :: alpha, a
631 character , intent(in) :: side, uplo, transa, diag
632 end subroutine
633 pure subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
634 use pm_kind , only: IK, TKG => RKD
635 implicit none
636 integer(IK) , intent(in) :: m, n, lda, ldb
637 real(TKG) , intent(inout) :: b
638 real(TKG) , intent(in) :: alpha, a
639 character , intent(in) :: side, uplo, transa, diag
640 end subroutine
641 pure subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
642 use pm_kind , only: IK, TKG => RKS
643 implicit none
644 integer(IK) , intent(in) :: m, n, lda, ldb
645 real(TKG) , intent(inout) :: b
646 real(TKG) , intent(in) :: alpha, a
647 character , intent(in) :: side, uplo, transa, diag
648 end subroutine
649 end interface
650
651!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
652
672 interface blasSYRK
673 pure subroutine zsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
674 use pm_kind , only: IK, TKG => RKD
675 implicit none
676 character , intent(in) :: uplo, trans
677 integer(IK) , intent(in) :: n, k, lda, ldc
678 complex(TKG), intent(in) :: a, alpha, beta
679 complex(TKG), intent(inout) :: c
680 end subroutine
681 pure subroutine csyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
682 use pm_kind , only: IK, TKG => RKS
683 implicit none
684 character , intent(in) :: uplo, trans
685 integer(IK) , intent(in) :: n, k, lda, ldc
686 complex(TKG), intent(in) :: a, alpha, beta
687 complex(TKG), intent(inout) :: c
688 end subroutine
689 pure subroutine dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
690 use pm_kind , only: IK, TKG => RKD
691 implicit none
692 character , intent(in) :: uplo, trans
693 integer(IK) , intent(in) :: n, k, lda, ldc
694 real(TKG) , intent(in) :: a, alpha, beta
695 real(TKG) , intent(inout) :: c
696 end subroutine
697 pure subroutine ssyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
698 use pm_kind , only: IK, TKG => RKS
699 implicit none
700 character , intent(in) :: uplo, trans
701 integer(IK) , intent(in) :: n, k, lda, ldc
702 real(TKG) , intent(in) :: a, alpha, beta
703 real(TKG) , intent(inout) :: c
704 end subroutine
705 end interface
706
707!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
708
728 interface blasHERK
729 pure subroutine zherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
730 use pm_kind , only: IK, TKG => RKD
731 implicit none
732 character , intent(in) :: uplo, trans
733 integer(IK) , intent(in) :: n, k, lda, ldc
734 real(TKG) , intent(in) :: alpha, beta
735 complex(TKG), intent(in) :: a
736 complex(TKG), intent(inout) :: c
737 end subroutine
738 pure subroutine cherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
739 use pm_kind , only: IK, TKG => RKS
740 implicit none
741 character , intent(in) :: uplo, trans
742 integer(IK) , intent(in) :: n, k, lda, ldc
743 real(TKG) , intent(in) :: alpha, beta
744 complex(TKG), intent(in) :: a
745 complex(TKG), intent(inout) :: c
746 end subroutine
747 end interface
748
749!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
750
751end module pm_blas ! LCOV_EXCL_LINE
Return the result of the multiplication of two input General matrices matA and matB added to a third ...
Definition: pm_blas.F90:396
Return the result of the multiplication of a General matrix matA in Rectangular Default packing forma...
Definition: pm_blas.F90:102
Return the result of the multiplication of two input Hermitian-General/General-Hermitian matrices mat...
Definition: pm_blas.F90:360
Return the result of the multiplication of a Hermitian matrix matA in Rectangular Default packing for...
Definition: pm_blas.F90:232
Return the Rank-K Update of a Complex Hermitian matrix.
Definition: pm_blas.F90:728
Return the result of the multiplication of a Hermitian matrix matA in Linear Full Packed format with ...
Definition: pm_blas.F90:156
Return the result of the multiplication of a Symmetric matrix matA in Linear Full Packed format with ...
Definition: pm_blas.F90:194
Interchange the elements of vectors x and y.
Definition: pm_blas.F90:60
Return the result of the multiplication of two input Symmetric-General/General-Symmetric matrices mat...
Definition: pm_blas.F90:307
Return the result of the multiplication of a Symmetric matrix matA in Rectangular Default packing for...
Definition: pm_blas.F90:270
Return the Rank-K Update of a Real or Complex Symmetric or a Complex Hermitian matrix.
Definition: pm_blas.F90:672
Return the result of the triangular/general matrix multiplication, using scalar alpha,...
Definition: pm_blas.F90:563
Return the following matrix-vector products, using the vector x and triangular matrix A or its transp...
Definition: pm_blas.F90:453
Return the solution to a triangular system of equations with multiple right-hand sides,...
Definition: pm_blas.F90:616
Return the following matrix-vector products, using the vector x and triangular matrix A or its transp...
Definition: pm_blas.F90:510
This module contains a set of generic interfaces to the BLAS routines.
Definition: pm_blas.F90:35
character(*, SK), parameter MODULE_NAME
Definition: pm_blas.F90:41
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567