openGPMP
Open Source Mathematics Package
mtx_arr_f90.cpp
Go to the documentation of this file.
1 /*************************************************************************
2  *
3  * Project
4  * _____ _____ __ __ _____
5  * / ____| __ \| \/ | __ \
6  * ___ _ __ ___ _ __ | | __| |__) | \ / | |__) |
7  * / _ \| '_ \ / _ \ '_ \| | |_ | ___/| |\/| | ___/
8  *| (_) | |_) | __/ | | | |__| | | | | | | |
9  * \___/| .__/ \___|_| |_|\_____|_| |_| |_|_|
10  * | |
11  * |_|
12  *
13  * Copyright (C) Akiel Aries, <akiel@akiel.org>, et al.
14  *
15  * This software is licensed as described in the file LICENSE, which
16  * you should have received as part of this distribution. The terms
17  * among other details are referenced in the official documentation
18  * seen here : https://akielaries.github.io/openGPMP/ along with
19  * important files seen in this project.
20  *
21  * You may opt to use, copy, modify, merge, publish, distribute
22  * and/or sell copies of the Software, and permit persons to whom
23  * the Software is furnished to do so, under the terms of the
24  * LICENSE file. As this is an Open Source effort, all implementations
25  * must be of the same methodology.
26  *
27  *
28  *
29  * This software is distributed on an AS IS basis, WITHOUT
30  * WARRANTY OF ANY KIND, either express or implied.
31  *
32  ************************************************************************/
33 #include <cassert>
34 #include <cstddef>
35 #include <cstdint>
36 #include <iostream>
37 #include <openGPMP/linalg/mtx.hpp>
38 #include <vector>
39 
40 /************************************************************************
41  *
42  * Matrix Operations on Arrays using Fortran Subroutines
43  * in `mtx_routines.f90`
44  *
45  ************************************************************************/
46 
47 // if source C++ API is being compiled, include these Fortran refs and wrappers
48 #if defined(__GPMP_CPP_API__)
49 
50 extern "C" {
51 // Matrix add routine (FLOAT)
52 void mtx_add_routine_float_(float *A,
53  float *B,
54  float *C,
55  std::size_t *mtx_size);
56 
57 // Matrix add routine (INT)
58 void mtx_add_routine_int_(int *A, int *B, int *C, std::size_t *mtx_size);
59 
60 void mtx_mult_routine_int_(int *A,
61  int *B,
62  int *C,
63  std::size_t *rows_a,
64  std::size_t *cols_a,
65  std::size_t *cols_b);
66 }
67 
68 // C++ wrapper for Fortran mtx addition subroutine FLOAT
70  float *B,
71  float *C,
72  std::size_t mtx_size) {
73  mtx_add_routine_float_(A, B, C, &mtx_size);
74 }
75 
76 // C++ wrapper for Fortran mtx addition subroutine INT
78  int *B,
79  int *C,
80  std::size_t mtx_size) {
81  mtx_add_routine_int_(A, B, C, &mtx_size);
82 }
83 
84 // C++ wrapper for Fortran mtx multiplication subroutine INT
85 void gpmp::linalg::Mtx::mtx_mult_f90(int *A,
86  int *B,
87  int *C,
88  std::size_t rows_a,
89  std::size_t cols_a,
90  std::size_t cols_b) {
91  mtx_mult_routine_int_(A, B, C, &rows_a, &cols_a, &cols_b);
92 }
93 
94 #endif
subroutine mtx_add_routine_float_(A, B, C, mtx_size)
FORTRAN Subroutine for Matrix Addition on flattened matrices as arrays of type float32....
subroutine mtx_mult_routine_int_(A, B, C, rows_a, cols_a, cols_b)
FORTRAN Subroutine for Matrix Multiplication using Fortran intrinsics. Contains C++ wrapper function.
subroutine mtx_add_routine_int_(A, B, C, mtx_size)
FORTRAN Subroutine for Matrix Addition on flattened matrices as arrays of type int32....
list C
Definition: linalg.py:24
list A
Definition: linalg.py:22
list B
Definition: linalg.py:23
void mtx_add_f90()
Definition: mtx.cpp:333