Mimetic Operators Library Enhanced 4.0
Loading...
Searching...
No Matches
Utils Class Reference

Utility Functions. More...

#include <utils.h>

Public Member Functions

void meshgrid (const vec &x, const vec &y, mat &X, mat &Y)
 An analog to the MATLAB 2D meshgrid operation.
 
void meshgrid (const vec &x, const vec &y, const vec &z, cube &X, cube &Y, cube &Z)
 An analog to the MATLAB 3D meshgrid operation.
 

Static Public Member Functions

static sp_mat spkron (const sp_mat &A, const sp_mat &B)
 A wrappper for implementing a sparse Kroenecker product.
 
static sp_mat spjoin_rows (const sp_mat &A, const sp_mat &B)
 An in place operation for joining two matrices by rows.
 
static sp_mat spjoin_cols (const sp_mat &A, const sp_mat &B)
 An in place operation for joining two matrices by columns.
 
static vec spsolve_eigen (const sp_mat &A, const vec &b)
 A wrappper for implementing a sparse solve using Eigen from SuperLU.
 

Detailed Description

Utility Functions.

Definition at line 28 of file utils.h.

Member Function Documentation

◆ meshgrid() [1/2]

void Utils::meshgrid ( const vec & x,
const vec & y,
const vec & z,
cube & X,
cube & Y,
cube & Z )

An analog to the MATLAB 3D meshgrid operation.

meshgrid(x,y,z,X,Y,Z) returns 3-D grid coordinates defined by the vectors x, y, and z. The grid represented by X, Y, and Z has size length(y)-by-length(x)-by-length(z).

Parameters
xa vector of x-indices
ya vector of y-indices
za vector of z-indices
Xa sparse matrix, will be filled by the function
Ya sparse matrix, will be filled by the function
Za sparse matrix, will be filled by the function

Definition at line 196 of file utils.cpp.

197 {
198 int m = x.n_elem;
199 int n = y.n_elem;
200 int o = z.n_elem;
201
202 assert(m > 0);
203 assert(n > 0);
204 assert(o > 0);
205
206 // Temporary Holder of sheet of cube
207 mat sheet(m, n, fill::zeros);
208
209 // Build X
210 vec t(n, fill::ones);
211
212 X.zeros(m, n, o);
213 Y.zeros(m, n, o);
214 Z.zeros(m, n, o);
215
216 // Sheet that repeats each slice
217 for (int ii = 0; ii < m; ++ii) {
218 sheet.row(ii) = x(ii) * t.t();
219 t.ones();
220 }
221
222 for (int kk = 0; kk < o; ++kk)
223 X.slice(kk) = sheet;
224
225 // Y Cube, repeats same sheet as well
226 for (int ii = 0; ii < m; ++ii)
227 sheet.row(ii) = y.t();
228
229 for (int kk = 0; kk < o; ++kk)
230 Y.slice(kk) = sheet;
231
232 // Z cube goes by slices each with same value
233 for (int kk = 0; kk < o; ++kk)
234 Z.slice(kk).fill(z(kk));
235}

◆ meshgrid() [2/2]

void Utils::meshgrid ( const vec & x,
const vec & y,
mat & X,
mat & Y )

An analog to the MATLAB 2D meshgrid operation.

returns 2-D grid coordinates based on the coordinates contained in vectors x and y. X is a matrix where each row is a copy of x, and Y is a matrix where each column is a copy of y. The grid represented by the coordinates X and Y has length(y) rows and length(x) columns. Key here is the rows is the y-coordinate, and the columns are the x-coordinate.

Parameters
xa vector of x-indices
ya vector of y-indices
Xa sparse matrix, will be filled by the function
Ya sparse matrix, will be filled by the function

Definition at line 173 of file utils.cpp.

173 {
174 int m = x.n_elem;
175 int n = y.n_elem;
176
177 assert(m > 0);
178 assert(n > 0);
179
180 // Build X
181 vec t(n, fill::ones);
182
183 X.zeros(n, m);
184 Y.zeros(n, m);
185
186 for (int ii = 0; ii < m; ++ii) {
187 X.col(ii) = x(ii) * t;
188 t.ones();
189 }
190
191 for (int ii = 0; ii < m; ++ii)
192 Y.col(ii) = y;
193}

◆ spjoin_cols()

sp_mat Utils::spjoin_cols ( const sp_mat & A,
const sp_mat & B )
static

An in place operation for joining two matrices by columns.

Parameters
Aa sparse matrix
Ba sparse matrix
Note
This is available in Armadillo >=8.5

Definition at line 138 of file utils.cpp.

138 {
139 sp_mat::const_iterator itA = A.begin();
140 sp_mat::const_iterator endA = A.end();
141 sp_mat::const_iterator itB = B.begin();
142 sp_mat::const_iterator endB = B.end();
143 u32 j = 0;
144
145 vec a = nonzeros(A);
146 vec b = nonzeros(B);
147
148 umat locations(2, a.n_elem + b.n_elem);
149 vec values(a.n_elem + b.n_elem);
150
151 while (itA != endA) {
152 locations(0, j) = itA.row();
153 locations(1, j) = itA.col();
154 values(j) = (*itA);
155 ++itA;
156 ++j;
157 }
158
159 while (itB != endB) {
160 locations(0, j) = itB.row() + A.n_rows;
161 locations(1, j) = itB.col();
162 values(j) = (*itB);
163 ++itB;
164 ++j;
165 }
166
167 sp_mat result(locations, values, A.n_rows + B.n_rows, A.n_cols, true);
168
169 return result;
170}

◆ spjoin_rows()

sp_mat Utils::spjoin_rows ( const sp_mat & A,
const sp_mat & B )
static

An in place operation for joining two matrices by rows.

Parameters
Aa sparse matrix
Ba sparse matrix
Note
This is available in Armadillo >8.0

Definition at line 103 of file utils.cpp.

103 {
104 sp_mat::const_iterator itA = A.begin();
105 sp_mat::const_iterator endA = A.end();
106 sp_mat::const_iterator itB = B.begin();
107 sp_mat::const_iterator endB = B.end();
108 u32 j = 0;
109
110 vec a = nonzeros(A);
111 vec b = nonzeros(B);
112
113 umat locations(2, a.n_elem + b.n_elem);
114 vec values(a.n_elem + b.n_elem);
115
116 while (itA != endA) {
117 locations(0, j) = itA.row();
118 locations(1, j) = itA.col();
119 values(j) = (*itA);
120 ++itA;
121 ++j;
122 }
123
124 while (itB != endB) {
125 locations(0, j) = itB.row();
126 locations(1, j) = itB.col() + A.n_cols;
127 values(j) = (*itB);
128 ++itB;
129 ++j;
130 }
131
132 sp_mat result(locations, values, A.n_rows, A.n_cols + B.n_cols, true);
133
134 return result;
135}

◆ spkron()

sp_mat Utils::spkron ( const sp_mat & A,
const sp_mat & B )
static

A wrappper for implementing a sparse Kroenecker product.

Parameters
Aa sparse matrix
Ba sparse matrix
Note
This is available in Armadillo >8.0

Definition at line 70 of file utils.cpp.

70 {
71 sp_mat::const_iterator itA = A.begin();
72 sp_mat::const_iterator endA = A.end();
73 sp_mat::const_iterator itB = B.begin();
74 sp_mat::const_iterator endB = B.end();
75 u32 j = 0;
76
77 vec a = nonzeros(A);
78 vec b = nonzeros(B);
79
80 umat locations(2, a.n_elem * b.n_elem);
81 vec values(a.n_elem * b.n_elem);
82
83 while (itA != endA) {
84 while (itB != endB) {
85 locations(0, j) = itA.row() * B.n_rows + itB.row();
86 locations(1, j) = itA.col() * B.n_cols + itB.col();
87 values(j) = (*itA) * (*itB);
88 ++j;
89 ++itB;
90 }
91
92 ++itA;
93 itB = B.begin();
94 }
95
96 sp_mat result(locations, values, A.n_rows * B.n_rows, A.n_cols * B.n_cols,
97 true);
98
99 return result;
100}

◆ spsolve_eigen()

static vec Utils::spsolve_eigen ( const sp_mat & A,
const vec & b )
static

A wrappper for implementing a sparse solve using Eigen from SuperLU.

Parameters
Aa sparse matrix LHS of Ax=b
ba vector for the RHS of Ax=b
Note
This function requires the EIGEN to be used when Armadillo is built

The documentation for this class was generated from the following files: